29 #ifndef OOMPH_FOURIER_DECOMPOSED_TIME_HARMONIC_LINEAR_ELASTICITY_TRACTION_ELEMENTS_HEADER
30 #define OOMPH_FOURIER_DECOMPOSED_TIME_HARMONIC_LINEAR_ELASTICITY_TRACTION_ELEMENTS_HEADER
34 #include <oomph-lib-config.h>
48 namespace TimeHarmonicFourierDecomposedLinearElasticityTractionElementHelper
57 unsigned n_dim =
load.size();
58 for (
unsigned i = 0;
i < n_dim;
i++)
60 load[
i] = std::complex<double>(0.0, 0.0);
75 template<
class ELEMENT>
132 ELEMENT* cast_element_pt =
dynamic_cast<ELEMENT*
>(element_pt);
135 for (
unsigned i = 0;
i < n_dim + 1;
i++)
148 ->U_index_time_harmonic_fourier_decomposed_linear_elasticity_traction
151 ->u_index_time_harmonic_fourier_decomposed_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)
240 template<
class ELEMENT>
241 void TimeHarmonicFourierDecomposedLinearElasticityTractionElement<
243 Vector<std::complex<double>>& traction)
245 unsigned n_dim = this->nodal_dimension();
249 interpolated_x(
s,
x);
253 outer_unit_normal(
s, unit_normal);
259 get_traction(ipt,
x, unit_normal, traction);
267 template<
class ELEMENT>
273 unsigned n_node = nnode();
277 unsigned n_position_type = this->nnodal_position_type();
278 if (n_position_type != 1)
280 throw OomphLibError(
"TimeHarmonicFourierDecomposedLinearElasticity is "
281 "not yet implemented for more than one position type",
288 const unsigned n_dim = this->nodal_dimension();
291 std::vector<std::complex<unsigned>> u_nodal_index(n_dim + 1);
292 for (
unsigned i = 0;
i < n_dim + 1;
i++)
304 ->U_index_time_harmonic_fourier_decomposed_linear_elasticity_traction
315 DShape dpsids(n_node, n_dim - 1);
318 unsigned n_intpt = integral_pt()->nweight();
321 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
324 double w = integral_pt()->weight(ipt);
327 dshape_local_at_knot(ipt, psi, dpsids);
336 for (
unsigned l = 0; l < n_node; l++)
339 for (
unsigned i = 0;
i < n_dim;
i++)
342 const double x_local = nodal_position(l,
i);
343 interpolated_x[
i] += x_local * psi(l);
346 for (
unsigned j = 0;
j < n_dim - 1;
j++)
348 interpolated_A(
j,
i) += x_local * dpsids(l,
j);
355 for (
unsigned i = 0;
i < n_dim - 1;
i++)
357 for (
unsigned j = 0;
j < n_dim - 1;
j++)
363 for (
unsigned k = 0;
k < n_dim;
k++)
365 A(
i,
j) += interpolated_A(
i,
k) * interpolated_A(
j,
k);
372 outer_unit_normal(ipt, interpolated_normal);
382 Adet =
A(0, 0) *
A(1, 1) -
A(0, 1) *
A(1, 0);
386 "Wrong dimension in "
387 "TimeHarmonicFourierDecomposedLinearElasticityTractionElement",
388 "TimeHarmonicFourierDecomposedLinearElasticityTractionElement::"
389 "fill_in_contribution_to_residuals()",
395 double W =
w *
sqrt(Adet);
399 get_traction(ipt, interpolated_x, interpolated_normal, traction);
402 for (
unsigned l = 0; l < n_node; l++)
405 for (
unsigned i = 0;
i < n_dim + 1;
i++)
408 local_eqn = this->nodal_local_eqn(l, u_nodal_index[
i].
real());
413 residuals[local_eqn] -=
414 traction[
i].real() * psi(l) * interpolated_x[0] *
W;
418 local_eqn = this->nodal_local_eqn(l, u_nodal_index[
i].
imag());
423 residuals[local_eqn] -=
424 traction[
i].imag() * psi(l) * interpolated_x[0] *
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
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
Definition: oomph_definitions.h:222
Definition: time_harmonic_fourier_decomposed_linear_elasticity_traction_elements.h:79
void traction(const Vector< double > &s, Vector< std::complex< double >> &traction)
Definition: time_harmonic_fourier_decomposed_linear_elasticity_traction_elements.h:242
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function.
Definition: time_harmonic_fourier_decomposed_linear_elasticity_traction_elements.h:218
TimeHarmonicFourierDecomposedLinearElasticityTractionElement(FiniteElement *const &element_pt, const int &face_index)
Definition: time_harmonic_fourier_decomposed_linear_elasticity_traction_elements.h:120
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
Definition: time_harmonic_fourier_decomposed_linear_elasticity_traction_elements.h:179
void output(std::ostream &outfile)
Output function.
Definition: time_harmonic_fourier_decomposed_linear_elasticity_traction_elements.h:200
Vector< std::complex< unsigned > > U_index_time_harmonic_fourier_decomposed_linear_elasticity_traction
Index at which the i-th displacement component is stored.
Definition: time_harmonic_fourier_decomposed_linear_elasticity_traction_elements.h:83
void(*&)(const Vector< double > &x, const Vector< double > &n, Vector< std::complex< double >> &traction) traction_fct_pt()
Reference to the traction function pointer.
Definition: time_harmonic_fourier_decomposed_linear_elasticity_traction_elements.h:162
void fill_in_contribution_to_residuals_time_harmonic_fourier_decomposed_linear_elasticity_traction(Vector< double > &residuals)
Helper function that actually calculates the residuals.
Definition: time_harmonic_fourier_decomposed_linear_elasticity_traction_elements.h:269
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: time_harmonic_fourier_decomposed_linear_elasticity_traction_elements.h:192
void output(FILE *file_pt)
C_style output function.
Definition: time_harmonic_fourier_decomposed_linear_elasticity_traction_elements.h:212
void output(std::ostream &outfile, const unsigned &n_plot)
Output function.
Definition: time_harmonic_fourier_decomposed_linear_elasticity_traction_elements.h:206
virtual void get_traction(const unsigned &intpt, const Vector< double > &x, const Vector< double > &n, Vector< std::complex< double >> &traction)
Definition: time_harmonic_fourier_decomposed_linear_elasticity_traction_elements.h:100
void(* Traction_fct_pt)(const Vector< double > &x, const Vector< double > &n, Vector< std::complex< double >> &result)
Definition: time_harmonic_fourier_decomposed_linear_elasticity_traction_elements.h:90
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Return the residuals.
Definition: time_harmonic_fourier_decomposed_linear_elasticity_traction_elements.h:171
@ 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
@ W
Definition: quadtree.h:63
void Zero_traction_fct(const Vector< double > &x, const Vector< double > &N, Vector< std::complex< double >> &load)
Default load function (zero traction)
Definition: time_harmonic_fourier_decomposed_linear_elasticity_traction_elements.h:53
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