32 #ifndef OOMPH_ADV_DIFF_FLUX_ELEMENTS_HEADER 33 #define OOMPH_ADV_DIFF_FLUX_ELEMENTS_HEADER 37 #include <oomph-lib-config.h> 41 #include "../generic/Qelements.h" 59 template <
class ELEMENT>
83 "Don't call empty constructor for AdvectionDiffusionFluxElement",
84 OOMPH_CURRENT_FUNCTION,
85 OOMPH_EXCEPTION_LOCATION);
129 const unsigned &
i)
const 138 void output(std::ostream &outfile,
const unsigned &nplot)
151 unsigned n_node =
nnode();
157 for(
unsigned i=0;
i<n_node;
i++) {test[
i] = psi[
i];}
172 unsigned n_node =
nnode();
178 for(
unsigned i=0;
i<n_node;
i++) {test[
i] = psi[
i];}
198 (*Flux_fct_pt)(x,flux);
238 template<
class ELEMENT>
252 ELEMENT* elem_pt =
dynamic_cast<ELEMENT*
>(bulk_el_pt);
254 if(elem_pt->dim()==3)
263 "This flux element will not work correctly if nodes are hanging\n",
264 OOMPH_CURRENT_FUNCTION,
265 OOMPH_EXCEPTION_LOCATION);
300 "Bulk element must inherit from AdvectionDiffusionEquations.";
302 "Nodes are one dimensional, but cannot cast the bulk element to\n";
303 error_string +=
"AdvectionDiffusionEquations<1>\n.";
305 "If you desire this functionality, you must implement it yourself\n";
309 OOMPH_CURRENT_FUNCTION,
310 OOMPH_EXCEPTION_LOCATION);
330 "Bulk element must inherit from AdvectionDiffusionEquations.";
332 "Nodes are two dimensional, but cannot cast the bulk element to\n";
333 error_string +=
"AdvectionDiffusionEquations<2>\n.";
335 "If you desire this functionality, you must implement it yourself\n";
339 OOMPH_CURRENT_FUNCTION,
340 OOMPH_EXCEPTION_LOCATION);
359 "Bulk element must inherit from AdvectionDiffusionEquations.";
361 "Nodes are three dimensional, but cannot cast the bulk element to\n";
362 error_string +=
"AdvectionDiffusionEquations<3>\n.";
364 "If you desire this functionality, you must implement it yourself\n";
368 OOMPH_CURRENT_FUNCTION,
369 OOMPH_EXCEPTION_LOCATION);
382 std::ostringstream error_stream;
383 error_stream <<
"Dimension of node is " <<
Dim 384 <<
". It should be 1,2, or 3!" << std::endl;
388 OOMPH_CURRENT_FUNCTION,
389 OOMPH_EXCEPTION_LOCATION);
398 template<
class ELEMENT>
406 const unsigned n_node =
nnode();
409 Shape psif(n_node), testf(n_node);
427 for(
unsigned ipt=0;ipt<n_intpt;ipt++)
447 for(
unsigned l=0;l<n_node;l++)
450 for(
unsigned i=0;
i<
Dim;
i++)
463 for(
unsigned l=0;l<n_node;l++)
471 residuals[local_eqn] += flux*testf[l]*
W;
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...
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.
AdvectionDiffusionFluxElement()
Broken empty constructor.
unsigned Dim
The spatial dimension of the problem.
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 ...
void output(std::ostream &outfile, const unsigned &nplot)
Output function – forward to broken version in FiniteElement until somebody decides what exactly the...
int & face_index()
Index of the face (a number that uniquely identifies the face in the element)
AdvectionDiffusionPrescribedFluxFctPt Flux_fct_pt
Function pointer to the (global) prescribed-flux function.
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 ...
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...
void get_flux(const Vector< double > &x, double &flux)
Function to calculate the prescribed flux at a given spatial position.
void(* AdvectionDiffusionPrescribedFluxFctPt)(const Vector< double > &x, double &flux)
Function pointer to the prescribed-flux function fct(x,f(x)) – x is a Vector!
void output(std::ostream &outfile)
Output function – forward to broken version in FiniteElement until somebody decides what exactly the...
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
A class for all elements that solve the Advection Diffusion equations using isoparametric elements...
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.
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
virtual void shape_at_knot(const unsigned &ipt, Shape &psi) const
Return the geometric shape function at the ipt-th integration point.
AdvectionDiffusionFluxElement(const AdvectionDiffusionFluxElement &dummy)
Broken copy constructor.
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector.
void fill_in_generic_residual_contribution_adv_diff_flux(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Add the element's contribution to its residual vector. flag=1(or 0): do (or don't) compute the Jacobi...
double shape_and_test_at_knot(const unsigned &ipt, Shape &psi, Shape &test) const
Function to compute the shape and test functions and to return the Jacobian of mapping between local ...
virtual unsigned u_index_adv_diff() const
Return the index at which the unknown value is stored. The default value, 0, is appropriate for singl...
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
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 ...
AdvectionDiffusionPrescribedFluxFctPt & flux_fct_pt()
Access function for the prescribed-flux function pointer.
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.
A class for elements that allow the imposition of an applied flux on the boundaries of Advection Diff...
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
void operator=(const AdvectionDiffusionFluxElement &)
Broken assignment operator.
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...
unsigned U_index_adv_diff
The index at which the unknown is stored at the nodes.
unsigned nnode() const
Return the number of nodes.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Add the element's contribution to its residual vector and its Jacobian matrix.
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...
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...