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
60 unsigned n_dim =
load.size();
61 for (
unsigned i = 0;
i < n_dim;
i++)
77 template<
class ELEMENT>
103 const unsigned& intpt,
135 ELEMENT* cast_element_pt =
dynamic_cast<ELEMENT*
>(element_pt);
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] =
140 cast_element_pt->u_index_axisymmetric_linear_elasticity(
i);
183 const unsigned&
i)
const
196 void output(std::ostream& outfile,
const unsigned& n_plot)
202 const unsigned n_node =
nnode();
219 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
232 for (
unsigned i = 0;
i < n_dim;
i++)
238 for (
unsigned l = 0; l < n_node; l++)
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++)
270 outfile << unit_normal[
i] <<
" ";
272 outfile << std::endl;
283 void output(FILE* file_pt,
const unsigned& n_plot)
306 template<
class ELEMENT>
310 unsigned n_dim = this->nodal_dimension();
314 interpolated_x(
s,
x);
318 outer_unit_normal(
s, unit_normal);
324 get_traction(time, ipt,
x, unit_normal, traction);
332 template<
class ELEMENT>
338 unsigned n_node = nnode();
341 double time = node_pt(0)->time_stepper_pt()->time_pt()->time();
345 unsigned n_position_type = this->nnodal_position_type();
346 if (n_position_type != 1)
348 throw OomphLibError(
"AxisymmetricLinearElasticity is not yet implemented "
349 "for more than one position type",
356 unsigned n_dim = this->nodal_dimension();
359 unsigned u_nodal_index[n_dim + 1];
360 for (
unsigned i = 0;
i < n_dim + 1;
i++)
363 this->U_index_axisymmetric_linear_elasticity_traction[
i];
373 DShape dpsids(n_node, n_dim - 1);
376 unsigned n_intpt = integral_pt()->nweight();
379 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
382 double w = integral_pt()->weight(ipt);
385 dshape_local_at_knot(ipt, psi, dpsids);
394 for (
unsigned l = 0; l < n_node; l++)
397 for (
unsigned i = 0;
i < n_dim;
i++)
400 const double x_local = nodal_position(l,
i);
401 interpolated_x[
i] += x_local * psi(l);
404 for (
unsigned j = 0;
j < n_dim - 1;
j++)
406 interpolated_A(
j,
i) += x_local * dpsids(l,
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++)
423 A(
i,
j) += interpolated_A(
i,
k) * interpolated_A(
j,
k);
430 outer_unit_normal(ipt, interpolated_normal);
440 Adet =
A(0, 0) *
A(1, 1) -
A(0, 1) *
A(1, 0);
444 "Wrong dimension in AxisymmetricLinearElasticityTractionElement",
445 "AxisymmetricLinearElasticityTractionElement::fill_in_contribution_"
452 double W =
w *
sqrt(Adet);
456 get_traction(time, ipt, interpolated_x, interpolated_normal, traction);
459 for (
unsigned l = 0; l < n_node; l++)
462 for (
unsigned i = 0;
i < n_dim + 1;
i++)
465 local_eqn = this->nodal_local_eqn(l, u_nodal_index[
i]);
470 residuals[local_eqn] -=
471 traction[
i] * psi(l) * interpolated_x[0] *
W;
492 template<
class ELASTICITY_BULK_ELEMENT,
class NAVIER_STOKES_BULK_ELEMENT>
539 ELASTICITY_BULK_ELEMENT* cast_element_pt =
540 dynamic_cast<ELASTICITY_BULK_ELEMENT*
>(element_pt);
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] =
545 cast_element_pt->u_index_axisymmetric_linear_elasticity(
i);
573 const double&
q()
const
596 void output(std::ostream& outfile,
const unsigned& n_plot)
602 const double q_value =
q();
613 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
627 NAVIER_STOKES_BULK_ELEMENT* ext_el_pt =
628 dynamic_cast<NAVIER_STOKES_BULK_ELEMENT*
>(
633 ext_el_pt->traction(s_ext, interpolated_normal, traction);
635 outfile << ext_el_pt->interpolated_x(s_ext, 0) <<
" "
636 << ext_el_pt->interpolated_x(s_ext, 1) <<
" "
637 << q_value * traction[0] <<
" " << q_value * traction[1] <<
" "
638 << q_value * traction[2] <<
" " << interpolated_normal[0] <<
" "
639 << interpolated_normal[1] <<
" " <<
zeta[0] << std::endl;
650 void output(FILE* file_pt,
const unsigned& n_plot)
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,
677 NAVIER_STOKES_BULK_ELEMENT>::
678 fill_in_contribution_to_residuals_axisym_fsi_traction(
682 unsigned n_node = nnode();
686 unsigned n_position_type = this->nnodal_position_type();
687 if (n_position_type != 1)
689 throw OomphLibError(
"LinearElasticity is not yet implemented for more "
690 "than one position type.",
697 unsigned n_dim = this->nodal_dimension();
700 unsigned u_nodal_index[n_dim + 1];
701 for (
unsigned i = 0;
i < n_dim + 1;
i++)
703 u_nodal_index[
i] = this->U_index_axisym_fsi_traction[
i];
711 DShape dpsids(n_node, n_dim - 1);
714 const double q_value =
q();
720 unsigned n_intpt = integral_pt()->nweight();
723 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
726 double w = integral_pt()->weight(ipt);
729 dshape_local_at_knot(ipt, psi, dpsids);
738 for (
unsigned l = 0; l < n_node; l++)
741 for (
unsigned i = 0;
i < n_dim;
i++)
744 const double x_local = nodal_position(l,
i);
745 interpolated_x[
i] += x_local * psi(l);
748 for (
unsigned j = 0;
j < n_dim - 1;
j++)
750 interpolated_A(
j,
i) += x_local * dpsids(l,
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++)
767 A(
i,
j) += interpolated_A(
i,
k) * interpolated_A(
j,
k);
774 outer_unit_normal(ipt, interpolated_normal);
784 Adet =
A(0, 0) *
A(1, 1) -
A(0, 1) *
A(1, 0);
788 "Wrong dimension in TimeHarmonicLinElastLoadedByPressureElement",
789 "TimeHarmonicLinElastLoadedByPressureElement::fill_in_contribution_"
796 double W =
w *
sqrt(Adet);
799 NAVIER_STOKES_BULK_ELEMENT* ext_el_pt =
800 dynamic_cast<NAVIER_STOKES_BULK_ELEMENT*
>(
801 this->external_element_pt(0, ipt));
805 ext_el_pt->traction(s_ext, interpolated_normal, traction);
808 for (
unsigned l = 0; l < n_node; l++)
811 for (
unsigned i = 0;
i < n_dim + 1;
i++)
814 local_eqn = this->nodal_local_eqn(l, u_nodal_index[
i]);
820 residuals[local_eqn] -=
821 q_value * traction[
i] * psi(l) * interpolated_x[0] *
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: axisym_linear_elasticity_traction_elements.h:81
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function.
Definition: axisym_linear_elasticity_traction_elements.h:283
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: axisym_linear_elasticity_traction_elements.h:181
void traction(const double &time, const Vector< double > &s, Vector< double > &traction)
Definition: axisym_linear_elasticity_traction_elements.h:307
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
Definition: axisym_linear_elasticity_traction_elements.h:168
Vector< unsigned > U_index_axisymmetric_linear_elasticity_traction
Index at which the i-th displacement component is stored.
Definition: axisym_linear_elasticity_traction_elements.h:84
void output(FILE *file_pt)
C_style output function.
Definition: axisym_linear_elasticity_traction_elements.h:277
void(* Traction_fct_pt)(const double &time, const Vector< double > &x, const Vector< double > &n, Vector< double > &result)
Definition: axisym_linear_elasticity_traction_elements.h:91
virtual void get_traction(const double &time, const unsigned &intpt, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction)
Definition: axisym_linear_elasticity_traction_elements.h:102
void output(std::ostream &outfile)
Output function.
Definition: axisym_linear_elasticity_traction_elements.h:189
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: axisym_linear_elasticity_traction_elements.h:150
AxisymmetricLinearElasticityTractionElement(FiniteElement *const &element_pt, const int &face_index)
Definition: axisym_linear_elasticity_traction_elements.h:123
void fill_in_contribution_to_residuals_axisymmetric_linear_elasticity_traction(Vector< double > &residuals)
Helper function that actually calculates the residuals.
Definition: axisym_linear_elasticity_traction_elements.h:334
void output(std::ostream &outfile, const unsigned &n_plot)
Output function.
Definition: axisym_linear_elasticity_traction_elements.h:196
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Return the residuals.
Definition: axisym_linear_elasticity_traction_elements.h:160
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
Definition: element_with_external_element.h:56
Vector< double > & external_element_local_coord(const unsigned &interaction_index, const unsigned &ipt)
Definition: element_with_external_element.h:136
void set_ninteraction(const unsigned &n_interaction)
Definition: element_with_external_element.h:178
FiniteElement *& external_element_pt(const unsigned &interaction_index, const unsigned &ipt)
Definition: element_with_external_element.h:107
void fill_in_jacobian_from_external_interaction_by_fd(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: element_with_external_element.h:345
Definition: axisym_linear_elasticity_traction_elements.h:497
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Return the residuals.
Definition: axisym_linear_elasticity_traction_elements.h:551
FSIAxisymmetricLinearElasticityTractionElement(FiniteElement *const &element_pt, const int &face_index)
Definition: axisym_linear_elasticity_traction_elements.h:521
double *& q_pt()
Definition: axisym_linear_elasticity_traction_elements.h:580
Vector< unsigned > U_index_axisym_fsi_traction
Index at which the i-th displacement component is stored.
Definition: axisym_linear_elasticity_traction_elements.h:509
void fill_in_contribution_to_residuals_axisym_fsi_traction(Vector< double > &residuals)
Helper function that actually calculates the residuals.
Definition: axisym_linear_elasticity_traction_elements.h:678
const double & q() const
Definition: axisym_linear_elasticity_traction_elements.h:573
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function.
Definition: axisym_linear_elasticity_traction_elements.h:650
void output(std::ostream &outfile, const unsigned &n_plot)
Definition: axisym_linear_elasticity_traction_elements.h:596
void output(FILE *file_pt)
C_style output function.
Definition: axisym_linear_elasticity_traction_elements.h:644
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
Definition: axisym_linear_elasticity_traction_elements.h:558
void output(std::ostream &outfile)
Output function.
Definition: axisym_linear_elasticity_traction_elements.h:587
static double Default_Q_Value
Definition: axisym_linear_elasticity_traction_elements.h:506
double * Q_pt
Definition: axisym_linear_elasticity_traction_elements.h:502
Definition: elements.h:4338
int & face_index()
Definition: elements.h:4626
void outer_unit_normal(const Vector< double > &s, Vector< double > &unit_normal) const
Compute outer unit normal at the specified local coordinate.
Definition: elements.cc:6006
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: elements.h:4497
double interpolated_x(const Vector< double > &s, const unsigned &i) const
Definition: elements.h:4528
Definition: elements.h:4998
Definition: elements.h:1313
void interpolated_zeta(const Vector< double > &s, Vector< double > &zeta) const
Definition: elements.cc:4675
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
virtual void output(std::ostream &outfile)
Definition: elements.h:3050
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Definition: elements.h:3161
virtual void build_face_element(const int &face_index, FaceElement *face_element_pt)
Definition: elements.cc:5132
virtual void shape(const Vector< double > &s, Shape &psi) const =0
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1963
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Definition: elements.h:3148
virtual unsigned nplot_points(const unsigned &nplot) const
Definition: elements.h:3186
double nodal_position(const unsigned &n, const unsigned &i) const
Definition: elements.h:2317
unsigned nodal_dimension() const
Return the required Eulerian dimension of the nodes in this element.
Definition: elements.h:2484
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.
Definition: oomph_definitions.h:222
Time *const & time_pt() const
Access function for the pointer to time (const version)
Definition: timesteppers.h:572
double & time()
Return the current value of the continuous time.
Definition: timesteppers.h:123
@ N
Definition: constructor.cpp:22
RealScalar s
Definition: level1_cplx_impl.h:130
char char char int int * k
Definition: level2_impl.h:374
EIGEN_DEVICE_FUNC const Scalar & q
Definition: SpecialFunctionsImpl.h:2019
EIGEN_STRONG_INLINE const Eigen::CwiseBinaryOp< Eigen::internal::scalar_zeta_op< typename DerivedX::Scalar >, const DerivedX, const DerivedQ > zeta(const Eigen::ArrayBase< DerivedX > &x, const Eigen::ArrayBase< DerivedQ > &q)
Definition: SpecialFunctionsArrayAPI.h:152
void Zero_traction_fct(const double &time, const Vector< double > &x, const Vector< double > &N, Vector< double > &load)
Default load function (zero traction)
Definition: axisym_linear_elasticity_traction_elements.h:55
@ 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
void output(std::ostream &outfile, const unsigned &nplot)
Overload output function.
Definition: overloaded_element_body.h:490
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2