31 #ifndef OOMPH_TPML_HELMHOLTZ_ELEMENTS_HEADER 32 #define OOMPH_TPML_HELMHOLTZ_ELEMENTS_HEADER 37 #include <oomph-lib-config.h> 42 #include "../generic/nodes.h" 43 #include "../generic/oomph_utilities.h" 44 #include "../generic/Telements.h" 45 #include "../generic/error_estimator.h" 65 template <
unsigned DIM,
unsigned NNODE_1D,
class PML_ELEMENT>
111 void output(std::ostream &outfile,
const unsigned &n_plot)
127 void output(FILE* file_pt,
const unsigned &n_plot)
135 void output_fct(std::ostream &outfile,
const unsigned &n_plot,
144 void output_fct(std::ostream &outfile,
const unsigned &n_plot,
187 for (
unsigned i=0;
i<DIM;
i++)
189 flux[count++]=complex_flux[
i].real();
190 flux[count++]=complex_flux[
i].imag();
212 template<
unsigned DIM,
unsigned NNODE_1D,
class PML_ELEMENT>
224 template<
unsigned DIM,
unsigned NNODE_1D,
class PML_ELEMENT>
234 unsigned n_node = this->
nnode();
241 for(
unsigned i=0;
i<n_node;
i++)
244 dtestdx(
i,0) = dpsidx(
i,0);
245 dtestdx(
i,1) = dpsidx(
i,1);
260 template<
unsigned DIM,
unsigned NNODE_1D,
class PML_ELEMENT>
290 template<
unsigned DIM,
unsigned NNODE_1D,
class PML_ELEMENT>
292 public virtual TElement<DIM-1,NNODE_1D>
307 template<
unsigned NNODE_1D,
class PML_ELEMENT>
329 template<
unsigned NNODE_1D,
class PML_ELEMENT>
349 template<
unsigned DIM,
unsigned NNODE_1D,
class PML_ELEMENT>
351 public virtual QElement<DIM-1,NNODE_1D>
370 template<
unsigned NNODE_1D,
class PML_ELEMENT>
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Standard flux from UnsteadyHeat equations.
void broken_copy(const std::string &class_name)
Issue error message and terminate execution.
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
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) ...
virtual double dshape_eulerian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidx) const
Return the geometric shape functions and also first derivatives w.r.t. global coordinates at the ipt-...
virtual void output(std::ostream &outfile)
Output the element data — typically the values at the nodes in a format suitable for post-processing...
PMLHelmholtz upgraded to become projectable.
double dshape_and_dtest_eulerian_helmholtz(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.
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.
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as ...
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...
TPMLHelmholtzElement()
Constructor: Call constructors for TElement and PMLHelmholtz equations.
void output(FILE *file_pt)
C-style output function: x,y,u or x,y,z,u.
unsigned nvertex_node() const
Number of vertex nodes in the element.
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as ...
unsigned required_nvalue(const unsigned &n) const
Broken assignment operator.
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Compute the geometric shape functions and also first derivatives w.r.t. global coordinates at local c...
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
EquivalentQElement()
Constructor: Call the constructor for the appropriate QElement.
virtual void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output an exact solution over the element.
static const unsigned Initial_Nvalue
Static unsigned that holds the (same) number of variables at every node.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional TElement. ...
EquivalentQElement()
Constructor: Call the constructor for the appropriate QElement.
void output(std::ostream &outfile)
Output function: x,y,u or x,y,z,u.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional TElement. ...
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional TElement. ...
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.
QPMLHelmholtzElement elements are linear/quadrilateral/ brick-shaped PMLHelmholtz elements with isopa...
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.
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 nnode() const
Return the number of nodes.
double dshape_and_dtest_eulerian_at_knot_helmholtz(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Shape, test functions & derivs. w.r.t. to global coords. Return Jacobian.