28 #ifndef OOMPH_REFINEABLE_PML_HELMHOLTZ_ELEMENTS_HEADER
29 #define OOMPH_REFINEABLE_PML_HELMHOLTZ_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>
95 for (
unsigned i = 0;
i <
DIM;
i++)
97 flux[count] = complex_flux[
i].real();
99 flux[count] = complex_flux[
i].imag();
116 unsigned n_node =
nnode();
132 for (
unsigned l = 0; l < n_node; l++)
134 values[0] += this->
nodal_value(l, u_nodal_index.real()) * psi[l];
135 values[1] += this->
nodal_value(l, u_nodal_index.imag()) * psi[l];
151 "Time-dependent version of get_interpolated_values() ";
152 error_message +=
"not implemented for this element \n";
155 "RefineablePMLHelmholtzEquations::get_interpolated_values()",
185 const unsigned& flag);
194 template<
unsigned DIM,
unsigned NNODE_1D>
244 return (NNODE_1D - 1);
265 template<
unsigned DIM,
unsigned NNODE_1D>
267 :
public virtual QElement<DIM - 1, NNODE_1D>
279 template<
unsigned NNODE_1D>
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Definition: error_estimator.h:79
FaceGeometry()
Definition: refineable_pml_helmholtz_elements.h:272
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
virtual Node * vertex_node_pt(const unsigned &j) const
Definition: elements.h:2500
Definition: oomph_definitions.h:222
Definition: pml_helmholtz_elements.h:62
PMLHelmholtzSourceFctPt Source_fct_pt
Pointer to source function:
Definition: pml_helmholtz_elements.h:727
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.
Definition: pml_helmholtz_elements.h:446
PMLHelmholtzSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
Definition: pml_helmholtz_elements.h:400
virtual std::complex< unsigned > u_index_helmholtz() const
Broken assignment operator.
Definition: pml_helmholtz_elements.h:91
PMLLayerElement()
Definition: refineable_pml_helmholtz_elements.h:286
Definition: pml_meshes.h:48
Definition: Qelements.h:459
Definition: pml_helmholtz_elements.h:754
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_pml_helmholtz_elements.h:60
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Definition: refineable_pml_helmholtz_elements.h:144
void further_build()
Further build: Copy source function pointer from father element.
Definition: refineable_pml_helmholtz_elements.h:169
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Definition: refineable_pml_helmholtz_elements.h:109
unsigned num_Z2_flux_terms()
Broken assignment operator.
Definition: refineable_pml_helmholtz_elements.h:83
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Definition: refineable_pml_helmholtz_elements.h:90
void fill_in_generic_residual_contribution_helmholtz(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Definition: refineable_pml_helmholtz_elements.cc:41
RefineablePMLHelmholtzEquations()
Constructor, simply call other constructors.
Definition: refineable_pml_helmholtz_elements.h:63
RefineablePMLHelmholtzEquations(const RefineablePMLHelmholtzEquations< DIM > &dummy)=delete
Broken copy constructor.
Definition: Qelements.h:2259
Definition: refineable_pml_helmholtz_elements.h:199
unsigned nrecovery_order()
Definition: refineable_pml_helmholtz_elements.h:242
unsigned ncont_interpolated_values() const
Broken assignment operator.
Definition: refineable_pml_helmholtz_elements.h:220
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: refineable_pml_helmholtz_elements.h:232
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: refineable_pml_helmholtz_elements.h:226
RefineableQPMLHelmholtzElement()
Constructor, simply call the other constructors.
Definition: refineable_pml_helmholtz_elements.h:202
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
Definition: refineable_pml_helmholtz_elements.h:238
RefineableQPMLHelmholtzElement(const RefineableQPMLHelmholtzElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void further_setup_hanging_nodes()
Definition: refineable_pml_helmholtz_elements.h:249
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
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