28#ifndef OOMPH_WAVE_FLUX_ELEMENTS_HEADER
29#define OOMPH_WAVE_FLUX_ELEMENTS_HEADER
32#include <oomph-lib-config.h>
50 template<
class ELEMENT>
71 "Don't call empty constructor for LinearWaveFluxElement",
116 const unsigned&
i)
const
188 (*Flux_fct_pt)(time, x, flux);
219 template<
class ELEMENT>
244 throw OomphLibError(
"This flux element will not work correctly if "
245 "nodes are hanging\n",
281 "Bulk element must inherit from LinearWaveEquations.";
283 "Nodes are one dimensional, but cannot cast the bulk element to\n";
285 error_string +=
"If you desire this functionality, you must "
286 "implement it yourself\n";
309 "Bulk element must inherit from LinearWaveEquations.";
311 "Nodes are two dimensional, but cannot cast the bulk element to\n";
313 error_string +=
"If you desire this functionality, you must "
314 "implement it yourself\n";
336 "Bulk element must inherit from LinearWaveEquations.";
337 error_string +=
"Nodes are three dimensional, but cannot cast the "
340 error_string +=
"If you desire this functionality, you must "
341 "implement it yourself\n";
358 <<
". It should be 1,2, or 3!" << std::endl;
370 template<
class ELEMENT>
394 const unsigned u_index_lin_wave = U_index_lin_wave;
401 for (
unsigned i = 0;
i < (Dim - 1);
i++)
420 for (
unsigned i = 0;
i < Dim;
i++)
429 for (
unsigned i = 0;
i < Dim;
i++)
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: .
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
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)
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
In a FaceElement, the "global" intrinsic coordinate of the element along the boundary,...
double J_eulerian(const Vector< double > &s) const
Return the Jacobian of mapping from local to global coordinates at local position s....
FaceGeometry class definition: This policy class is used to allow construction of face elements that ...
A general Finite Element class.
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
virtual void output(std::ostream &outfile)
Output the element data — typically the values at the nodes in a format suitable for post-processing.
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.
double nodal_position(const unsigned &n, const unsigned &i) const
Return the i-th coordinate at local node n. If the node is hanging, the appropriate interpolation is ...
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
virtual void build_face_element(const int &face_index, FaceElement *face_element_pt)
Function for building a lower dimensional FaceElement on the specified face of the FiniteElement....
bool has_hanging_nodes() const
Return boolean to indicate if any of the element's nodes are geometrically hanging.
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
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.
A class for elements that allow the imposition of an applied flux on the boundaries of LinearWave ele...
void output(FILE *file_pt, const unsigned &n_plot)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Compute the element residual vector.
void(* LinearWavePrescribedFluxFctPt)(const double &time, const Vector< double > &x, double &flux)
Function pointer to the prescribed-flux function fct(x,f(x)) – x is a Vector!
LinearWavePrescribedFluxFctPt & flux_fct_pt()
Access function for the prescribed-flux function pointer.
unsigned U_index_lin_wave
The index at which the unknown is stored at the nodes.
LinearWaveFluxElement(const LinearWaveFluxElement &dummy)=delete
Broken copy constructor.
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Specify the value of nodal zeta from the face geometry The "global" intrinsic coordinate of the eleme...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute the element's residual vector and its Jacobian matrix.
unsigned Dim
The spatial dimension of the problem.
void output(FILE *file_pt)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
void operator=(const LinearWaveFluxElement &)=delete
Broken assignment operator.
LinearWavePrescribedFluxFctPt Flux_fct_pt
Function pointer to the (global) prescribed-flux function.
void get_flux(const double &time, const Vector< double > &x, double &flux)
Function to calculate the prescribed flux at a given spatial position and at a gien time.
LinearWaveFluxElement()
Broken empty constructor.
void fill_in_generic_residual_contribution_lin_wave_flux(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Compute the element residual vector. flag=1(or 0): do (or don't) compute the Jacobian as well.
void output(std::ostream &outfile)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
void output(std::ostream &outfile, const unsigned &n_plot)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
double shape_and_test(const Vector< double > &s, Shape &psi, Shape &test) const
Function to compute the shape and test functions and to return the Jacobian of mapping between local ...
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
An OomphLibError object which should be thrown when an run-time error is encountered....
RefineableElements are FiniteElements that may be subdivided into children to provide a better local ...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
TAdvectionDiffusionReactionElement()
Constructor: Call constructors for TElement and AdvectionDiffusionReaction equations.
Time *const & time_pt() const
Access function for the pointer to time (const version)
double & time()
Return the current value of the continuous time.
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).