159 for (
unsigned i = 0;
i <
dim;
i++)
175 for (
unsigned i = 0;
i <
dim;
i++)
220 for (
unsigned i = 0;
i < n_flux;
i++)
428 <<
"Coupling data between elements not included in jacobian\n"
429 <<
"You should call DGMesh::setup_face_neighbour_info(true) to "
431 <<
"that this information is included in the jacobian\n";
622 <<
"Cannot use a discontinuous formulation for the mass matrix when\n"
623 <<
"there are external data.\n "
624 <<
"Do not call Problem::enable_discontinuous_formulation()\n";
640 for (
unsigned n = 0;
n <
n_dof;
n++)
676 const int& face_index,
705 ->neighbour_face_pt(0)
710 ->neighbour_face_pt(0)
726 error_stream <<
"Slope limiting is not implemented for this dimension: "
755 else if (
args[0] > 0.0)
765 double min = std::fabs(
args[0]);
768 for (
unsigned i = 1;
i <
n_arg;
i++)
810 if (std::fabs(
args[0]) < this->
M * h * h)
826 const double tol = 1.0e-16;
832 const double h =
x_r -
x_l;
833 const double x0 = 0.5 * (
x_l +
x_r);
A Base class for DGElements.
DGMesh * DG_mesh_pt
Pointer to Mesh, which will be responsible for the neighbour finding.
void pre_compute_mass_matrix()
Function that computes and stores the (inverse) mass matrix.
void slope_limit(SlopeLimiter *const &slope_limiter_pt)
Limit the slope within the element.
virtual unsigned required_nflux()
Set the number of flux components.
bool Mass_matrix_has_been_computed
Boolean flag to indicate whether the mass matrix has been computed.
DGFaceElement * face_element_pt(const unsigned &i)
Access function for the faces.
bool Mass_matrix_reuse_is_enabled
Boolean flag to indicate whether to reuse the mass matrix.
void get_neighbouring_face_and_local_coordinate(const int &face_index, const Vector< double > &s, FaceElement *&face_element_pt, Vector< double > &s_face)
Return the neighbour info.
virtual void get_inverse_mass_matrix_times_residuals(Vector< double > &minv_res)
Function that returns the current value of the residuals multiplied by the inverse mass matrix (virtu...
DenseDoubleMatrix * M_pt
Pointer to storage for a mass matrix that can be recycled if desired.
Base class for Discontinuous Galerkin Faces. These are responsible for calculating the normal fluxes ...
virtual void numerical_flux(const Vector< double > &n_out, const Vector< double > &u_int, const Vector< double > &u_ext, Vector< double > &flux)
Calculate the normal flux, which is the dot product of our approximation to the flux with the outer u...
virtual void get_interpolation_data(Vector< Data * > &interpolation_data)
Get the data that are used to interpolate the unkowns in the element. These must be returned in order...
void setup_neighbour_info(const bool &add_neighbour_data_to_bulk)
Setup information from the neighbouring face, return a set of nodes (as data) in the neighour....
void report_info()
Output information about the present element and its neighbour.
Vector< FaceElement * > Neighbour_face_pt
Vector of neighbouring face elements at the integration points.
void add_flux_contributions(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Add the contribution from integrating the numerical flux.
Vector< Vector< unsigned > > Neighbour_external_data
Vector of the vectors that will store the number of the bulk external data that correspond to the dof...
virtual void dnumerical_flux_du(const Vector< double > &n_out, const Vector< double > &u_int, const Vector< double > &u_ext, DenseMatrix< double > &dflux_du_int, DenseMatrix< double > &dflux_du_ext)
Calculate the derivative of the normal flux, which is the dot product of our approximation to the flu...
virtual void interpolated_u(const Vector< double > &s, Vector< double > &f)
Return the interpolated values of the unknown fluxes.
virtual void numerical_flux_at_knot(const unsigned &ipt, const Shape &psi, Vector< double > &flux, DenseMatrix< double > &dflux_du_int, DenseMatrix< double > &dflux_du_ext, unsigned flag)
Calculate the normal numerical flux at the integration point. This is the most general interface that...
Vector< Vector< double > > Neighbour_local_coordinate
Vector of neighbouring local coordinates at the integration points.
virtual unsigned flux_index(const unsigned &i) const
Return the index at which the i-th unknown flux is stored.
virtual unsigned required_nflux()
Set the number of flux components.
virtual void neighbour_finder(FiniteElement *const &bulk_element_pt, const int &face_index, const Vector< double > &s_bulk, FaceElement *&face_element_pt, Vector< double > &s_face)
static double FaceTolerance
void set_value(const unsigned &i, const double &value_)
Set the i-th stored data value to specified value. The only reason that we require an explicit set fu...
Class of matrices containing doubles, and stored as a DenseMatrix<double>, but with solving functiona...
virtual void lubksub(DoubleVector &rhs)
LU backsubstitution.
virtual void ludecompose()
LU decomposition using DenseLU (default linea solver)
void initialise(const T &val)
Initialize all values in the matrix to val.
void resize(const unsigned long &n)
Resize to a square nxn matrix; any values already present will be transfered.
FaceElements are elements that coincide with the faces of higher-dimensional "bulk" elements....
int & face_index()
Index of the face (a number that uniquely identifies the face in the element)
void outer_unit_normal(const Vector< double > &s, Vector< double > &unit_normal) const
Compute outer unit normal at the specified local coordinate.
double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s. Overloaded to get information from bulk...
unsigned & bulk_node_number(const unsigned &n)
Return the bulk node number that corresponds to the n-th local node number.
FiniteElement *& bulk_element_pt()
Pointer to higher-dimensional "bulk" element.
double J_eulerian_at_knot(const unsigned &ipt) const
Return the Jacobian of the mapping from local to global coordinates at the ipt-th integration point O...
void get_local_coordinate_in_bulk(const Vector< double > &s, Vector< double > &s_bulk) const
Calculate the vector of local coordinate in the bulk element given the local coordinates in this Face...
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
double nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n. Produces suitably interpolated values for hanging nodes...
double size() const
Calculate the size of the element (length, area, volume,...) in Eulerian computational coordinates....
virtual void shape(const Vector< double > &s, Shape &psi) const =0
Calculate the geometric shape functions at local coordinate s. This function must be overloaded for e...
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
unsigned nnode() const
Return the number of nodes.
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
unsigned nodal_dimension() const
Return the required Eulerian dimension of the nodes in this element.
virtual void shape_at_knot(const unsigned &ipt, Shape &psi) const
Return the geometric shape function at the ipt-th integration point.
unsigned nexternal_data() const
Return the number of external data objects.
virtual void fill_in_contribution_to_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &mass_matrix)
Add the elemental contribution to the mass matrix matrix. and the residuals vector....
static double Default_fd_jacobian_step
Double used for the default finite difference step in elemental jacobian calculations.
virtual void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the elemental contribution to the residuals vector. Note that this function will NOT initialise t...
unsigned ndof() const
Return the number of equations/dofs in the element.
int external_local_eqn(const unsigned &i, const unsigned &j)
Return the local equation number corresponding to the j-th value stored at the i-th external data.
unsigned add_external_data(Data *const &data_pt, const bool &fd=true)
Add a (pointer to an) external data object to the element and return its index (i....
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
double minmod(Vector< double > &args)
The basic minmod function.
bool MUSCL
Boolean flag to indicate a MUSCL or straight min-mod limiter.
double minmodB(Vector< double > &args, const double &h)
The modified minmod function.
void limit(const unsigned &i, const Vector< DGElement * > &required_element_pt)
The limit function.
double M
Constant that is used in the modified min-mod function to provide better properties at extrema.
double & x(const unsigned &i)
Return the i-th nodal coordinate.
double value(const unsigned &i) const
Return i-th value (dofs or pinned) at this node either directly or via hanging node representation....
An OomphLibError object which should be thrown when an run-time error is encountered....
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Base class for slope limiters.
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
OomphInfo oomph_info
Single (global) instantiation of the OomphInfo object – this is used throughout the library as a "rep...