27#ifndef OOMPH_STEADY_AXISYM_ADV_DIFF_ELEMENTS_HEADER
28#define OOMPH_STEADY_AXISYM_ADV_DIFF_ELEMENTS_HEADER
33#include <oomph-lib-config.h>
125 const unsigned&
nplot,
164 const double&
pe()
const
181 double& source)
const
191 (*Source_fct_pt)(x, source);
208 for (
unsigned i = 0;
i < 2;
i++)
216 (*Wind_fct_pt)(x, wind);
237 for (
unsigned j = 0;
j < 2;
j++)
246 for (
unsigned j = 0;
j < 2;
j++)
294 double interpolated_u = 0.0;
302 return (interpolated_u);
342 global_eqn_number.resize(
n_u_dof, 0);
423 template<
unsigned NNODE_1D>
425 :
public virtual QElement<2, NNODE_1D>,
525 template<
unsigned NNODE_1D>
526 double QSteadyAxisymAdvectionDiffusionElement<
541 for (
unsigned j = 0;
j < 2;
j++)
559 template<
unsigned NNODE_1D>
561 NNODE_1D>::dshape_and_dtest_eulerian_at_knot_adv_diff(
const unsigned&
ipt,
582 template<
unsigned NNODE_1D>
584 :
public virtual QElement<1, NNODE_1D>
605 template<
class ELEMENT>
632 "SteadyAxisymAdvectionDiffusionFluxElement",
686 const unsigned&
i)
const
767 (*Beta_fct_pt)(x, beta);
783 (*Alpha_fct_pt)(x, alpha);
816 template<
class ELEMENT>
819 const int& face_index)
843 throw OomphLibError(
"This flux element will not work correctly if "
844 "nodes are hanging\n",
871 std::string
error_string =
"Bulk element must inherit from "
872 "SteadyAxisymAdvectionDiffusionEquations.";
874 "Nodes are two dimensional, but cannot cast the bulk element to\n";
875 error_string +=
"SteadyAxisymAdvectionDiffusionEquations<2>\n.";
877 "If you desire this functionality, you must implement it yourself\n";
895 template<
class ELEMENT>
904 const unsigned u_index_axisym_adv_diff = U_index_adv_diff;
924 for (
unsigned i = 0;
i < 1;
i++)
941 double interpolated_u = 0.0;
951 for (
unsigned i = 0;
i < 2;
i++)
980 r * (beta - alpha * interpolated_u) *
testf(
l) * W;
984 if ((
flag) && (alpha != 0.0))
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
long & eqn_number(const unsigned &i)
Return the equation number of the i-th stored variable.
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...
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
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)
double nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n. Produces suitably interpolated values for hanging nodes...
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...
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.
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
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 ...
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Compute the geometric shape functions and also first derivatives w.r.t. global coordinates at local c...
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....
double raw_nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n but do NOT take hanging nodes into account.
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.
An OomphLibError object which should be thrown when an run-time error is encountered....
QSteadyAxisymAdvectionDiffusionElement elements are linear/quadrilateral/brick-shaped Axisymmetric Ad...
QSteadyAxisymAdvectionDiffusionElement(const QSteadyAxisymAdvectionDiffusionElement< NNODE_1D > &dummy)=delete
Broken copy constructor.
double dshape_and_dtest_eulerian_adv_diff(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Shape, test functions & derivs. w.r.t. to global coords. Return Jacobian.
unsigned required_nvalue(const unsigned &n) const
Broken assignment operator.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: r,z,u at n_plot^2 plot points.
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output function for an exact solution: r,z,u_exact at n_plot^2 plot points.
void output(FILE *file_pt)
C-style output function: r,z,u.
double dshape_and_dtest_eulerian_at_knot_adv_diff(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Shape, test functions & derivs. w.r.t. to global coords. at integration point ipt....
static const unsigned Initial_Nvalue
Static array of ints to hold number of variables at nodes: Initial_Nvalue[n].
void output(std::ostream &outfile)
Output function: r,z,u.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function: r,z,u at n_plot^2 plot points.
QSteadyAxisymAdvectionDiffusionElement()
Constructor: Call constructors for QElement and Advection Diffusion equations.
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...
A class for all elements that solve the Steady Axisymmetric Advection Diffusion equations using isopa...
virtual double dshape_and_dtest_eulerian_adv_diff(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Shape/test functions and derivs w.r.t. to global coords at local coord. s; return Jacobian of mapping...
SteadyAxisymAdvectionDiffusionSourceFctPt source_fct_pt() const
Access function: Pointer to source function. Const version.
SteadyAxisymAdvectionDiffusionEquations(const SteadyAxisymAdvectionDiffusionEquations &dummy)=delete
Broken copy constructor.
virtual void dinterpolated_u_adv_diff_ddata(const Vector< double > &s, Vector< double > &du_ddata, Vector< unsigned > &global_eqn_number)
Return derivative of u at point s with respect to all data that can affect its value....
SteadyAxisymAdvectionDiffusionSourceFctPt Source_fct_pt
Pointer to source function:
virtual double dshape_and_dtest_eulerian_at_knot_adv_diff(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of ...
void(* SteadyAxisymAdvectionDiffusionSourceFctPt)(const Vector< double > &x, double &f)
Function pointer to source function fct(x,f(x)) – x is a Vector!
static double Default_peclet_number
Static default value for the Peclet number.
virtual void get_wind_axisym_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Get wind at (Eulerian) position x and/or local coordinate s. This function is virtual to allow overlo...
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
unsigned self_test()
Self-test: Return 0 for OK.
virtual void get_source_axisym_adv_diff(const unsigned &ipt, const Vector< double > &x, double &source) const
Get source term at (Eulerian) position x. This function is virtual to allow overloading in multi-phys...
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: r,z,u_exact at nplot^2 plot points.
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: .
const double & pe() const
Peclet number.
SteadyAxisymAdvectionDiffusionWindFctPt Wind_fct_pt
Pointer to wind function:
SteadyAxisymAdvectionDiffusionEquations()
Constructor: Initialise the Source_fct_pt and Wind_fct_pt to null and set (pointer to) Peclet number ...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector (wrapper)
double * Pe_pt
Pointer to global Peclet number.
SteadyAxisymAdvectionDiffusionWindFctPt & wind_fct_pt()
Access function: Pointer to wind function.
void output(FILE *file_pt)
C_style output with default number of plot points.
SteadyAxisymAdvectionDiffusionWindFctPt wind_fct_pt() const
Access function: Pointer to wind function. Const version.
double *& pe_pt()
Pointer to Peclet number.
double interpolated_u_adv_diff(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
virtual void fill_in_generic_residual_contribution_adv_diff(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Add the element's contribution to its residual vector only (if flag=and/or element Jacobian matrix.
virtual unsigned u_index_axisym_adv_diff() const
Broken assignment operator.
void output(std::ostream &outfile)
Output with default number of plot points.
SteadyAxisymAdvectionDiffusionSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
void(* SteadyAxisymAdvectionDiffusionWindFctPt)(const Vector< double > &x, Vector< double > &wind)
Function pointer to wind function fct(x,w(x)) – x is a Vector!
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Add the element's contribution to its residual vector and the element Jacobian matrix (wrapper)
A class for elements that allow the imposition of an applied Robin boundary condition on the boundari...
void get_alpha(const Vector< double > &x, double &alpha)
Function to calculate the prescribed alpha at a given spatial position.
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_beta(const Vector< double > &x, double &beta)
Function to calculate the prescribed beta at a given spatial position.
SteadyAxisymAdvectionDiffusionPrescribedAlphaFctPt & alpha_fct_pt()
Access function for the prescribed-alpha function pointer.
unsigned U_index_adv_diff
The index at which the unknown is stored at the nodes.
void(* SteadyAxisymAdvectionDiffusionPrescribedAlphaFctPt)(const Vector< double > &x, double &alpha)
Function pointer to the prescribed-alpha function fct(x,alpha(x)) – x is a Vector!
SteadyAxisymAdvectionDiffusionPrescribedBetaFctPt & beta_fct_pt()
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 ...
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(std::ostream &outfile)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
SteadyAxisymAdvectionDiffusionFluxElement(const SteadyAxisymAdvectionDiffusionFluxElement &dummy)=delete
Broken copy constructor.
SteadyAxisymAdvectionDiffusionFluxElement()
Broken empty constructor.
SteadyAxisymAdvectionDiffusionPrescribedAlphaFctPt Alpha_fct_pt
Function pointer to the (global) prescribed-alpha function.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector.
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...
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...
SteadyAxisymAdvectionDiffusionPrescribedBetaFctPt Beta_fct_pt
Function pointer to the (global) prescribed-beta function.
void(* SteadyAxisymAdvectionDiffusionPrescribedBetaFctPt)(const Vector< double > &x, double &beta)
Function pointer to the prescribed-beta function fct(x,beta(x)) – x is a Vector!
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).