32 #ifndef OOMPH_REFINEABLE_HELMHOLTZ_ELEMENTS_HEADER 33 #define OOMPH_REFINEABLE_HELMHOLTZ_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/error_estimator.h" 62 template <
unsigned DIM>
98 for (
unsigned i=0;
i<DIM;
i++)
100 flux[count] =actual_flux[
i].real();
101 flux[count+1]=actual_flux[
i].imag();
114 values.resize(2,0.0);
117 unsigned n_node =
nnode();
126 for(
unsigned l=0;l<n_node;l++)
144 "Time-dependent version of get_interpolated_values() ";
145 error_message +=
"not implemented for this element \n";
148 "RefineableHelmholtzEquations::get_interpolated_values()",
149 OOMPH_EXCEPTION_LOCATION);
182 const unsigned& flag);
193 template <
unsigned DIM,
unsigned NNODE_1D>
259 template<
unsigned DIM,
unsigned NNODE_1D>
261 public virtual QElement<DIM-1,NNODE_1D>
HelmholtzSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
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: Same order as shape functions.
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
virtual RefineableElement * father_element_pt() const
Return a pointer to the father element.
RefineableQHelmholtzElement(const RefineableQHelmholtzElement< DIM, NNODE_1D > &dummy)
Broken copy constructor.
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Standard flux.from Helmholtz equations.
HelmholtzSourceFctPt Source_fct_pt
Pointer to source function:
unsigned nvertex_node() const
Number of vertex nodes in the element.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes...
void get_flux(const Vector< double > &s, Vector< std::complex< double > > &flux) const
Get flux: flux[i] = du/dx_i for real and imag part.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
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...
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
RefineableHelmholtzEquations()
Constructor, simply call other constructors.
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
virtual unsigned nvertex_node() const
double *& k_squared_pt()
Get pointer to square of wavenumber.
double * K_squared_pt
Pointer to square of wavenumber.
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...
RefineableHelmholtzEquations(const RefineableHelmholtzEquations< DIM > &dummy)
Broken copy constructor.
virtual std::complex< unsigned > u_index_helmholtz() const
Broken assignment operator.
unsigned num_Z2_flux_terms()
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 further_build()
Further build: Copy source function pointer from father element.
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...
void fill_in_generic_residual_contribution_helmholtz(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 nnode() const
Return the number of nodes.
RefineableQHelmholtzElement()
Constructor, simply call the other constructors.
unsigned ncont_interpolated_values() const
Broken assignment operator.
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...