29 #ifndef OOMPH_PML_TIME_HARMONIC_LINEAR_ELASTICITY_ELEMENTS_HEADER
30 #define OOMPH_PML_TIME_HARMONIC_LINEAR_ELASTICITY_ELEMENTS_HEADER
34 #include <oomph-lib-config.h>
46 #include "../generic/Qelements.h"
47 #include "../generic/mesh.h"
48 #include "../generic/hermite_elements.h"
50 #include "../generic/projection.h"
51 #include "../generic/nodes.h"
52 #include "../generic/oomph_utilities.h"
53 #include "../generic/pml_meshes.h"
54 #include "../generic/pml_mapping_functions.h"
77 template<
unsigned DIM>
99 const unsigned i)
const
101 return std::complex<unsigned>(
i,
i + DIM);
109 unsigned n_node =
nnode();
117 for (
unsigned i = 0;
i < DIM;
i++)
120 std::complex<unsigned> u_nodal_index =
124 disp[
i] = std::complex<double>(0.0, 0.0);
127 for (
unsigned l = 0; l < n_node; l++)
129 const std::complex<double> u_value(
133 disp[
i] += u_value * psi[l];
143 unsigned n_node =
nnode();
152 std::complex<unsigned> u_nodal_index =
156 std::complex<double> interpolated_u(0.0, 0.0);
159 for (
unsigned l = 0; l < n_node; l++)
161 const std::complex<double> u_value(
165 interpolated_u += u_value * psi[l];
168 return (interpolated_u);
185 inline std::complex<double>
E(
const unsigned&
i,
188 const unsigned& l)
const
194 inline double nu()
const
228 DenseMatrix<std::complex<double>>& sigma)
const = 0;
237 Vector<std::complex<double>>& b)
const
244 for (
unsigned i = 0;
i < n;
i++)
246 b[
i] = std::complex<double>(0.0, 0.0);
252 (*Body_force_fct_pt)(x, b);
263 values_to_pin.resize(DIM * 2);
264 for (
unsigned j = 0; j < DIM * 2; j++)
266 values_to_pin[j] = j;
285 std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list)
const
289 std::pair<unsigned long, unsigned> dof_lookup;
292 const unsigned n_node = this->
nnode();
295 int local_unknown = 0;
298 for (
unsigned n = 0; n < n_node; n++)
301 for (
unsigned i = 0;
i < 2 * DIM;
i++)
307 if (local_unknown >= 0)
311 dof_lookup.first = this->
eqn_number(local_unknown);
312 dof_lookup.second = 0;
315 dof_lookup_list.push_front(dof_lookup);
334 Vector<std::complex<double>>& pml_inverse_jacobian_diagonal,
335 std::complex<double>& pml_jacobian_det)
339 for (
unsigned k = 0; k < DIM; k++)
341 pml_inverse_jacobian_diagonal[k] = std::complex<double>(1.0, 0.0);
343 pml_jacobian_det = std::complex<double>(1.0, 0.0);
350 for (
unsigned k = 0; k < DIM; k++)
358 for (
unsigned k = 0; k < DIM; k++)
368 std::ostringstream error_message_stream;
369 error_message_stream <<
"Pml_mapping_pt needs to be set "
373 OOMPH_CURRENT_FUNCTION,
374 OOMPH_EXCEPTION_LOCATION);
381 double wavenumber_squared = 2.0 * (1.0 + this->
nu()) * this->
omega_sq();
383 for (
unsigned k = 0; k < DIM; k++)
388 std::complex<double> pml_gamma =
393 pml_inverse_jacobian_diagonal[k] = 1.0 / pml_gamma;
395 pml_jacobian_det *= pml_gamma;
445 template<
unsigned DIM>
475 ->fill_in_generic_contribution_to_residuals_time_harmonic_linear_elasticity(
476 residuals, jacobian, 1);
486 const unsigned& nplot,
497 void output(std::ostream& outfile,
const unsigned& n_plot);
508 void output(FILE* file_pt,
const unsigned& n_plot);
516 std::ostream& outfile,
519 const unsigned& nplot);
529 const unsigned& n_plot);
537 const unsigned& n_plot);
557 std::ostringstream error_stream;
558 error_stream <<
"There is no time-dependent compute_error() " << std::endl
559 <<
"for pml time harmonic linear elasticity elements"
562 error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
581 template<
unsigned DIM,
unsigned NNODE_1D>
583 :
public virtual QElement<DIM, NNODE_1D>,
601 void output(std::ostream& outfile,
const unsigned& n_plot)
614 void output(FILE* file_pt,
const unsigned& n_plot)
711 template<
class TIME_HARMONIC_LINEAR_ELAST_ELEMENT>
733 unsigned nnod = this->
nnode();
734 for (
unsigned j = 0; j < nnod; j++)
737 data_values.push_back(std::make_pair(this->
node_pt(j), fld));
748 return 2 * this->
dim();
758 std::stringstream error_stream;
759 error_stream <<
"Elements only store four fields so fld can't be"
760 <<
" " << fld << std::endl;
762 error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
782 unsigned n_dim = this->
dim();
783 unsigned n_node = this->
nnode();
784 DShape dpsidx(n_node, n_dim);
798 unsigned n_node = this->
nnode();
811 double interpolated_u = 0.0;
814 for (
unsigned l = 0; l < n_node; l++)
818 if (nvalue != 2 * n_dim)
820 std::stringstream error_stream;
822 <<
"Current implementation only works for non-resized nodes\n"
823 <<
"but nvalue= " << nvalue <<
"!= 2 dim = " << 2 * n_dim
826 OOMPH_CURRENT_FUNCTION,
827 OOMPH_EXCEPTION_LOCATION);
830 interpolated_u += this->
nodal_value(t, l, fld) * psi[l];
832 return interpolated_u;
839 return this->
nnode();
849 if (nvalue != 2 * n_dim)
851 std::stringstream error_stream;
853 <<
"Current implementation only works for non-resized nodes\n"
854 <<
"but nvalue= " << nvalue <<
"!= 2 dim = " << 2 * n_dim
857 error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
869 template<
class ELEMENT>
882 template<
class ELEMENT>
901 template<
unsigned NNODE_1D>
A mapping function proposed by Bermudez et al, similar to the one above but is continuous across the ...
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
unsigned nvalue() const
Return number of values stored in data object (incl pinned ones).
unsigned ntstorage() const
Return total number of doubles stored per value to record time history of each value (one for steady ...
//////////////////////////////////////////////////////////////////////////// ////////////////////////...
FaceGeometry()
Constructor must call the constructor of the underlying solid element.
FaceGeometry()
Constructor must call the constructor of the underlying element.
FaceGeometry()
Constructor must call the constructor of the underlying element.
FaceGeometry()
Constructor must call the constructor of the underlying element.
FaceGeometry()
Constructor must call the constructor of the underlying element.
FaceGeometry()
Constructor must call the constructor of the underlying element.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
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 .
unsigned long eqn_number(const unsigned &ieqn_local) const
Return the global equation number corresponding to the ieqn_local-th local equation number.
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
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 ...
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....
Class to hold the mapping function (gamma) for the Pml which defines how the coordinates are transfor...
virtual std::complex< double > gamma(const double &nu_i, const double &pml_width_i, const double &wavenumber_squared, const double &alpha_shift=0.0)=0
Pure virtual to return Pml mapping gamma, where gamma is the as function of where where h is the v...
An isotropic elasticity tensor defined in terms of Young's modulus and Poisson's ratio....
double nu() const
Poisson's ratio.
A base class for elements that solve the equations of time-harmonic linear elasticity in Cartesian co...
double *& omega_sq_pt()
Access function for square of non-dim frequency.
BodyForceFctPt body_force_fct_pt() const
Access function: Pointer to body force function (const version)
unsigned ndof_types() const
The number of "DOF types" that degrees of freedom in this element are sub-divided into: for now lump ...
void body_force(const Vector< double > &x, Vector< std::complex< double >> &b) const
Evaluate body force at Eulerian coordinate x (returns zero vector if no body force function pointer h...
BodyForceFctPt & body_force_fct_pt()
Access function: Pointer to body force function.
std::complex< double > E(const unsigned &i, const unsigned &j, const unsigned &k, const unsigned &l) const
Access function to the entries in the elasticity tensor.
PMLTimeHarmonicIsotropicElasticityTensor *& elasticity_tensor_pt()
Return the pointer to the elasticity_tensor.
virtual void get_stress(const Vector< double > &s, DenseMatrix< std::complex< double >> &sigma) const =0
Return the Cauchy stress tensor, as calculated from the elasticity tensor at specified local coordina...
double nu() const
Access function to nu in the elasticity tensor.
static ContinuousBermudezPMLMapping Default_pml_mapping
Static so that the class doesn't need to instantiate a new default everytime it uses it.
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...
static double Default_omega_sq_value
Static default value for square of frequency.
PMLMapping * Pml_mapping_pt
Pointer to class which holds the pml mapping function, also known as gamma.
BodyForceFctPt Body_force_fct_pt
Pointer to body force function.
virtual std::complex< unsigned > u_index_time_harmonic_linear_elasticity(const unsigned i) const
Return the index at which the i-th real or imag unknown displacement component is stored....
const double & omega_sq() const
Access function for square of non-dim frequency.
PMLTimeHarmonicIsotropicElasticityTensor * Elasticity_tensor_pt
Pointer to the elasticity tensor.
std::complex< double > interpolated_u_time_harmonic_linear_elasticity(const Vector< double > &s, const unsigned &i) const
Return FE interpolated displacement u[i] at local coordinate s.
PMLTimeHarmonicLinearElasticityEquationsBase()
Constructor: Set null pointers for constitutive law and for isotropic growth function....
PMLMapping *const & pml_mapping_pt() const
Return a pointer to the PML Mapping object (const version)
void get_strain(const Vector< double > &s, DenseMatrix< std::complex< double >> &strain) const
Return the strain tensor.
PMLMapping *& pml_mapping_pt()
Return a pointer to the PML Mapping object.
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contain...
void interpolated_u_time_harmonic_linear_elasticity(const Vector< double > &s, Vector< std::complex< double >> &disp) const
Compute vector of FE interpolated displacement u at local coordinate s.
void(* BodyForceFctPt)(const Vector< double > &x, Vector< std::complex< double >> &b)
Function pointer to function that specifies the body force as a function of the Cartesian coordinates...
void compute_pml_coefficients(const unsigned &ipt, const Vector< double > &x, Vector< std::complex< double >> &pml_inverse_jacobian_diagonal, std::complex< double > &pml_jacobian_det)
Compute pml coefficients at position x and integration point ipt. pml_inverse_jacobian_diagonal conta...
double * Omega_sq_pt
Square of nondim frequency.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
The jacobian is calculated by finite differences by default, We need only to take finite differences ...
void output_imag(std::ostream &outfile, const double &phi, const unsigned &n_plot)
Output function for imaginary part of full time-dependent solution u = Im( (u_r +i u_i) exp(-i omega ...
void get_stress(const Vector< double > &s, DenseMatrix< std::complex< double >> &sigma) const
Return the Cauchy stress tensor, as calculated from the elasticity tensor at specified local coordina...
void compute_norm(double &norm)
Compute norm of solution: square of the L2 norm.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Return the residuals for the solid equations (the discretised principle of virtual displacements)
unsigned required_nvalue(const unsigned &n) const
Number of values required at node n.
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Dummy, time dependent error checker.
void output(std::ostream &outfile)
Output: x,y,[z],u_r,v_r,[w_r],u_i,v_i,[w_i].
virtual void fill_in_generic_contribution_to_residuals_time_harmonic_linear_elasticity(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Private helper function to compute residuals and (if requested via flag) also the Jacobian matrix.
void output_total_real(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt incoming_wave_fct_pt, const double &phi, const unsigned &nplot)
Output function for real part of full time-dependent solution constructed by adding the scattered fie...
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact solution x,y,[z],u_r,v_r,[w_r],u_i,v_i,[w_i].
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)) a...
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
void output(FILE *file_pt)
C-style output: x,y,[z],u_r,v_r,[w_r],u_i,v_i,[w_i].
PMLTimeHarmonicLinearElasticityEquations()
Constructor.
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. (includes present value!)
unsigned nvalue_of_field(const unsigned &fld)
Return number of values 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...
unsigned nfields_for_projection()
Number of fields to be projected: 2*dim, corresponding to real and imag parts of the displacement com...
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 ...
ProjectablePMLTimeHarmonicLinearElasticityElement()
Constructor [this was only required explicitly from gcc 4.5.2 onwards...].
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
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.
unsigned nhistory_values_for_coordinate_projection()
Number of positional history values: Read out from positional timestepper (Note: count includes curre...
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
void output(std::ostream &outfile)
Output function.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function.
void output(FILE *file_pt)
C-style output function.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function.
QPMLTimeHarmonicLinearElasticityElement()
Constructor.
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)
//////////////////////////////////////////////////////////////////// ////////////////////////////////...