32 #ifndef OOMPH_WAVE_FLUX_ELEMENTS_HEADER 33 #define OOMPH_WAVE_FLUX_ELEMENTS_HEADER 36 #include <oomph-lib-config.h> 40 #include "../generic/Qelements.h" 55 template <
class ELEMENT>
77 "Don't call empty constructor for LinearWaveFluxElement",
78 OOMPH_CURRENT_FUNCTION,
79 OOMPH_EXCEPTION_LOCATION);
123 const unsigned &
i)
const 132 void output(std::ostream &outfile,
const unsigned &n_plot)
142 void output(FILE* file_pt,
const unsigned &n_plot)
155 unsigned n_node =
nnode();
161 for(
unsigned i=0;
i<n_node;
i++) {test[
i] = psi[
i];}
181 (*Flux_fct_pt)(time,x,flux);
219 template<
class ELEMENT>
233 ELEMENT* elem_pt =
dynamic_cast<ELEMENT*
>(bulk_el_pt);
236 if(elem_pt->dim()==3)
245 "This flux element will not work correctly if nodes are hanging\n",
246 OOMPH_CURRENT_FUNCTION,
247 OOMPH_EXCEPTION_LOCATION);
281 "Bulk element must inherit from LinearWaveEquations.";
283 "Nodes are one dimensional, but cannot cast the bulk element to\n";
284 error_string +=
"LinearWaveEquations<1>\n.";
286 "If you desire this functionality, you must implement it yourself\n";
290 OOMPH_CURRENT_FUNCTION,
291 OOMPH_EXCEPTION_LOCATION);
311 "Bulk element must inherit from LinearWaveEquations.";
313 "Nodes are two dimensional, but cannot cast the bulk element to\n";
314 error_string +=
"LinearWaveEquations<2>\n.";
316 "If you desire this functionality, you must implement it yourself\n";
320 OOMPH_CURRENT_FUNCTION,
321 OOMPH_EXCEPTION_LOCATION);
340 "Bulk element must inherit from LinearWaveEquations.";
342 "Nodes are three dimensional, but cannot cast the bulk element to\n";
343 error_string +=
"LinearWaveEquations<3>\n.";
345 "If you desire this functionality, you must implement it yourself\n";
349 OOMPH_CURRENT_FUNCTION,
350 OOMPH_EXCEPTION_LOCATION);
363 std::ostringstream error_stream;
364 error_stream <<
"Dimension of node is " <<
Dim 365 <<
". It should be 1,2, or 3!" << std::endl;
369 OOMPH_CURRENT_FUNCTION,
370 OOMPH_EXCEPTION_LOCATION);
380 template<
class ELEMENT>
387 const unsigned n_node =
nnode();
393 Shape psif(n_node), testf(n_node);
409 for(
unsigned ipt=0;ipt<n_intpt;ipt++)
413 for(
unsigned i=0;
i<(
Dim-1);
i++)
432 for(
unsigned i=0;
i<
Dim;
i++)
434 interpolated_x[
i] = 0.0;
438 for(
unsigned l=0;l<n_node;l++)
441 for(
unsigned i=0;
i<
Dim;
i++)
454 for(
unsigned l=0;l<n_node;l++)
461 residuals[local_eqn] -= flux*testf[l]*
W;
LinearWavePrescribedFluxFctPt Flux_fct_pt
Function pointer to the (global) prescribed-flux function.
bool has_hanging_nodes() const
Return boolean to indicate if any of the element's nodes are geometrically hanging.
void broken_copy(const std::string &class_name)
Issue error message and terminate execution.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Compute the element residual vector.
virtual void output(std::ostream &outfile)
Output the element data — typically the values at the nodes in a format suitable for post-processing...
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 ...
int & face_index()
Index of the face (a number that uniquely identifies the face in the element)
A general Finite Element class.
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, when viewed as part of a compound geometric object is specified using the boundary coordinate defined by the mesh. Note: Boundary coordinates will have been set up when creating the underlying mesh, and their values will have been stored at the nodes.
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...
A class for elements that allow the imposition of an applied flux on the boundaries of LinearWave ele...
virtual unsigned u_index_lin_wave() const
Return the index at which the unknown value is stored. The default value, 0, is appropriate for singl...
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
double J_eulerian(const Vector< double > &s) const
Return the Jacobian of mapping from local to global coordinates at local position s...
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute the element's residual vector and its Jacobian matrix.
void operator=(const LinearWaveFluxElement &)
Broken assignment operator.
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
LinearWaveFluxElement(const LinearWaveFluxElement &dummy)
Broken copy constructor.
unsigned U_index_lin_wave
The index at which the unknown is stored at the nodes.
void get_flux(const double &time, const Vector< double > &x, double &flux)
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
LinearWaveFluxElement()
Broken empty constructor.
void output(FILE *file_pt)
double & time()
Return the current value of the continuous time.
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
void output(std::ostream &outfile)
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(FILE *file_pt, const unsigned &n_plot)
Output function – forward to broken version in FiniteElement until somebody decides what exactly the...
unsigned Dim
The spatial dimension of the problem.
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 ...
void broken_assign(const std::string &class_name)
Issue error message and terminate execution.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
Time *const & time_pt() const
Access function for the pointer to time (const version)
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
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...
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!
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
The "global" intrinsic coordinate of the element when viewed as part of a geometric object should be ...
unsigned nnode() const
Return the number of nodes.
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...
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
void output(std::ostream &outfile, const unsigned &n_plot)
Output function – forward to broken version in FiniteElement until somebody decides what exactly the...
LinearWavePrescribedFluxFctPt & flux_fct_pt()
Access function for the prescribed-flux function pointer.
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...