27 #ifndef OOMPH_PML_FOURIER_DECOMPOSED_HELMHOLTZ_ELEMENTS_HEADER
28 #define OOMPH_PML_FOURIER_DECOMPOSED_HELMHOLTZ_ELEMENTS_HEADER
33 #include <oomph-lib-config.h>
41 #include "../generic/projection.h"
42 #include "../generic/nodes.h"
43 #include "../generic/Qelements.h"
44 #include "../generic/oomph_utilities.h"
45 #include "../generic/pml_meshes.h"
46 #include "../generic/projection.h"
47 #include "../generic/oomph_definitions.h"
55 namespace Legendre_functions_helper
58 extern double factorial(
const unsigned& l);
61 extern double plgndr1(
const unsigned& n,
const double& x);
64 extern double plgndr2(
const unsigned& l,
88 virtual std::complex<double>
gamma(
const double& nu_i,
89 const double& pml_width_i,
90 const double& k_squared) = 0;
97 const double& pml_width_i,
98 const double& pml_inner_boundary,
99 const double& k_squared) = 0;
117 std::complex<double>
gamma(
const double& nu_i,
118 const double& pml_width_i,
119 const double& k_squared)
124 std::complex<double>(1.0 / sqrt(k_squared), 0) *
125 std::complex<double>(0.0, 1.0 / (std::fabs(pml_width_i - nu_i)));
132 const double& pml_width_i,
133 const double& pml_inner_boundary,
134 const double& k_squared)
137 double log_arg = 1.0 - std::fabs(nu_i / pml_width_i);
138 return std::complex<double>(pml_inner_boundary + nu_i,
139 -log(log_arg) / sqrt(k_squared));
200 return std::complex<unsigned>(0, 1);
218 "Please set pointer to k_squared using access fct to pointer!",
219 OOMPH_CURRENT_FUNCTION,
220 OOMPH_EXCEPTION_LOCATION);
262 const unsigned n_plot = 5;
268 void output(std::ostream& outfile,
const unsigned& n_plot);
277 const unsigned& n_plot);
282 const unsigned n_plot = 5;
288 void output(FILE* file_pt,
const unsigned& n_plot);
293 const unsigned& n_plot,
299 std::ostream& outfile,
300 const unsigned& n_plot,
304 throw OomphLibError(
"There is no time-dependent output_fct() for "
305 "PMLFourierDecomposedHelmholtz elements ",
306 OOMPH_CURRENT_FUNCTION,
307 OOMPH_EXCEPTION_LOCATION);
318 const unsigned& n_plot,
336 throw OomphLibError(
"There is no time-dependent compute_error() for "
337 "PMLFourierDecomposedHelmholtz elements",
338 OOMPH_CURRENT_FUNCTION,
339 OOMPH_EXCEPTION_LOCATION);
367 std::complex<double>& source)
const
372 source = std::complex<double>(0.0, 0.0);
377 (*Source_fct_pt)(x, source);
387 values_to_pin.resize(2);
388 for (
unsigned j = 0; j < 2; j++)
390 values_to_pin[j] = j;
397 Vector<std::complex<double>>& flux)
const
400 const unsigned n_node =
nnode();
410 const std::complex<double> zero(0.0, 0.0);
411 for (
unsigned j = 0; j < 2; j++)
417 for (
unsigned l = 0; l < n_node; l++)
420 const std::complex<double> u_value(
427 for (
unsigned j = 0; j < 2; j++)
429 flux[j] += u_value * dpsidx(l, j);
452 residuals, jacobian, 1);
462 const unsigned n_node =
nnode();
471 std::complex<double> interpolated_u(0.0, 0.0);
474 const unsigned u_nodal_index_real =
476 const unsigned u_nodal_index_imag =
480 for (
unsigned l = 0; l < n_node; l++)
483 const std::complex<double> u_value(
487 interpolated_u += u_value * psi[l];
489 return interpolated_u;
505 Vector<std::complex<double>>& pml_laplace_factor,
506 std::complex<double>& pml_k_squared_factor)
510 for (
unsigned k = 0; k < 2; k++)
518 for (
unsigned k = 0; k < 2; k++)
532 for (
unsigned k = 0; k < 2; k++)
538 nu[k], pml_width[k], k_squared_local);
550 pml_laplace_factor[0] = pml_gamma[1] / pml_gamma[0];
551 pml_laplace_factor[1] = pml_gamma[0] / pml_gamma[1];
552 pml_k_squared_factor = pml_gamma[0] * pml_gamma[1];
558 for (
unsigned k = 0; k < 2; k++)
560 pml_laplace_factor[k] = std::complex<double>(1.0, 0.0);
563 pml_k_squared_factor = std::complex<double>(1.0, 0.0);
572 std::complex<double>& complex_r)
599 complex_r = std::complex<double>(r, 0.0);
630 DShape& dtestdx)
const = 0;
640 DShape& dtestdx)
const = 0;
647 const unsigned& flag);
681 template<
unsigned NNODE_1D>
683 :
public virtual QElement<2, NNODE_1D>,
723 void output(std::ostream& outfile,
const unsigned& n_plot)
735 const unsigned& n_plot)
748 void output(FILE* file_pt,
const unsigned& n_plot)
756 const unsigned& n_plot,
760 outfile, n_plot, exact_soln_pt);
770 const unsigned& n_plot,
774 outfile, phi, n_plot, exact_soln_pt);
782 const unsigned& n_plot,
787 outfile, n_plot, time, exact_soln_pt);
821 template<
unsigned NNODE_1D>
831 const double J = this->dshape_eulerian(
s, psi, dpsidx);
848 template<
unsigned NNODE_1D>
858 const double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
880 template<
unsigned NNODE_1D>
882 :
public virtual QElement<1, NNODE_1D>
899 template<
class FOURIER_DECOMPOSED_HELMHOLTZ_ELEMENT>
916 std::stringstream error_stream;
917 error_stream <<
"Fourier decomposed Helmholtz elements only store 2 "
919 << fld <<
" is illegal \n";
921 error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
926 unsigned nnod = this->
nnode();
930 for (
unsigned j = 0; j < nnod; j++)
933 data_values[j] = std::make_pair(this->
node_pt(j), fld);
953 std::stringstream error_stream;
954 error_stream <<
"Helmholtz elements only store two fields so fld = "
955 << fld <<
" is illegal\n";
957 error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
979 std::stringstream error_stream;
980 error_stream <<
"Helmholtz elements only store two fields so fld = "
981 << fld <<
" is illegal.\n";
983 error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
986 unsigned n_dim = this->
dim();
987 unsigned n_node = this->
nnode();
989 DShape dpsidx(n_node, n_dim), dtestdx(n_node, n_dim);
991 this->dshape_and_dtest_eulerian_pml_fourier_decomposed_helmholtz(
992 s, psi, dpsidx, test, dtestdx);
1000 const unsigned& fld,
1006 std::stringstream error_stream;
1007 error_stream <<
"Helmholtz elements only store two fields so fld = "
1008 << fld <<
" is illegal\n";
1010 error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1014 std::complex<unsigned> complex_u_nodal_index =
1015 this->u_index_pml_fourier_decomposed_helmholtz();
1016 unsigned u_nodal_index = 0;
1019 u_nodal_index = complex_u_nodal_index.real();
1023 u_nodal_index = complex_u_nodal_index.imag();
1028 unsigned n_node = this->
nnode();
1032 this->
shape(s, psi);
1035 double interpolated_u = 0.0;
1038 for (
unsigned l = 0; l < n_node; l++)
1040 interpolated_u += this->
nodal_value(t, l, u_nodal_index) * psi[l];
1042 return interpolated_u;
1052 std::stringstream error_stream;
1053 error_stream <<
"Helmholtz elements only store two fields so fld = "
1054 << fld <<
" is illegal\n";
1056 error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1059 return this->
nnode();
1069 std::stringstream error_stream;
1070 error_stream <<
"Helmholtz elements only store two fields so fld = "
1071 << fld <<
" is illegal\n";
1073 error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1076 std::complex<unsigned> complex_u_nodal_index =
1077 this->u_index_pml_fourier_decomposed_helmholtz();
1078 unsigned u_nodal_index = 0;
1081 u_nodal_index = complex_u_nodal_index.real();
1085 u_nodal_index = complex_u_nodal_index.imag();
1093 void output(std::ostream& outfile,
const unsigned& nplot)
1104 template<
class ELEMENT>
1117 template<
class ELEMENT>
1135 template<
unsigned NNODE_1D>
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
unsigned ntstorage() const
Return total number of doubles stored per value to record time history of each value (one for steady ...
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
A general Finite Element class.
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
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...
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...
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
unsigned nnode() const
Return the number of nodes.
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
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...
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as .
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
TimeStepper *& position_time_stepper_pt()
Return a pointer to the position timestepper.
An OomphLibError object which should be thrown when an run-time error is encountered....
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
bool Pml_is_enabled
Boolean indicating if element is used in pml mode.
std::vector< bool > Pml_direction_active
Coordinate direction along which pml boundary is constant; alternatively: coordinate direction in whi...
Vector< double > Pml_outer_boundary
Coordinate of outer pml boundary (Storage is provided for any coordinate direction; only the entries ...
Vector< double > Pml_inner_boundary
Coordinate of inner pml boundary (Storage is provided for any coordinate direction; only the entries ...
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
std::complex< double > interpolated_u_pml_fourier_decomposed_helmholtz(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
int pml_fourier_wavenumber()
Get the Fourier wavenumber.
virtual void fill_in_generic_residual_contribution_pml_fourier_decomposed_helmholtz(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Compute element residual Vector only (if flag=and/or element Jacobian matrix.
PMLMappingAndTransformedCoordinate *const & pml_mapping_and_transformed_coordinate_pt() const
Return a pointer to the PML Mapping object (const version)
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
double * K_squared_pt
Pointer to k^2 (wavenumber squared)
void output_real(std::ostream &outfile, const double &phi, const unsigned &n_plot)
Output function for real part of full time-dependent solution u = Re( (u_r +i u_i) exp(-i omega t) at...
double *& alpha_pt()
Get pointer to complex shift.
void(* PMLFourierDecomposedHelmholtzSourceFctPt)(const Vector< double > &x, std::complex< double > &f)
Function pointer to source function fct(x,f(x)) – x is a Vector!
void output(std::ostream &outfile)
Output with default number of plot points.
virtual double dshape_and_dtest_eulerian_pml_fourier_decomposed_helmholtz(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Shape/test functions and derivs w.r.t. to global coords at local coord. s; return Jacobian of mapping...
void compute_complex_r(const unsigned &ipt, const Vector< double > &x, std::complex< double > &complex_r)
Compute complex variable r at position x[0] and integration point ipt.
static double Default_Physical_Constant_Value
Static default value for the physical constants (initialised to zero)
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Dummy, time dependent error checker.
double * Alpha_pt
Pointer to wavenumber complex shift.
unsigned self_test()
Self-test: Return 0 for OK.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Add the element's contribution to its residual vector and element Jacobian matrix (wrapper)
virtual void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: (dummy time-dependent version to keep intel compiler happy)
PMLFourierDecomposedHelmholtzEquations(const PMLFourierDecomposedHelmholtzEquations &dummy)=delete
Broken copy constructor.
double k_squared()
Get k squared.
virtual void get_source_pml_fourier_decomposed_helmholtz(const unsigned &ipt, const Vector< double > &x, std::complex< double > &source) const
Get source term at (Eulerian) position x. This function is virtual to allow overloading in multi-phys...
int * N_pml_fourier_pt
Pointer to Fourier wave number.
PMLFourierDecomposedHelmholtzSourceFctPt Source_fct_pt
Pointer to source function:
PMLFourierDecomposedHelmholtzSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
PMLFourierDecomposedHelmholtzSourceFctPt source_fct_pt() const
Access function: Pointer to source function. Const version.
int *& pml_fourier_wavenumber_pt()
Get pointer to Fourier wavenumber.
void compute_norm(double &norm)
Compute norm of fe solution.
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 fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector (wrapper)
void values_to_be_pinned_on_outer_pml_boundary(Vector< unsigned > &values_to_pin)
Pure virtual function in which we specify the values to be pinned (and set to zero) on the outer edge...
void compute_pml_coefficients(const unsigned &ipt, const Vector< double > &x, Vector< std::complex< double >> &pml_laplace_factor, std::complex< double > &pml_k_squared_factor)
Compute pml coefficients at position x and integration point ipt. pml_laplace_factor is used in the r...
void output(FILE *file_pt)
C_style output with default number of plot points.
double *& k_squared_pt()
Get pointer to frequency.
void output_real_fct(std::ostream &outfile, const double &phi, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output function for real part of full time-dependent fct u = Re( (u_r +i u_i) exp(-i omega t) at phas...
double alpha()
Get complex shift.
virtual double dshape_and_dtest_eulerian_at_knot_pml_fourier_decomposed_helmholtz(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of ...
PMLMappingAndTransformedCoordinate *& pml_mapping_and_transformed_coordinate_pt()
Return a pointer to the PML Mapping object.
static BermudezPMLMappingAndTransformedCoordinate Default_pml_mapping_and_transformed_coordinate
Static so that the class doesn't need to instantiate a new default everytime it uses it.
PMLFourierDecomposedHelmholtzEquations()
Constructor.
PMLMappingAndTransformedCoordinate * Pml_mapping_and_transformed_coordinate_pt
Pointer to class which holds the pml mapping function (also known as gamma) and the associated transf...
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.
virtual std::complex< unsigned > u_index_pml_fourier_decomposed_helmholtz() const
Broken assignment operator.
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....
Wrapper class for projectable elements. Adds "projectability" to the underlying ELEMENT.
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
unsigned nhistory_values_for_projection(const unsigned &fld)
Number of history values to be stored for fld-th field. (Note: count includes current value!...
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Return interpolated field fld at local coordinate s, at time level t (t=0: present; t>0: history valu...
ProjectablePMLFourierDecomposedHelmholtzElement()
Constructor [this was only required explicitly from gcc 4.5.2 onwards...].
unsigned nfields_for_projection()
Number of fields to be projected: 2 (real and imag part)
unsigned nhistory_values_for_coordinate_projection()
Number of positional history values (Note: count includes current value!)
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Return Jacobian of mapping and shape functions of field fld at local coordinate s.
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Specify the values associated with field fld. The information is returned in a vector of pairs which ...
void output(std::ostream &outfile, const unsigned &nplot)
Output FE representation of soln: x,y,u or x,y,z,u at n_plot^DIM plot points.
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld: One per node.
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
//////////////////////////////////////////////////////////////////////// ////////////////////////////...
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(std::ostream &outfile, const unsigned &n_plot)
Output function: r,z,u at n_plot^2 plot points.
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. r,z,u_exact at n_plot^2 plot points (Calls the s...
static const unsigned Initial_Nvalue
Static int that holds the number of variables at nodes: always the same.
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 at n_plot^2 plot points.
void output_real_fct(std::ostream &outfile, const double &phi, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output function for real part of full time-dependent fct u = Re( (u_r +i u_i) exp(-i omega t) at phas...
QPMLFourierDecomposedHelmholtzElement(const QPMLFourierDecomposedHelmholtzElement< NNODE_1D > &dummy)=delete
Broken copy constructor.
void output(FILE *file_pt)
C-style output function: r,z,u.
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. at integration point ipt....
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function: r,z,u at n_plot^2 plot points.
QPMLFourierDecomposedHelmholtzElement()
Constructor: Call constructors for QElement and PMLFourierDecomposedHelmholtz equations.
void output_real(std::ostream &outfile, const double &phi, const unsigned &n_plot)
Output function for real part of full time-dependent solution u = Re( (u_r +i u_i) exp(-i omega t) at...
void output(std::ostream &outfile)
Output function: r,z,u.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
unsigned ntstorage() const
Return the number of doubles required to represent history (one for steady)
void output()
Doc the command line arguments.
double factorial(const unsigned &l)
Factorial.
double plgndr1(const unsigned &n, const double &x)
Legendre polynomials depending on one parameter.
double plgndr2(const unsigned &l, const unsigned &m, const double &x)
Legendre polynomials depending on two parameters.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...