32 #ifndef OOMPH_REFINEABLE_POISSON_ELEMENTS_HEADER 33 #define OOMPH_REFINEABLE_POISSON_ELEMENTS_HEADER 37 #include <oomph-lib-config.h> 42 #include "../generic/refineable_quad_element.h" 43 #include "../generic/refineable_brick_element.h" 44 #include "../generic/hp_refineable_elements.h" 45 #include "../generic/error_estimator.h" 61 template <
unsigned DIM>
94 std::ostream &outfile,
96 double& error,
double& norm);
108 unsigned n_node =
nnode();
123 for(
unsigned l=0;l<n_node;l++)
125 values[0] += this->
nodal_value(l,u_nodal_index)*psi[l];
140 "Time-dependent version of get_interpolated_values() ";
141 error_message +=
"not implemented for this element \n";
144 "RefineablePoissonEquations::get_interpolated_values()",
145 OOMPH_EXCEPTION_LOCATION);
173 const unsigned& flag);
180 dresidual_dnodal_coordinates);
190 template <
unsigned DIM,
unsigned NNODE_1D>
248 template<
unsigned DIM>
310 void compute_energy_error(std::ostream &outfile,
312 double& error,
double& norm);
329 template<
unsigned DIM,
unsigned NNODE_1D>
331 public virtual QElement<DIM-1,NNODE_1D>
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
void broken_copy(const std::string &class_name)
Issue error message and terminate execution.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation:
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
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. Overwrites default implementation in FiniteElement base class. dresidual_dnodal_coordinates(l,i,j) = d res(l) / dX_{ij}.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
virtual void set_integration_scheme(Integral *const &integral_pt)
Set the spatial integration scheme.
unsigned nvertex_node() const
Number of vertex nodes in the element.
p-refineable version of 2D QPoissonElement elements
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
RefineablePoissonEquations()
Constructor, simply call other constructors.
RefineableQPoissonElement()
Constructor, simply call the other constructors.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
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...
PRefineableQPoissonElement(const PRefineableQPoissonElement< DIM > &dummy)
Broken copy constructor.
void operator=(const RefineableQPoissonElement< DIM, NNODE_1D > &)
Broken assignment operator.
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as ...
virtual unsigned u_index_poisson() const
Return the index at which the unknown value is stored. The default value, 0, is appropriate for singl...
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
unsigned nvertex_node() const
Number of vertex nodes in the element.
virtual unsigned nvertex_node() const
RefineableQPoissonElement(const RefineableQPoissonElement< DIM, NNODE_1D > &dummy)
Broken copy constructor.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 1.
~PRefineableQPoissonElement()
Destructor (to avoid memory leaks)
RefineablePoissonEquations(const RefineablePoissonEquations< DIM > &dummy)
Broken copy constructor.
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
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...
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: flux[i] = du/dx_i.
void operator=(const PRefineableQPoissonElement< DIM > &)
Broken assignment operator.
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...
void broken_assign(const std::string &class_name)
Issue error message and terminate execution.
PoissonSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
void operator=(const RefineablePoissonEquations< DIM > &)
Broken assignment operator.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes...
PRefineableQPoissonElement()
Constructor, simply call the other constructors.
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
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...
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement. ...
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...
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 1.
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Standard flux.from Poisson equations.
unsigned nnode() const
Return the number of nodes.
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.
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...
void further_build()
Further build: Copy source function pointer from father element.
PoissonSourceFctPt Source_fct_pt
Pointer to source function: