29 #ifndef OOMPH_PML_TIME_HARMONIC_LINEAR_ELASTICITY_TRACTION_ELEMENTS_HEADER
30 #define OOMPH_PML_TIME_HARMONIC_LINEAR_ELASTICITY_TRACTION_ELEMENTS_HEADER
34 #include <oomph-lib-config.h>
39 #include "../generic/Qelements.h"
47 namespace PMLTimeHarmonicLinearElasticityTractionElementHelper
56 unsigned n_dim =
load.size();
57 for (
unsigned i = 0;
i < n_dim;
i++)
59 load[
i] = std::complex<double>(0.0, 0.0);
73 template<
class ELEMENT>
130 ELEMENT* cast_element_pt =
dynamic_cast<ELEMENT*
>(element_pt);
132 for (
unsigned i = 0;
i < n_dim;
i++)
135 cast_element_pt->u_index_time_harmonic_linear_elasticity(
i);
145 ELEMENT* elem_pt =
dynamic_cast<ELEMENT*
>(element_pt);
147 if (elem_pt->dim() == 3)
156 throw OomphLibError(
"This flux element will not work correctly "
157 "if nodes are hanging\n",
201 const unsigned&
i)
const
214 void output(std::ostream& outfile,
const unsigned& nplot)
226 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
236 for (
unsigned i = 0;
i <
ndim + 1;
i++)
238 outfile <<
x[
i] <<
" ";
242 for (
unsigned i = 0;
i <
ndim + 1;
i++)
248 for (
unsigned i = 0;
i <
ndim + 1;
i++)
253 outfile << std::endl;
267 void output(FILE* file_pt,
const unsigned& n_plot)
289 template<
class ELEMENT>
293 unsigned n_dim = this->nodal_dimension();
297 interpolated_x(
s,
x);
301 outer_unit_normal(
s, unit_normal);
307 get_traction(ipt,
x, unit_normal, traction);
315 template<
class ELEMENT>
321 unsigned n_node = nnode();
325 unsigned n_position_type = this->nnodal_position_type();
326 if (n_position_type != 1)
328 throw OomphLibError(
"PMLTimeHarmonicLinearElasticity is not yet "
329 "implemented for more than one position type",
336 unsigned n_dim = this->nodal_dimension();
339 std::vector<std::complex<unsigned>> u_nodal_index(n_dim);
340 for (
unsigned i = 0;
i < n_dim;
i++)
349 this->U_index_time_harmonic_linear_elasticity_traction[
i];
359 DShape dpsids(n_node, n_dim - 1);
362 unsigned n_intpt = integral_pt()->nweight();
365 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
368 double w = integral_pt()->weight(ipt);
371 dshape_local_at_knot(ipt, psi, dpsids);
380 for (
unsigned l = 0; l < n_node; l++)
383 for (
unsigned i = 0;
i < n_dim;
i++)
386 const double x_local = nodal_position(l,
i);
387 interpolated_x[
i] += x_local * psi(l);
390 for (
unsigned j = 0;
j < n_dim - 1;
j++)
392 interpolated_A(
j,
i) += x_local * dpsids(l,
j);
399 for (
unsigned i = 0;
i < n_dim - 1;
i++)
401 for (
unsigned j = 0;
j < n_dim - 1;
j++)
407 for (
unsigned k = 0;
k < n_dim;
k++)
409 A(
i,
j) += interpolated_A(
i,
k) * interpolated_A(
j,
k);
416 outer_unit_normal(ipt, interpolated_normal);
426 Adet =
A(0, 0) *
A(1, 1) -
A(0, 1) *
A(1, 0);
430 "Wrong dimension in PMLTimeHarmonicLinearElasticityTractionElement",
431 "PMLTimeHarmonicLinearElasticityTractionElement::fill_in_"
432 "contribution_to_residuals()",
438 double W =
w *
sqrt(Adet);
442 get_traction(ipt, interpolated_x, interpolated_normal, traction);
445 for (
unsigned l = 0; l < n_node; l++)
448 for (
unsigned i = 0;
i < n_dim;
i++)
451 local_eqn = this->nodal_local_eqn(l, u_nodal_index[
i].
real());
456 residuals[local_eqn] -= traction[
i].real() * psi(l) *
W;
460 local_eqn = this->nodal_local_eqn(l, u_nodal_index[
i].
imag());
465 residuals[local_eqn] -= traction[
i].imag() * psi(l) *
W;
AnnoyingScalar imag(const AnnoyingScalar &)
Definition: AnnoyingScalar.h:132
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
double interpolated_x(const Vector< double > &s, const unsigned &i) const
Definition: elements.h:4528
Definition: elements.h:4998
Definition: elements.h:1313
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
unsigned dim() const
Definition: elements.h:2611
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
unsigned nodal_dimension() const
Return the required Eulerian dimension of the nodes in this element.
Definition: elements.h:2484
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Definition: elements.h:3174
bool has_hanging_nodes() const
Definition: elements.h:2470
unsigned ndim() const
Access function to # of Eulerian coordinates.
Definition: geom_objects.h:177
Definition: oomph_definitions.h:222
Definition: pml_time_harmonic_linear_elasticity_traction_elements.h:77
void(* Traction_fct_pt)(const Vector< double > &x, const Vector< double > &n, Vector< std::complex< double >> &result)
Definition: pml_time_harmonic_linear_elasticity_traction_elements.h:88
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: pml_time_harmonic_linear_elasticity_traction_elements.h:199
void traction(const Vector< double > &s, Vector< std::complex< double >> &traction)
Definition: pml_time_harmonic_linear_elasticity_traction_elements.h:290
PMLTimeHarmonicLinearElasticityTractionElement(FiniteElement *const &element_pt, const int &face_index)
Definition: pml_time_harmonic_linear_elasticity_traction_elements.h:118
void fill_in_contribution_to_residuals_time_harmonic_linear_elasticity_traction(Vector< double > &residuals)
Helper function that actually calculates the residuals.
Definition: pml_time_harmonic_linear_elasticity_traction_elements.h:317
void output(std::ostream &outfile)
Output function.
Definition: pml_time_harmonic_linear_elasticity_traction_elements.h:207
void(*&)(const Vector< double > &x, const Vector< double > &n, Vector< std::complex< double >> &traction) traction_fct_pt()
Reference to the traction function pointer.
Definition: pml_time_harmonic_linear_elasticity_traction_elements.h:169
Vector< std::complex< unsigned > > U_index_time_harmonic_linear_elasticity_traction
Index at which the i-th displacement component is stored.
Definition: pml_time_harmonic_linear_elasticity_traction_elements.h:81
void output(FILE *file_pt)
C_style output function.
Definition: pml_time_harmonic_linear_elasticity_traction_elements.h:261
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function.
Definition: pml_time_harmonic_linear_elasticity_traction_elements.h:267
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
Definition: pml_time_harmonic_linear_elasticity_traction_elements.h:186
void output(std::ostream &outfile, const unsigned &nplot)
Output function.
Definition: pml_time_harmonic_linear_elasticity_traction_elements.h:214
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Return the residuals.
Definition: pml_time_harmonic_linear_elasticity_traction_elements.h:178
virtual void get_traction(const unsigned &intpt, const Vector< double > &x, const Vector< double > &n, Vector< std::complex< double >> &traction)
Definition: pml_time_harmonic_linear_elasticity_traction_elements.h:98
Definition: refineable_elements.h:97
@ N
Definition: constructor.cpp:22
float real
Definition: datatypes.h:10
RealScalar s
Definition: level1_cplx_impl.h:130
char char char int int * k
Definition: level2_impl.h:374
void Zero_traction_fct(const Vector< double > &x, const Vector< double > &N, Vector< std::complex< double >> &load)
Default load function (zero traction)
Definition: pml_time_harmonic_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