28 #ifndef OOMPH_REFINEABLE_POISSON_ELEMENTS_HEADER
29 #define OOMPH_REFINEABLE_POISSON_ELEMENTS_HEADER
33 #include <oomph-lib-config.h>
38 #include "../generic/refineable_quad_element.h"
39 #include "../generic/refineable_brick_element.h"
40 #include "../generic/hp_refineable_elements.h"
41 #include "../generic/error_estimator.h"
55 template<
unsigned DIM>
90 std::ostream& outfile,
106 unsigned n_node =
nnode();
121 for (
unsigned l = 0; l < n_node; l++)
123 values[0] += this->
nodal_value(l, u_nodal_index) * psi[l];
139 "Time-dependent version of get_interpolated_values() ";
140 error_message +=
"not implemented for this element \n";
143 "RefineablePoissonEquations::get_interpolated_values()",
144 OOMPH_EXCEPTION_LOCATION);
172 const unsigned& flag);
188 template<
unsigned DIM,
unsigned NNODE_1D>
237 return (NNODE_1D - 1);
249 template<
unsigned DIM>
318 std::ostream& outfile,
336 template<
unsigned DIM,
unsigned NNODE_1D>
338 :
public virtual QElement<DIM - 1, NNODE_1D>
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 QElement.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
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 unsigned nvertex_node() const
Return the number of vertex nodes in this element. Broken virtual function in "pure" finite elements.
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.
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
virtual Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element. Broken virtual function in "pure" finite elements.
virtual void set_integration_scheme(Integral *const &integral_pt)
Set the spatial integration scheme.
Class for multidimensional Gauss Lobatto Legendre integration rules empty - just establishes template...
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
An OomphLibError object which should be thrown when an run-time error is encountered....
A class that is used to template the p-refineable Q elements by dimension. It's really nothing more t...
p-refineable version of 2D QPoissonElement elements
void operator=(const PRefineableQPoissonElement< DIM > &)=delete
Broken assignment operator.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 1.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation:
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
void further_build()
Further build: Copy source function pointer from father element.
unsigned nvertex_node() const
Number of vertex nodes in the element.
PRefineableQPoissonElement(const PRefineableQPoissonElement< DIM > &dummy)=delete
Broken copy constructor.
~PRefineableQPoissonElement()
Destructor (to avoid memory leaks)
void compute_energy_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_grad_pt, double &error, double &norm)
Get error against and norm of exact solution.
PRefineableQPoissonElement()
Constructor, simply call the other constructors.
A class for all isoparametric elements that solve the Poisson equations.
PoissonSourceFctPt Source_fct_pt
Pointer to source function:
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: flux[i] = du/dx_i.
virtual unsigned u_index_poisson() const
Return the index at which the unknown value is stored. The default value, 0, is appropriate for singl...
PoissonSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
//////////////////////////////////////////////////////////////////////// ////////////////////////////...
////////////////////////////////////////////////////////////////// //////////////////////////////////...
RefineableElements are FiniteElements that may be subdivided into children to provide a better local ...
virtual RefineableElement * father_element_pt() const
Return a pointer to the father element.
////////////////////////////////////////////////////////////////////////
virtual void get_dresidual_dnodal_coordinates(RankThreeTensor< double > &dresidual_dnodal_coordinates)
Compute derivatives of elemental residual vector with respect to nodal coordinates....
void fill_in_generic_residual_contribution_poisson(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Add element's contribution to elemental residual vector and/or Jacobian matrix flag=1: compute both f...
void operator=(const RefineablePoissonEquations< DIM > &)=delete
Broken assignment operator.
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Get the function value u in Vector. Note: Given the generality of the interface (this function is usu...
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Standard flux.from Poisson equations.
void compute_exact_Z2_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_flux_pt, double &error, double &norm)
Get error against and norm of exact flux.
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
RefineablePoissonEquations()
Constructor, simply call other constructors.
RefineablePoissonEquations(const RefineablePoissonEquations< DIM > &dummy)=delete
Broken copy constructor.
void further_build()
Further build: Copy source function pointer from father element.
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Get the function value u in Vector. Note: Given the generality of the interface (this function is usu...
A class that is used to template the refineable Q elements by dimension. It's really nothing more tha...
Refineable version of 2D QPoissonElement elements.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 1.
RefineableQPoissonElement()
Constructor, simply call the other constructors.
unsigned nvertex_node() const
Number of vertex nodes in the element.
void operator=(const RefineableQPoissonElement< DIM, NNODE_1D > &)=delete
Broken assignment operator.
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
RefineableQPoissonElement(const RefineableQPoissonElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes....
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.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...