28#ifndef OOMPH_AXISYMMETRIC_LINEAR_ELASTICITY_ELEMENTS_HEADER
29#define OOMPH_AXISYMMETRIC_LINEAR_ELASTICITY_ELEMENTS_HEADER
33#include <oomph-lib-config.h>
61 const unsigned&
i)
const
68 const unsigned&
i)
const
101 const unsigned&
i)
const
143 for (
unsigned i = 0;
i < 3;
i++)
179 double interpolated_u = 0.0;
189 return (interpolated_u);
207 for (
unsigned i = 0;
i < 3;
i++)
234 for (
unsigned i = 0;
i < 3;
i++)
283 std::ostringstream error_message;
284 error_message <<
"No pointer to Poisson's ratio set. Please set one!\n";
334 for (
unsigned i = 0;
i <
n;
i++)
341 (*Body_force_fct_pt)(time, x, b);
364 std::pair<unsigned long, unsigned>
dof_lookup;
376 for (
unsigned i = 0;
i < 3;
i++)
462 ->fill_in_generic_contribution_to_residuals_axisymmetric_linear_elasticity(
472 const unsigned&
nplot,
478 const unsigned&
nplot,
537 template<
unsigned NNODE_1D>
539 :
public virtual QElement<2, NNODE_1D>,
580 template<
unsigned NNODE_1D>
582 :
public virtual QElement<1, NNODE_1D>
658 template<
class AXISYM_LINEAR_ELAST_ELEMENT>
680 for (
unsigned j = 0;
j <
nnod;
j++)
705 error_stream <<
"Elements only store two fields so fld can't be"
706 <<
" " <<
fld << std::endl;
752 double interpolated_u = 0.0;
759 return interpolated_u;
766 return this->
nnode();
782 template<
class ELEMENT>
795 template<
class ELEMENT>
A base class for elements that solve the axisymmetric (in cylindrical polars) equations of linear ela...
double * Youngs_modulus_pt
Pointer to the Young's modulus.
double d2u_dt2_axisymmetric_linear_elasticity(const unsigned &n, const unsigned &i) const
d^2u/dt^2 at local node n
const double & lambda_sq() const
Access function for timescale ratio (nondim density)
unsigned ndof_types() const
The number of "DOF types" that degrees of freedom in this element are sub-divided into: for now lump ...
double *& nu_pt()
Access function for pointer to Poisson's ratio.
void interpolated_d2u_dt2_axisymmetric_linear_elasticity(const Vector< double > &s, Vector< double > &d2u_dt2) const
Compute vector of FE interpolated accel d2u/dt2 at local coordinate s.
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()
Access function: Pointer to body force function.
virtual unsigned u_index_axisymmetric_linear_elasticity(const unsigned &i) const
Return the index at which the i-th (i=0: r, i=1: z; i=2: theta) unknown displacement component is sto...
double * Nu_pt
Pointer to Poisson's ratio.
double youngs_modulus() const
Access function to Young's modulus.
double & nu() const
Access function for Poisson's ratio.
double *& lambda_sq_pt()
Access function for pointer to timescale ratio (nondim density)
static double Default_lambda_sq_value
Static default value for timescale ratio (1.0 for natural scaling)
void interpolated_du_dt_axisymmetric_linear_elasticity(const Vector< double > &s, Vector< double > &du_dt) const
Compute vector of FE interpolated velocity du/dt at local coordinate s.
double du_dt_axisymmetric_linear_elasticity(const unsigned &n, const unsigned &i) const
du/dt at local node n
AxisymmetricLinearElasticityEquationsBase()
Constructor: Set null pointers for constitutive law. Set physical parameter values to default values,...
double *& youngs_modulus_pt()
Return the pointer to Young's modulus.
static double Default_youngs_modulus_value
Static default value for Young's modulus (1.0 – for natural scaling, i.e. all stresses have been non-...
void body_force(const double &time, const Vector< double > &x, Vector< double > &b) const
Evaluate body force at Eulerian coordinate x at present time (returns zero vector if no body force fu...
BodyForceFctPt body_force_fct_pt() const
Access function: Pointer to body force function (const version)
BodyForceFctPt Body_force_fct_pt
Pointer to body force function.
double * Lambda_sq_pt
Timescale ratio (non-dim. density)
double interpolated_u_axisymmetric_linear_elasticity(const Vector< double > &s, const unsigned &i) const
Return FE interpolated displacement u[i] (i=0: r, i=1: z; i=2: theta) at local coordinate s.
void(* BodyForceFctPt)(const double &time, const Vector< double > &x, Vector< double > &b)
Function pointer to function that specifies the body force as a function of the Cartesian coordinates...
void interpolated_u_axisymmetric_linear_elasticity(const Vector< double > &s, Vector< double > &disp) const
Compute vector of FE interpolated displacement u at local coordinate s.
A class for elements that solve the axisymmetric (in cylindrical polars) equations of linear elastici...
void get_strain(const Vector< double > &s, DenseMatrix< double > &strain)
Get strain (3x3 entries; r, z, phi)
virtual void fill_in_generic_contribution_to_residuals_axisymmetric_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(std::ostream &outfile)
Output: r,z, u_r, u_z, u_theta.
AxisymmetricLinearElasticityEquations()
Constructor.
void output(FILE *file_pt)
C-style output: r,z, u_r, u_z, u_theta.
unsigned required_nvalue(const unsigned &n) const
Number of values required at node n.
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 compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Validate against exact solution. Solution is provided via function pointer. Plot at a given number of...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Return the residuals for the equations (the discretised principle of virtual displacements)
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact solution: r,z, u_r, u_z, u_theta.
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 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 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.
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...
TimeStepper *& time_stepper_pt()
Access function for pointer to time stepper: Null if object is not time-dependent.
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....
Axisym linear elasticity upgraded to become projectable.
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...
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 ...
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
unsigned nhistory_values_for_coordinate_projection()
Number of positional history values: Read out from positional timestepper (Note: count includes curre...
ProjectableAxisymLinearElasticityElement()
Constructor [this was only required explicitly from gcc 4.5.2 onwards...].
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_projection(const unsigned &fld)
Number of history values to be stored for fld-th field. (includes present value!)
unsigned nfields_for_projection()
Number of fields to be projected: 3, corresponding to the displacement components.
Wrapper class for projectable elements. Adds "projectability" to the underlying ELEMENT.
An Element that solves the equations of axisymmetric (in cylindrical polars) linear elasticity,...
void output(std::ostream &outfile, const unsigned &n_plot)
Output function.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function.
QAxisymmetricLinearElasticityElement()
Constructor.
void output(FILE *file_pt)
C-style output function.
void output(std::ostream &outfile)
Output function.
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...
Base class for time-stepping schemes. Timestepper provides an approximation of the temporal derivativ...
unsigned ntstorage() const
Return the number of doubles required to represent history (one for steady)
virtual double weight(const unsigned &i, const unsigned &j) const
Access function for j-th weight for the i-th derivative.
bool is_steady() const
Flag to indicate if a timestepper has been made steady (possibly temporarily to switch off time-depen...
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).