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()",
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>
Definition: error_estimator.h:79
FaceGeometry()
Definition: refineable_poisson_elements.h:343
Definition: elements.h:4998
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
virtual unsigned nvertex_node() const
Definition: elements.h:2491
virtual void shape(const Vector< double > &s, Shape &psi) const =0
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Definition: elements.h:1759
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1963
virtual Node * vertex_node_pt(const unsigned &j) const
Definition: elements.h:2500
virtual void set_integration_scheme(Integral *const &integral_pt)
Set the spatial integration scheme.
Definition: elements.cc:3210
Definition: integral.h:1281
Definition: oomph_definitions.h:222
Definition: Qelements.h:2274
p-refineable version of 2D QPoissonElement elements
Definition: refineable_poisson_elements.h:254
void operator=(const PRefineableQPoissonElement< DIM > &)=delete
Broken assignment operator.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 1.
Definition: refineable_poisson_elements.h:286
unsigned nrecovery_order()
Definition: refineable_poisson_elements.h:312
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: refineable_poisson_elements.h:298
void further_build()
Further build: Copy source function pointer from father element.
Definition: refineable_poisson_elements.cc:615
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: refineable_poisson_elements.h:292
PRefineableQPoissonElement(const PRefineableQPoissonElement< DIM > &dummy)=delete
Broken copy constructor.
~PRefineableQPoissonElement()
Destructor (to avoid memory leaks)
Definition: refineable_poisson_elements.h:270
void compute_energy_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_grad_pt, double &error, double &norm)
Get error against and norm of exact solution.
Definition: refineable_poisson_elements.cc:529
PRefineableQPoissonElement()
Constructor, simply call the other constructors.
Definition: refineable_poisson_elements.h:257
Definition: poisson_elements.h:56
PoissonSourceFctPt Source_fct_pt
Pointer to source function:
Definition: poisson_elements.h:523
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: flux[i] = du/dx_i.
Definition: poisson_elements.h:364
virtual unsigned u_index_poisson() const
Definition: poisson_elements.h:85
PoissonSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
Definition: poisson_elements.h:281
Definition: Qelements.h:459
Definition: poisson_elements.h:542
A Rank 3 Tensor class.
Definition: matrices.h:1370
Definition: refineable_elements.h:97
virtual RefineableElement * father_element_pt() const
Return a pointer to the father element.
Definition: refineable_elements.h:539
Definition: refineable_poisson_elements.h:59
virtual void get_dresidual_dnodal_coordinates(RankThreeTensor< double > &dresidual_dnodal_coordinates)
Definition: refineable_poisson_elements.cc:345
void fill_in_generic_residual_contribution_poisson(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Definition: refineable_poisson_elements.cc:140
void operator=(const RefineablePoissonEquations< DIM > &)=delete
Broken assignment operator.
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Definition: refineable_poisson_elements.h:99
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Standard flux.from Poisson equations.
Definition: refineable_poisson_elements.h:83
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.
Definition: refineable_poisson_elements.cc:37
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
Definition: refineable_poisson_elements.h:77
RefineablePoissonEquations()
Constructor, simply call other constructors.
Definition: refineable_poisson_elements.h:62
RefineablePoissonEquations(const RefineablePoissonEquations< DIM > &dummy)=delete
Broken copy constructor.
void further_build()
Further build: Copy source function pointer from father element.
Definition: refineable_poisson_elements.h:156
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Definition: refineable_poisson_elements.h:132
Definition: Qelements.h:2259
Definition: refineable_poisson_elements.h:193
unsigned nrecovery_order()
Definition: refineable_poisson_elements.h:235
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: refineable_poisson_elements.h:225
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 1.
Definition: refineable_poisson_elements.h:213
RefineableQPoissonElement()
Constructor, simply call the other constructors.
Definition: refineable_poisson_elements.h:196
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: refineable_poisson_elements.h:219
void operator=(const RefineableQPoissonElement< DIM, NNODE_1D > &)=delete
Broken assignment operator.
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
Definition: refineable_poisson_elements.h:231
RefineableQPoissonElement(const RefineableQPoissonElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void further_setup_hanging_nodes()
Definition: refineable_poisson_elements.h:242
RealScalar s
Definition: level1_cplx_impl.h:130
#define DIM
Definition: linearised_navier_stokes_elements.h:44
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
int error
Definition: calibrate.py:297
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
t
Definition: plotPSD.py:36
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2