28 #ifndef OOMPH_ADV_DIFF_FLUX_ELEMENTS_HEADER
29 #define OOMPH_ADV_DIFF_FLUX_ELEMENTS_HEADER
33 #include <oomph-lib-config.h>
37 #include "../generic/Qelements.h"
52 template<
class ELEMENT>
73 "Don't call empty constructor for AdvectionDiffusionFluxElement",
74 OOMPH_CURRENT_FUNCTION,
75 OOMPH_EXCEPTION_LOCATION);
109 residuals, jacobian, 1);
119 const unsigned&
i)
const
133 void output(std::ostream& outfile,
const unsigned& nplot)
148 unsigned n_node =
nnode();
154 for (
unsigned i = 0;
i < n_node;
i++)
172 unsigned n_node =
nnode();
178 for (
unsigned i = 0;
i < n_node;
i++)
200 (*Flux_fct_pt)(x, flux);
231 template<
class ELEMENT>
244 ELEMENT* elem_pt =
dynamic_cast<ELEMENT*
>(bulk_el_pt);
246 if (elem_pt->dim() == 3)
255 throw OomphLibError(
"This flux element will not work correctly if "
256 "nodes are hanging\n",
257 OOMPH_CURRENT_FUNCTION,
258 OOMPH_EXCEPTION_LOCATION);
293 "Bulk element must inherit from AdvectionDiffusionEquations.";
295 "Nodes are one dimensional, but cannot cast the bulk element to\n";
296 error_string +=
"AdvectionDiffusionEquations<1>\n.";
297 error_string +=
"If you desire this functionality, you must "
298 "implement it yourself\n";
301 error_string, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
321 "Bulk element must inherit from AdvectionDiffusionEquations.";
323 "Nodes are two dimensional, but cannot cast the bulk element to\n";
324 error_string +=
"AdvectionDiffusionEquations<2>\n.";
325 error_string +=
"If you desire this functionality, you must "
326 "implement it yourself\n";
329 error_string, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
348 "Bulk element must inherit from AdvectionDiffusionEquations.";
349 error_string +=
"Nodes are three dimensional, but cannot cast the "
351 error_string +=
"AdvectionDiffusionEquations<3>\n.";
352 error_string +=
"If you desire this functionality, you must "
353 "implement it yourself\n";
356 error_string, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
368 std::ostringstream error_stream;
369 error_stream <<
"Dimension of node is " <<
Dim
370 <<
". It should be 1,2, or 3!" << std::endl;
373 error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
382 template<
class ELEMENT>
388 const unsigned n_node = nnode();
391 Shape psif(n_node), testf(n_node);
394 const unsigned n_intpt = integral_pt()->nweight();
404 const unsigned u_index_adv_diff = U_index_adv_diff;
409 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
412 for (
unsigned i = 0;
i < (
Dim - 1);
i++)
414 s[
i] = integral_pt()->knot(ipt,
i);
418 double w = integral_pt()->weight(ipt);
422 double J = shape_and_test(
s, psif, testf);
431 for (
unsigned l = 0; l < n_node; l++)
434 for (
unsigned i = 0;
i <
Dim;
i++)
436 interpolated_x[
i] += nodal_position(l,
i) * psif[l];
442 get_flux(interpolated_x, flux);
447 for (
unsigned l = 0; l < n_node; l++)
450 local_eqn = nodal_local_eqn(l, u_index_adv_diff);
455 residuals[local_eqn] += flux * testf[l] *
W;
A class for all elements that solve the Advection Diffusion equations using isoparametric elements.
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...
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
void output(std::ostream &outfile)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
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 ...
AdvectionDiffusionFluxElement(const AdvectionDiffusionFluxElement &dummy)=delete
Broken copy constructor.
void operator=(const AdvectionDiffusionFluxElement &)=delete
Broken assignment operator.
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.
void output(std::ostream &outfile, const unsigned &nplot)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
unsigned U_index_adv_diff
The index at which the unknown is stored at the nodes.
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 fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector.
AdvectionDiffusionPrescribedFluxFctPt Flux_fct_pt
Function pointer to the (global) prescribed-flux function.
AdvectionDiffusionPrescribedFluxFctPt & flux_fct_pt()
Access function for the prescribed-flux function pointer.
unsigned Dim
The spatial dimension of the problem.
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 get_flux(const Vector< double > &x, double &flux)
Function to calculate the prescribed flux at a given spatial position.
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...
void(* AdvectionDiffusionPrescribedFluxFctPt)(const Vector< double > &x, double &flux)
Function pointer to the prescribed-flux function fct(x,f(x)) – x is a Vector!
AdvectionDiffusionFluxElement()
Broken empty constructor.
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....
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...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
A general Finite Element class.
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
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...
unsigned nnode() const
Return the number of nodes.
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....
virtual void shape_at_knot(const unsigned &ipt, Shape &psi) const
Return the geometric shape function at the ipt-th integration point.
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...
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...
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
unsigned Dim
Dimension of zeta tuples (set by get_dim_helper) – needed because we store the scalar coordinates in ...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...