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
57 unsigned n_dim =
load.size();
58 for (
unsigned i = 0;
i < n_dim;
i++)
74 template<
class ELEMENT>
99 const unsigned& intpt,
130 ELEMENT* elem_pt =
dynamic_cast<ELEMENT*
>(element_pt);
132 if (elem_pt->dim() == 3)
141 throw OomphLibError(
"This flux element will not work correctly "
142 "if nodes are hanging\n",
155 ELEMENT* cast_element_pt =
dynamic_cast<ELEMENT*
>(element_pt);
156 this->U_index_linear_elasticity_traction.resize(n_dim);
157 for (
unsigned i = 0;
i < n_dim;
i++)
159 this->U_index_linear_elasticity_traction[
i] =
160 cast_element_pt->u_index_linear_elasticity(
i);
201 const unsigned&
i)
const
213 void output(std::ostream& outfile,
const unsigned& n_plot)
225 void output(FILE* file_pt,
const unsigned& n_plot)
248 template<
class ELEMENT>
252 unsigned n_dim = this->nodal_dimension();
256 interpolated_x(
s,
x);
260 outer_unit_normal(
s, unit_normal);
266 get_traction(time, ipt,
x, unit_normal, traction);
273 template<
class ELEMENT>
279 unsigned n_node = nnode();
282 double time = node_pt(0)->time_stepper_pt()->time_pt()->time();
286 unsigned n_position_type = this->nnodal_position_type();
287 if (n_position_type != 1)
289 throw OomphLibError(
"LinearElasticity is not yet implemented for more "
290 "than one position type",
297 unsigned n_dim = this->nodal_dimension();
300 unsigned u_nodal_index[n_dim];
301 for (
unsigned i = 0;
i < n_dim;
i++)
303 u_nodal_index[
i] = this->U_index_linear_elasticity_traction[
i];
313 DShape dpsids(n_node, n_dim - 1);
316 unsigned n_intpt = integral_pt()->nweight();
319 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
322 double w = integral_pt()->weight(ipt);
325 dshape_local_at_knot(ipt, psi, dpsids);
334 for (
unsigned l = 0; l < n_node; l++)
337 for (
unsigned i = 0;
i < n_dim;
i++)
340 const double x_local = nodal_position(l,
i);
341 interpolated_x[
i] += x_local * psi(l);
344 for (
unsigned j = 0;
j < n_dim - 1;
j++)
346 interpolated_A(
j,
i) += x_local * dpsids(l,
j);
353 for (
unsigned i = 0;
i < n_dim - 1;
i++)
355 for (
unsigned j = 0;
j < n_dim - 1;
j++)
361 for (
unsigned k = 0;
k < n_dim;
k++)
363 A(
i,
j) += interpolated_A(
i,
k) * interpolated_A(
j,
k);
370 outer_unit_normal(ipt, interpolated_normal);
380 Adet =
A(0, 0) *
A(1, 1) -
A(0, 1) *
A(1, 0);
384 "Wrong dimension in LinearElasticityTractionElement",
385 "LinearElasticityTractionElement::fill_in_contribution_to_"
392 double W =
w *
sqrt(Adet);
396 get_traction(time, ipt, interpolated_x, interpolated_normal, traction);
399 for (
unsigned l = 0; l < n_node; l++)
402 for (
unsigned i = 0;
i < n_dim;
i++)
404 local_eqn = this->nodal_local_eqn(l, u_nodal_index[
i]);
409 residuals[local_eqn] -= traction[
i] * psi(l) *
W;
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
RowVector3d w
Definition: Matrix_resize_int.cpp:3
void load(Archive &ar, ParticleHandler &handl)
Definition: Particles.h:21
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
Definition: elements.h:4338
int & face_index()
Definition: elements.h:4626
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: elements.h:4497
Definition: elements.h:4998
Definition: elements.h:1313
virtual void output(std::ostream &outfile)
Definition: elements.h:3050
virtual void build_face_element(const int &face_index, FaceElement *face_element_pt)
Definition: elements.cc:5132
unsigned nodal_dimension() const
Return the required Eulerian dimension of the nodes in this element.
Definition: elements.h:2484
bool has_hanging_nodes() const
Definition: elements.h:2470
Definition: linear_elasticity_traction_elements.h:77
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function.
Definition: linear_elasticity_traction_elements.h:225
void output(std::ostream &outfile)
Output function.
Definition: linear_elasticity_traction_elements.h:207
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: linear_elasticity_traction_elements.h:199
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Return the residuals.
Definition: linear_elasticity_traction_elements.h:180
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
Definition: linear_elasticity_traction_elements.h:187
void output(std::ostream &outfile, const unsigned &n_plot)
Output function.
Definition: linear_elasticity_traction_elements.h:213
void fill_in_contribution_to_residuals_linear_elasticity_traction(Vector< double > &residuals)
Helper function that actually calculates the residuals.
Definition: linear_elasticity_traction_elements.h:275
void(* Traction_fct_pt)(const double &time, const Vector< double > &x, const Vector< double > &n, Vector< double > &result)
Definition: linear_elasticity_traction_elements.h:87
void output(FILE *file_pt)
C_style output function.
Definition: linear_elasticity_traction_elements.h:219
void(*&)(const double &time, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction) traction_fct_pt()
Reference to the traction function pointer.
Definition: linear_elasticity_traction_elements.h:170
LinearElasticityTractionElement(FiniteElement *const &element_pt, const int &face_index)
Definition: linear_elasticity_traction_elements.h:119
virtual void get_traction(const double &time, const unsigned &intpt, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction)
Definition: linear_elasticity_traction_elements.h:98
Vector< unsigned > U_index_linear_elasticity_traction
Index at which the i-th displacement component is stored.
Definition: linear_elasticity_traction_elements.h:80
void traction(const double &time, const Vector< double > &s, Vector< double > &traction)
Definition: linear_elasticity_traction_elements.h:249
Definition: oomph_definitions.h:222
Definition: refineable_elements.h:97
@ N
Definition: constructor.cpp:22
RealScalar s
Definition: level1_cplx_impl.h:130
char char char int int * k
Definition: level2_impl.h:374
void Zero_traction_fct(const double &time, const Vector< double > &x, const Vector< double > &N, Vector< double > &load)
Default load function (zero traction)
Definition: linear_elasticity_traction_elements.h:52
@ W
Definition: quadtree.h:63
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
list x
Definition: plotDoE.py:28
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2