28#ifndef OOMPH_POISSON_FLUX_ELEMENTS_HEADER
29#define OOMPH_POISSON_FLUX_ELEMENTS_HEADER
34#include <oomph-lib-config.h>
48 template<
class ELEMENT>
65 throw OomphLibError(
"Don't call empty constructor for PoissonFluxElement",
83 const unsigned&
i)
const
118 const unsigned n_plot = 5;
142 for (
unsigned i = 0;
i <
el_dim + 1;
i++)
234 (*Flux_fct_pt)(x, flux);
245 const unsigned&
flag);
270 template<
class ELEMENT>
294 throw OomphLibError(
"This flux element will not work correctly if "
295 "nodes are hanging\n",
331 "Bulk element must inherit from PoissonEquations.";
333 "Nodes are one dimensional, but cannot cast the bulk element to\n";
335 error_string +=
"If you desire this functionality, you must "
336 "implement it yourself\n";
359 "Bulk element must inherit from PoissonEquations.";
361 "Nodes are two dimensional, but cannot cast the bulk element to\n";
363 error_string +=
"If you desire this functionality, you must "
364 "implement it yourself\n";
386 "Bulk element must inherit from PoissonEquations.";
387 error_string +=
"Nodes are three dimensional, but cannot cast the "
390 error_string +=
"If you desire this functionality, you must "
391 "implement it yourself\n";
408 <<
". It should be 1,2, or 3!" << std::endl;
420 template<
class ELEMENT>
425 const unsigned&
flag)
443 const unsigned u_index_poisson = U_index_poisson;
450 for (
unsigned i = 0;
i < (Dim - 1);
i++)
472 for (
unsigned i = 0;
i < Dim;
i++)
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: .
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 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...
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...
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 std::string tecplot_zone_string(const unsigned &nplot) const
Return string for tecplot zone header (when plotting nplot points in each "coordinate direction")
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.
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Get cector of local coordinates of plot point i (when plotting nplot points in each "coordinate direc...
virtual unsigned nplot_points(const unsigned &nplot) const
Return total number of plot points (when plotting nplot points in each "coordinate direction")
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....
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Add tecplot zone "footer" to output stream (when plotting nplot points in each "coordinate direction"...
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...
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.
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
An OomphLibError object which should be thrown when an run-time error is encountered....
A class for elements that allow the imposition of an applied flux on the boundaries of Poisson elemen...
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 output(FILE *file_pt)
C-style output function – forward to broken version in FiniteElement until somebody decides what exac...
PoissonPrescribedFluxFctPt & flux_fct_pt()
Access function for the prescribed-flux function pointer.
PoissonFluxElement()
Broken empty constructor.
PoissonFluxElement(const PoissonFluxElement &dummy)=delete
Broken copy constructor.
unsigned Dim
The spatial dimension of the problem.
void operator=(const PoissonFluxElement &)=delete
Broken assignment operator.
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 ...
void get_flux(const Vector< double > &x, double &flux)
Function to calculate the prescribed flux at a given spatial position.
void(* PoissonPrescribedFluxFctPt)(const Vector< double > &x, double &flux)
Function pointer to the prescribed-flux function fct(x,f(x)) – x is a Vector!
void fill_in_generic_residual_contribution_poisson_flux(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Add the element's contribution to its residual vector. flag=1(or 0): do (or don't) compute the contri...
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.
unsigned U_index_poisson
The index at which the unknown is stored at the nodes.
void output(std::ostream &outfile)
Output function.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function – forward to broken version in FiniteElement until somebody decides what exac...
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...
PoissonPrescribedFluxFctPt Flux_fct_pt
Function pointer to the (global) prescribed-flux function.
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.
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).