29#ifndef OOMPH_AXISYMMETRIC_LINEAR_ELASTICITY_TRACTION_ELEMENTS_HEADER
30#define OOMPH_AXISYMMETRIC_LINEAR_ELASTICITY_TRACTION_ELEMENTS_HEADER
34#include <oomph-lib-config.h>
50 namespace AxisymmetricLinearElasticityTractionElementHelper
77 template<
class ELEMENT>
103 const unsigned&
intpt,
136 this->U_index_axisymmetric_linear_elasticity_traction.resize(
n_dim + 1);
137 for (
unsigned i = 0;
i <
n_dim + 1;
i++)
139 this->U_index_axisymmetric_linear_elasticity_traction[
i] =
183 const unsigned&
i)
const
232 for (
unsigned i = 0;
i <
n_dim;
i++)
241 for (
unsigned i = 0;
i <
n_dim;
i++)
256 for (
unsigned i = 0;
i <
n_dim;
i++)
262 for (
unsigned i = 0;
i <
n_dim + 1;
i++)
268 for (
unsigned i = 0;
i <
n_dim;
i++)
306 template<
class ELEMENT>
332 template<
class ELEMENT>
346 if (n_position_type != 1)
348 throw OomphLibError(
"AxisymmetricLinearElasticity is not yet implemented "
349 "for more than one position type",
360 for (
unsigned i = 0;
i <
n_dim + 1;
i++)
363 this->U_index_axisymmetric_linear_elasticity_traction[
i];
397 for (
unsigned i = 0;
i <
n_dim;
i++)
404 for (
unsigned j = 0;
j <
n_dim - 1;
j++)
413 for (
unsigned i = 0;
i <
n_dim - 1;
i++)
415 for (
unsigned j = 0;
j <
n_dim - 1;
j++)
421 for (
unsigned k = 0;
k <
n_dim;
k++)
440 Adet = A(0, 0) * A(1, 1) - A(0, 1) * A(1, 0);
444 "Wrong dimension in AxisymmetricLinearElasticityTractionElement",
445 "AxisymmetricLinearElasticityTractionElement::fill_in_contribution_"
462 for (
unsigned i = 0;
i <
n_dim + 1;
i++)
492 template<
class ELASTICITY_BULK_ELEMENT,
class NAVIER_STOKES_BULK_ELEMENT>
541 this->U_index_axisym_fsi_traction.resize(
n_dim + 1);
542 for (
unsigned i = 0;
i <
n_dim + 1;
i++)
544 this->U_index_axisym_fsi_traction[
i] =
573 const double&
q()
const
665 template<
class ELASTICITY_BULK_ELEMENT,
class NAVIER_STOKES_BULK_ELEMENT>
666 double FSIAxisymmetricLinearElasticityTractionElement<
667 ELASTICITY_BULK_ELEMENT,
668 NAVIER_STOKES_BULK_ELEMENT>::Default_Q_Value = 1.0;
674 template<
class ELASTICITY_BULK_ELEMENT,
class NAVIER_STOKES_BULK_ELEMENT>
675 void FSIAxisymmetricLinearElasticityTractionElement<
676 ELASTICITY_BULK_ELEMENT,
678 fill_in_contribution_to_residuals_axisym_fsi_traction(
682 unsigned n_node = nnode();
689 throw OomphLibError(
"LinearElasticity is not yet implemented for more "
690 "than one position type.",
697 unsigned n_dim = this->nodal_dimension();
701 for (
unsigned i = 0;
i <
n_dim + 1;
i++)
720 unsigned n_intpt = integral_pt()->nweight();
726 double w = integral_pt()->weight(
ipt);
741 for (
unsigned i = 0;
i <
n_dim;
i++)
744 const double x_local = nodal_position(
l,
i);
748 for (
unsigned j = 0;
j <
n_dim - 1;
j++)
757 for (
unsigned i = 0;
i <
n_dim - 1;
i++)
759 for (
unsigned j = 0;
j <
n_dim - 1;
j++)
765 for (
unsigned k = 0;
k <
n_dim;
k++)
784 Adet = A(0, 0) * A(1, 1) - A(0, 1) * A(1, 0);
788 "Wrong dimension in TimeHarmonicLinElastLoadedByPressureElement",
789 "TimeHarmonicLinElastLoadedByPressureElement::fill_in_contribution_"
801 this->external_element_pt(0,
ipt));
811 for (
unsigned i = 0;
i <
n_dim + 1;
i++)
A class for elements that allow the imposition of an applied traction in the equations of axisymmetri...
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function.
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Specify the value of nodal zeta from the face geometry The "global" intrinsic coordinate of the eleme...
void traction(const double &time, const Vector< double > &s, Vector< double > &traction)
Compute traction vector at specified local coordinate Should only be used for post-processing; ignore...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
Vector< unsigned > U_index_axisymmetric_linear_elasticity_traction
Index at which the i-th displacement component is stored.
void output(FILE *file_pt)
C_style output function.
void(* Traction_fct_pt)(const double &time, const Vector< double > &x, const Vector< double > &n, Vector< double > &result)
Pointer to an imposed traction function. Arguments: Eulerian coordinate; outer unit normal; applied t...
virtual void get_traction(const double &time, const unsigned &intpt, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction)
Get the traction vector: Pass number of integration point (dummy), Eulerian coordinate and normal vec...
void(*&)(const double &time, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction) traction_fct_pt()
Reference to the traction function pointer.
void output(std::ostream &outfile)
Output function.
AxisymmetricLinearElasticityTractionElement(FiniteElement *const &element_pt, const int &face_index)
Constructor, which takes a "bulk" element and the value of the index and its limit.
void fill_in_contribution_to_residuals_axisymmetric_linear_elasticity_traction(Vector< double > &residuals)
Helper function that actually calculates the residuals.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Return the residuals.
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.
This is a base class for all elements that require external sources (e.g. FSI, multi-domain problems ...
Vector< double > & external_element_local_coord(const unsigned &interaction_index, const unsigned &ipt)
Access function to get source element's local coords for specified interaction index at specified int...
void set_ninteraction(const unsigned &n_interaction)
Set the number of interactions in the element This function is usually called in the specific element...
void fill_in_jacobian_from_external_interaction_by_fd(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Calculate the contributions to the jacobian from all external interaction degrees of freedom (geometr...
FiniteElement *& external_element_pt(const unsigned &interaction_index, const unsigned &ipt)
Access function to source element for specified interaction index at specified integration point.
A class for elements that allow the imposition of an applied traction in the equations of axisymmetri...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Return the residuals.
FSIAxisymmetricLinearElasticityTractionElement(FiniteElement *const &element_pt, const int &face_index)
Constructor, which takes a "bulk" element and the value of the index and its limit.
const double & q() const
Return the ratio of the stress scales used to non-dimensionalise the fluid and elasticity equations....
Vector< unsigned > U_index_axisym_fsi_traction
Index at which the i-th displacement component is stored.
void fill_in_contribution_to_residuals_axisym_fsi_traction(Vector< double > &residuals)
Helper function that actually calculates the residuals.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: Plot traction etc at Gauss points nplot is ignored.
void output(FILE *file_pt)
C_style output function.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
double *& q_pt()
Return a pointer the ratio of stress scales used to non-dimensionalise the fluid and solid equations.
void output(std::ostream &outfile)
Output function.
static double Default_Q_Value
Static default value for the ratio of stress scales used in the fluid and solid equations (default is...
double * Q_pt
Pointer to the ratio, , of the stress used to non-dimensionalise the fluid stresses to the stress us...
FaceElements are elements that coincide with the faces of higher-dimensional "bulk" elements....
int & face_index()
Index of the face (a number that uniquely identifies the face in the element)
void outer_unit_normal(const Vector< double > &s, Vector< double > &unit_normal) const
Compute outer unit normal at the specified local coordinate.
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
In a FaceElement, the "global" intrinsic coordinate of the element along the boundary,...
double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s. Overloaded to get information from bulk...
FaceGeometry class definition: This policy class is used to allow construction of face elements that ...
A general Finite Element class.
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
virtual void dshape_local_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsids) const
Return the geometric shape function and its derivative w.r.t. the local coordinates at the ipt-th int...
virtual void output(std::ostream &outfile)
Output the element data — typically the values at the nodes in a format suitable for post-processing.
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Return string for tecplot zone header (when plotting nplot points in each "coordinate direction")
unsigned nnodal_position_type() const
Return the number of coordinate types that the element requires to interpolate the geometry between t...
void interpolated_zeta(const Vector< double > &s, Vector< double > &zeta) const
Calculate the interpolated value of zeta, the intrinsic coordinate of the element when viewed as a co...
double size() const
Calculate the size of the element (length, area, volume,...) in Eulerian computational coordinates....
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...
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
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 nnode() const
Return the number of nodes.
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Get cector of local coordinates of plot point i (when plotting nplot points in each "coordinate direc...
virtual unsigned nplot_points(const unsigned &nplot) const
Return total number of plot points (when plotting nplot points in each "coordinate direction")
double nodal_position(const unsigned &n, const unsigned &i) const
Return the i-th coordinate at local node n. If the node is hanging, the appropriate interpolation is ...
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
virtual void build_face_element(const int &face_index, FaceElement *face_element_pt)
Function for building a lower dimensional FaceElement on the specified face of the FiniteElement....
unsigned nodal_dimension() const
Return the required Eulerian dimension of the nodes in this element.
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
An OomphLibError object which should be thrown when an run-time error is encountered....
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...
TAdvectionDiffusionReactionElement()
Constructor: Call constructors for TElement and AdvectionDiffusionReaction equations.
Time *const & time_pt() const
Access function for the pointer to time (const version)
double & time()
Return the current value of the continuous time.
void Zero_traction_fct(const double &time, const Vector< double > &x, const Vector< double > &N, Vector< double > &load)
Default load function (zero traction)
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).