27 #ifndef OOMPH_TPOISSON_ELEMENTS_HEADER
28 #define OOMPH_TPOISSON_ELEMENTS_HEADER
33 #include <oomph-lib-config.h>
38 #include "../generic/nodes.h"
39 #include "../generic/oomph_utilities.h"
40 #include "../generic/Telements.h"
41 #include "../generic/error_estimator.h"
58 template<
unsigned DIM,
unsigned NNODE_1D>
91 void output(std::ostream& outfile,
const unsigned& n_plot)
107 void output(FILE* file_pt,
const unsigned& n_plot)
116 const unsigned& n_plot,
126 const unsigned& n_plot,
169 return (NNODE_1D - 1);
211 template<
unsigned DIM,
unsigned NNODE_1D>
219 unsigned n_node = this->nnode();
222 double J = this->dshape_eulerian(
s, psi, dpsidx);
226 for (
unsigned i = 0;
i < n_node;
i++)
229 dtestdx(
i, 0) = dpsidx(
i, 0);
230 dtestdx(
i, 1) = dpsidx(
i, 1);
244 template<
unsigned DIM,
unsigned NNODE_1D>
253 double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
272 template<
unsigned DIM,
unsigned NNODE_1D>
285 const double J = this->dshape_eulerian_at_knot(
286 ipt, psi, dpsidx, djacobian_dX, d_dpsidx_dX);
291 d_dtestdx_dX = d_dpsidx_dX;
304 template<
unsigned DIM,
unsigned NNODE_1D>
306 :
public virtual TElement<DIM - 1, NNODE_1D>
317 template<
unsigned NNODE_1D>
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional TElement.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional TElement.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as .
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
A class for all isoparametric elements that solve the Poisson equations.
void output(std::ostream &outfile)
Output with default number of plot points.
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: flux[i] = du/dx_i.
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: x,y,u_exact or x,y,z,u_exact at n_plot^DIM plot points.
////////////////////////////////////////////////////////////////// //////////////////////////////////...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
//////////////////////////////////////////////////////////////////////
void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Output function for a time-dependent exact solution. x,y,u_exact (calls the steady version)
void output(FILE *file_pt)
C-style output function: x,y,u or x,y,z,u.
void output(std::ostream &outfile)
Output function: x,y,u or x,y,z,u.
TPoissonElement()
Constructor: Call constructors for TElement and Poisson equations.
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output function for an exact solution: x,y,u_exact.
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
unsigned nvertex_node() const
Number of vertex nodes in the element.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: x,y,u or x,y,z,u at n_plot^DIM plot points.
double dshape_and_dtest_eulerian_at_knot_poisson(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Shape, test functions & derivs. w.r.t. to global coords. Return Jacobian.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function: x,y,u or x,y,z,u at n_plot^DIM plot points.
double dshape_and_dtest_eulerian_poisson(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
Access function for Nvalue: # of ‘values’ (pinned or dofs) at node n (always returns the same value a...
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
void operator=(const TPoissonElement< DIM, NNODE_1D > &)=delete
Broken assignment operator.
TPoissonElement(const TPoissonElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Standard flux.from Poisson equations.
static const unsigned Initial_Nvalue
Static unsigned that holds the (same) number of variables at every node.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...