30#ifndef OOMPH_UNSTEADY_HEAT_FLUX_ELEMENTS_HEADER
31#define OOMPH_UNSTEADY_HEAT_FLUX_ELEMENTS_HEADER
35#include <oomph-lib-config.h>
53 template<
class ELEMENT>
114 const unsigned&
i)
const
184 (*Flux_fct_pt)(time, x, flux);
214 template<
class ELEMENT>
239 throw OomphLibError(
"This flux element will not work correctly if "
240 "nodes are hanging\n",
276 "Bulk element must inherit from UnsteadyHeatEquations.";
278 "Nodes are one dimensional, but cannot cast the bulk element to\n";
280 error_string +=
"If you desire this functionality, you must "
281 "implement it yourself\n";
304 "Bulk element must inherit from UnsteadyHeatEquations.";
306 "Nodes are two dimensional, but cannot cast the bulk element to\n";
308 error_string +=
"If you desire this functionality, you must "
309 "implement it yourself\n";
331 "Bulk element must inherit from UnsteadyHeatEquations.";
332 error_string +=
"Nodes are three dimensional, but cannot cast the "
335 error_string +=
"If you desire this functionality, you must "
336 "implement it yourself\n";
353 <<
". It should be 1,2, or 3!" << std::endl;
365 template<
class ELEMENT>
389 const unsigned u_index_ust_heat = U_index_ust_heat;
396 for (
unsigned i = 0;
i < (Dim - 1);
i++)
415 for (
unsigned i = 0;
i < Dim;
i++)
424 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 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.
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.
A class for elements that allow the imposition of an applied flux on the boundaries of UnsteadyHeat e...
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 ...
UnsteadyHeatPrescribedFluxFctPt & flux_fct_pt()
Broken assignment operator.
unsigned Dim
The spatial dimension of the problem.
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.
UnsteadyHeatPrescribedFluxFctPt Flux_fct_pt
Function pointer to the (global) prescribed-flux function.
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 elem...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Compute the element residual vector.
void output(FILE *file_pt)
C-style output function – forward to broken version in FiniteElement until somebody decides what exac...
UnsteadyHeatFluxElement(FiniteElement *const &bulk_el_pt, const int &face_index)
Constructor, takes the pointer to the "bulk" element and the index of the face to be created.
unsigned U_index_ust_heat
The index at which the unknown is stored at the nodes.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function – forward to broken version in FiniteElement until somebody decides what exac...
void(* UnsteadyHeatPrescribedFluxFctPt)(const double &time, const Vector< double > &x, double &flux)
Function pointer to the prescribed-flux function fct(x,f(x)) – x is a Vector!
UnsteadyHeatFluxElement(const UnsteadyHeatFluxElement &dummy)=delete
Broken copy constructor.
void output(std::ostream &outfile)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
void fill_in_generic_residual_contribution_ust_heat_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, const unsigned &n_plot)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute the element's residual vector and its Jacobian matrix.
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).