29#ifndef OOMPH_TIME_HARMONIC_LINEAR_ELASTICITY_ELEMENTS_HEADER
30#define OOMPH_TIME_HARMONIC_LINEAR_ELASTICITY_ELEMENTS_HEADER
34#include <oomph-lib-config.h>
62 template<
unsigned DIM>
70 const unsigned i)
const
72 return std::complex<unsigned>(
i,
i +
DIM);
88 for (
unsigned i = 0;
i <
DIM;
i++)
95 disp[
i] = std::complex<double>(0.0, 0.0);
100 const std::complex<double>
u_value(
127 std::complex<double> interpolated_u(0.0, 0.0);
132 const std::complex<double>
u_value(
139 return (interpolated_u);
167 inline double E(
const unsigned&
i,
170 const unsigned&
l)
const
204 DenseMatrix<std::complex<double>>& sigma)
const = 0;
213 Vector<std::complex<double>>& b)
const
220 for (
unsigned i = 0;
i <
n;
i++)
222 b[
i] = std::complex<double>(0.0, 0.0);
233 (*Body_force_fct_pt)(time, x, b);
256 std::pair<unsigned long, unsigned>
dof_lookup;
268 for (
unsigned i = 0;
i < 2 *
DIM;
i++)
313 template<
unsigned DIM>
344 ->fill_in_generic_contribution_to_residuals_time_harmonic_linear_elasticity(
355 const unsigned&
nplot,
400 template<
unsigned DIM,
unsigned NNODE_1D>
402 :
public virtual QElement<DIM, NNODE_1D>,
525 template<
class TIME_HARMONIC_LINEAR_ELAST_ELEMENT>
548 for (
unsigned j = 0;
j <
nnod;
j++)
562 return 2 * this->
dim();
573 error_stream <<
"Elements only store four fields so fld can't be"
574 <<
" " <<
fld << std::endl;
625 double interpolated_u = 0.0;
632 if (nvalue != 2 *
n_dim)
636 <<
"Current implementation only works for non-resized nodes\n"
637 <<
"but nvalue= " << nvalue <<
"!= 2 dim = " << 2 *
n_dim
646 return interpolated_u;
653 return this->
nnode();
663 if (nvalue != 2 *
n_dim)
667 <<
"Current implementation only works for non-resized nodes\n"
668 <<
"but nvalue= " << nvalue <<
"!= 2 dim = " << 2 *
n_dim
683 template<
class ELEMENT>
696 template<
class ELEMENT>
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
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 ...
Class for dense matrices, storing all the values of the matrix as a pointer to a pointer with assorte...
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.
FaceGeometry class definition: This policy class is used to allow construction of face elements that ...
A general Finite Element class.
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...
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
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...
TimeStepper *& position_time_stepper_pt()
Return a pointer to the position timestepper.
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
An OomphLibError object which should be thrown when an run-time error is encountered....
Wrapper class for projectable elements. Adds "projectability" to the underlying ELEMENT.
Time-harmonic linear elasticity upgraded to become projectable.
unsigned nfields_for_projection()
Number of fields to be projected: 2*dim, corresponding to real and imag parts of the displacement com...
unsigned nhistory_values_for_projection(const unsigned &fld)
Number of history values to be stored for fld-th field. (includes present value!)
unsigned nhistory_values_for_coordinate_projection()
Number of positional history values: Read out from positional timestepper (Note: count includes curre...
ProjectableTimeHarmonicLinearElasticityElement()
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 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 nvalue_of_field(const unsigned &fld)
Return number of values in field fld.
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 ...
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.
An Element that solves the equations of linear elasticity in Cartesian coordinates,...
void output(FILE *file_pt)
C-style output function.
void output(std::ostream &outfile)
Output function.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function.
QTimeHarmonicLinearElasticityElement()
Constructor.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
A base class that represents the fourth-rank elasticity tensor defined such that.
A base class for elements that solve the equations of time-harmonic linear elasticity in Cartesian co...
BodyForceFctPt & body_force_fct_pt()
Access function: Pointer to body force function.
double * Omega_sq_pt
Square of nondim frequency.
TimeHarmonicElasticityTensor *& elasticity_tensor_pt()
Return the pointer to the elasticity_tensor.
void body_force(const Vector< double > &x, Vector< std::complex< double > > &b) const
Evaluate body force at Eulerian coordinate x at present time (returns zero vector if no body force fu...
const double & omega_sq() const
Access function for square of non-dim frequency.
double *& omega_sq_pt()
Access function for square of non-dim frequency.
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...
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.
TimeHarmonicElasticityTensor * Elasticity_tensor_pt
Pointer to the elasticity tensor.
static double Default_omega_sq_value
Static default value for square of frequency.
unsigned ndof_types() const
The number of "DOF types" that degrees of freedom in this element are sub-divided into: for now lump ...
void(* BodyForceFctPt)(const double &t, 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...
TimeHarmonicLinearElasticityEquationsBase()
Constructor: Set null pointers for constitutive law and for isotropic growth function....
BodyForceFctPt Body_force_fct_pt
Pointer to body force function.
double E(const unsigned &i, const unsigned &j, const unsigned &k, const unsigned &l) const
Access function to the entries in the elasticity tensor.
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...
BodyForceFctPt body_force_fct_pt() const
Access function: Pointer to body force function (const version)
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....
void get_strain(const Vector< double > &s, DenseMatrix< std::complex< double > > &strain) const
Return the strain 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.
A class for elements that solve the equations of linear elasticity in cartesian coordinates.
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 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_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
The jacobian is calculated by finite differences by default, We need only to take finite differences ...
unsigned required_nvalue(const unsigned &n) const
Number of values required at node n.
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(FILE *file_pt)
C-style output: x,y,[z],u_r,v_r,[w_r],u_i,v_i,[w_i].
void output(std::ostream &outfile)
Output: x,y,[z],u_r,v_r,[w_r],u_i,v_i,[w_i].
TimeHarmonicLinearElasticityEquations()
Constructor.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Return the residuals for the solid equations (the discretised principle of virtual displacements)
unsigned ntstorage() const
Return the number of doubles required to represent history (one for steady)
Time *const & time_pt() const
Access function for the pointer to time (const version)
double & time()
Return the current value of the continuous time.
A slight extension to the standard template vector class so that we can include "graceful" array rang...
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).