29 #ifndef OOMPH_LINEAR_ELASTICITY_TRACTION_ELEMENTS_HEADER
30 #define OOMPH_LINEAR_ELASTICITY_TRACTION_ELEMENTS_HEADER
34 #include <oomph-lib-config.h>
39 #include "../generic/Qelements.h"
47 namespace LinearElasticityTractionElementHelper
53 const Vector<double>& x,
54 const Vector<double>&
N,
55 Vector<double>& load);
67 template<
class ELEMENT>
92 const unsigned& intpt,
123 ELEMENT* elem_pt =
dynamic_cast<ELEMENT*
>(element_pt);
125 if (elem_pt->dim() == 3)
134 throw OomphLibError(
"This flux element will not work correctly "
135 "if nodes are hanging\n",
136 OOMPH_CURRENT_FUNCTION,
137 OOMPH_EXCEPTION_LOCATION);
148 ELEMENT* cast_element_pt =
dynamic_cast<ELEMENT*
>(element_pt);
149 this->U_index_linear_elasticity_traction.resize(n_dim);
150 for (
unsigned i = 0;
i < n_dim;
i++)
152 this->U_index_linear_elasticity_traction[
i] =
153 cast_element_pt->u_index_linear_elasticity(
i);
194 const unsigned&
i)
const
206 void output(std::ostream& outfile,
const unsigned& n_plot)
218 void output(FILE* file_pt,
const unsigned& n_plot)
241 template<
class ELEMENT>
245 unsigned n_dim = this->nodal_dimension();
249 interpolated_x(
s, x);
253 outer_unit_normal(
s, unit_normal);
259 get_traction(time, ipt, x, unit_normal, traction);
266 template<
class ELEMENT>
272 unsigned n_node = nnode();
275 double time = node_pt(0)->time_stepper_pt()->time_pt()->time();
279 unsigned n_position_type = this->nnodal_position_type();
280 if (n_position_type != 1)
282 throw OomphLibError(
"LinearElasticity is not yet implemented for more "
283 "than one position type",
284 OOMPH_CURRENT_FUNCTION,
285 OOMPH_EXCEPTION_LOCATION);
290 unsigned n_dim = this->nodal_dimension();
293 unsigned u_nodal_index[n_dim];
294 for (
unsigned i = 0;
i < n_dim;
i++)
296 u_nodal_index[
i] = this->U_index_linear_elasticity_traction[
i];
306 DShape dpsids(n_node, n_dim - 1);
309 unsigned n_intpt = integral_pt()->nweight();
312 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
315 double w = integral_pt()->weight(ipt);
318 dshape_local_at_knot(ipt, psi, dpsids);
327 for (
unsigned l = 0; l < n_node; l++)
330 for (
unsigned i = 0;
i < n_dim;
i++)
333 const double x_local = nodal_position(l,
i);
334 interpolated_x[
i] += x_local * psi(l);
337 for (
unsigned j = 0; j < n_dim - 1; j++)
339 interpolated_A(j,
i) += x_local * dpsids(l, j);
346 for (
unsigned i = 0;
i < n_dim - 1;
i++)
348 for (
unsigned j = 0; j < n_dim - 1; j++)
354 for (
unsigned k = 0; k < n_dim; k++)
356 A(
i, j) += interpolated_A(
i, k) * interpolated_A(j, k);
363 outer_unit_normal(ipt, interpolated_normal);
373 Adet = A(0, 0) * A(1, 1) - A(0, 1) * A(1, 0);
377 "Wrong dimension in LinearElasticityTractionElement",
378 "LinearElasticityTractionElement::fill_in_contribution_to_"
380 OOMPH_EXCEPTION_LOCATION);
385 double W = w * sqrt(Adet);
389 get_traction(time, ipt, interpolated_x, interpolated_normal, traction);
392 for (
unsigned l = 0; l < n_node; l++)
395 for (
unsigned i = 0;
i < n_dim;
i++)
397 local_eqn = this->nodal_local_eqn(l, u_nodal_index[
i]);
402 residuals[local_eqn] -= traction[
i] * psi(l) *
W;
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
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)
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,...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
A general Finite Element class.
virtual void output(std::ostream &outfile)
Output the element data — typically the values at the nodes in a format suitable for post-processing.
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.
bool has_hanging_nodes() const
Return boolean to indicate if any of the element's nodes are geometrically hanging.
A class for elements that allow the imposition of an applied traction in the equations of linear elas...
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function.
void output(std::ostream &outfile)
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 fill_in_contribution_to_residuals(Vector< double > &residuals)
Return the residuals.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function.
void fill_in_contribution_to_residuals_linear_elasticity_traction(Vector< double > &residuals)
Helper function that actually calculates the residuals.
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...
void output(FILE *file_pt)
C_style output function.
void(*&)(const double &time, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction) traction_fct_pt()
Reference to the traction function pointer.
LinearElasticityTractionElement(FiniteElement *const &element_pt, const int &face_index)
Constructor, which takes a "bulk" element and the value of the index and its limit.
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), Eulerlian coordinate and normal ve...
Vector< unsigned > U_index_linear_elasticity_traction
Index at which the i-th displacement component is stored.
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...
An OomphLibError object which should be thrown when an run-time error is encountered....
RefineableElements are FiniteElements that may be subdivided into children to provide a better local ...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
void Zero_traction_fct(const double &time, const Vector< double > &x, const Vector< double > &N, Vector< double > &load)
Default load function (zero traction)
//////////////////////////////////////////////////////////////////// ////////////////////////////////...