27 #ifndef OOMPH_TPML_FOURIER_DECOMPOSED_HELMHOLTZ_ELEMENTS_HEADER
28 #define OOMPH_TPML_FOURIER_DECOMPOSED_HELMHOLTZ_ELEMENTS_HEADER
33 #include <oomph-lib-config.h>
38 #include "../generic/nodes.h"
39 #include "../generic/oomph_utilities.h"
40 #include "../generic/Telements.h"
41 #include "../generic/error_estimator.h"
60 template<
unsigned NNODE_1D>
104 void output(std::ostream& outfile,
const unsigned& n_plot)
120 void output(FILE* file_pt,
const unsigned& n_plot)
129 const unsigned& n_plot,
133 outfile, n_plot, exact_soln_pt);
140 const unsigned& n_plot,
145 outfile, n_plot, time, exact_soln_pt);
173 return (NNODE_1D - 1);
189 for (
unsigned i = 0;
i < 2;
i++)
191 flux[count++] = complex_flux[
i].real();
192 flux[count++] = complex_flux[
i].imag();
223 template<
unsigned NNODE_1D>
232 unsigned n_node = this->nnode();
235 double J = this->dshape_eulerian(
s, psi, dpsidx);
239 for (
unsigned i = 0;
i < n_node;
i++)
242 dtestdx(
i, 0) = dpsidx(
i, 0);
243 dtestdx(
i, 1) = dpsidx(
i, 1);
257 template<
unsigned NNODE_1D>
267 double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
285 template<
unsigned NNODE_1D>
287 :
public virtual TElement<1, NNODE_1D>
304 template<
unsigned NNODE_1D>
323 template<
unsigned NNODE_1D>
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
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 TElement.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as .
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
void output(std::ostream &outfile)
Output with default number of plot points.
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: r,z,u_re_exact,u_im_exact at n_plot^2 plot points.
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.
PMLLayerElement()
Constructor: Call the constructor for the appropriate QElement.
PMLLayerElement()
Constructor: Call the constructor for the appropriate QElement.
General definition of policy class defining the elements to be used in the actual PML layers....
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
//////////////////////////////////////////////////////////////////////// ////////////////////////////...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
//////////////////////////////////////////////////////////////////////
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: r,z,u n_plot^2 plot points.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
TPMLFourierDecomposedHelmholtzElement(const TPMLFourierDecomposedHelmholtzElement< NNODE_1D > &dummy)=delete
Broken copy constructor.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
unsigned nvertex_node() const
Number of vertex nodes in the element.
void output(std::ostream &outfile)
Output function: r,z,u.
static const unsigned Initial_Nvalue
Static unsigned that holds the (same) number of variables at every node.
void output(FILE *file_pt)
C-style output function: r,z,u or x,y,z,u.
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output function for an exact solution: r,z,u_exact.
TPMLFourierDecomposedHelmholtzElement()
Constructor: Call constructors for TElement and PMLFourierDecomposedHelmholtz equations.
unsigned required_nvalue(const unsigned &n) const
Broken assignment operator.
double dshape_and_dtest_eulerian_at_knot_pml_fourier_decomposed_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.
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
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)
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Standard flux from UnsteadyHeat equations.
double dshape_and_dtest_eulerian_pml_fourier_decomposed_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(FILE *file_pt, const unsigned &n_plot)
C-style output function: r,z,u at n_plot^2 plot points.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...