26#ifndef OOMPH_FD_HH_TIME_HARMONIC_LIN_ELAST_HEADER
27#define OOMPH_FD_HH_TIME_HARMONIC_LIN_ELAST_HEADER
32#include "fourier_decomposed_helmholtz.h"
33#include "time_harmonic_fourier_decomposed_linear_elasticity.h"
46 template<
class ELASTICITY_BULK_ELEMENT,
class HELMHOLTZ_BULK_ELEMENT>
101 throw OomphLibError(
"This flux element will not work correctly "
102 "if nodes are hanging\n",
122 ->U_index_fourier_decomposed_time_harmonic_linear_elasticity_helmholtz_traction
124 for (
unsigned i = 0;
i <
n_dim + 1;
i++)
127 ->U_index_fourier_decomposed_time_harmonic_linear_elasticity_helmholtz_traction
130 ->u_index_time_harmonic_fourier_decomposed_linear_elasticity(
i);
158 const double&
q()
const
221 for (
unsigned i = 0;
i <
n_dim;
i++)
230 <<
zeta[0] << std::endl;
256 template<
class ELASTICITY_BULK_ELEMENT,
class HELMHOLTZ_BULK_ELEMENT>
258 FourierDecomposedTimeHarmonicLinElastLoadedByHelmholtzPressureBCElement<
259 ELASTICITY_BULK_ELEMENT,
260 HELMHOLTZ_BULK_ELEMENT>::Default_Q_Value = 1.0;
266 template<
class ELASTICITY_BULK_ELEMENT,
class HELMHOLTZ_BULK_ELEMENT>
267 void FourierDecomposedTimeHarmonicLinElastLoadedByHelmholtzPressureBCElement<
268 ELASTICITY_BULK_ELEMENT,
270 fill_in_contribution_to_residuals_helmholtz_traction(
274 unsigned n_node = nnode();
281 throw OomphLibError(
"LinearElasticity is not yet implemented for more "
282 "than one position type.",
289 const unsigned n_dim = this->nodal_dimension();
293 for (
unsigned i = 0;
i <
n_dim + 1;
i++)
297 ->U_index_fourier_decomposed_time_harmonic_linear_elasticity_helmholtz_traction
315 unsigned n_intpt = integral_pt()->nweight();
321 double w = integral_pt()->weight(
ipt);
336 for (
unsigned i = 0;
i <
n_dim;
i++)
339 const double x_local = nodal_position(
l,
i);
343 for (
unsigned j = 0;
j <
n_dim - 1;
j++)
351 double r = interpolated_x[0];
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++)
382 Adet = A(0, 0) * A(1, 1) - A(0, 1) * A(1, 0);
386 "Wrong dimension in TimeHarmonicLinElastLoadedByPressureElement",
387 "TimeHarmonicLinElastLoadedByPressureElement::fill_in_contribution_"
404 for (
unsigned i = 0;
i <
n_dim;
i++)
412 for (
unsigned i = 0;
i <
n_dim + 1;
i++)
450 template<
class HELMHOLTZ_BULK_ELEMENT,
class ELASTICITY_BULK_ELEMENT>
524 for (
unsigned i = 0;
i < (
Dim - 1);
i++)
541 ->interpolated_u_time_harmonic_fourier_decomposed_linear_elasticity(
556 <<
" " << flux <<
" " <<
zeta[0] << std::endl;
630 const unsigned&
flag);
649 template<
class HELMHOLTZ_BULK_ELEMENT,
class ELASTICITY_BULK_ELEMENT>
653 FourierDecomposedHelmholtzFluxFromNormalDisplacementBCElement(
677 throw OomphLibError(
"This flux element will not work correctly if "
678 "nodes are hanging\n",
714 "Bulk element must inherit from HelmholtzEquations.";
716 "Nodes are three dimensional, but cannot cast the bulk element to\n";
719 "If you desire this functionality, you must implement it yourself\n";
728 eqn_pt->u_index_fourier_decomposed_helmholtz();
737 template<
class HELMHOLTZ_BULK_ELEMENT,
class ELASTICITY_BULK_ELEMENT>
741 fill_in_generic_residual_contribution_helmholtz_flux_from_displacement(
744 const unsigned&
flag)
747 const unsigned n_node = nnode();
753 const unsigned n_intpt = integral_pt()->nweight();
762 const std::complex<unsigned> u_index_helmholtz =
763 U_index_helmholtz_from_displacement;
770 for (
unsigned i = 0;
i < (Dim - 1);
i++)
772 s[
i] = integral_pt()->knot(
ipt,
i);
776 double w = integral_pt()->weight(
ipt);
792 for (
unsigned i = 0;
i < Dim;
i++)
794 interpolated_x[
i] += nodal_position(
l,
i) *
psif[
l];
799 double r = interpolated_x[0];
810 ->interpolated_u_time_harmonic_fourier_decomposed_linear_elasticity(
816 std::complex<double> flux =
824 local_eqn = nodal_local_eqn(
l, u_index_helmholtz.real());
836 local_eqn = nodal_local_eqn(
l, u_index_helmholtz.imag());
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
This is a base class for all elements that require external sources (e.g. FSI, multi-domain problems ...
Vector< double > & external_element_local_coord(const unsigned &interaction_index, const unsigned &ipt)
Access function to get source element's local coords for specified interaction index at specified int...
void set_ninteraction(const unsigned &n_interaction)
Set the number of interactions in the element This function is usually called in the specific element...
void fill_in_jacobian_from_external_interaction_by_fd(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Calculate the contributions to the jacobian from all external interaction degrees of freedom (geometr...
FiniteElement *& external_element_pt(const unsigned &interaction_index, const unsigned &ipt)
Access function to source element for specified interaction index at specified integration point.
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)
void outer_unit_normal(const Vector< double > &s, Vector< double > &unit_normal) const
Compute outer unit normal at the specified local coordinate.
double J_eulerian(const Vector< double > &s) const
Return the Jacobian of mapping from local to global coordinates at local position s....
double J_eulerian_at_knot(const unsigned &ipt) const
Return the Jacobian of the mapping from local to global coordinates at the ipt-th integration point O...
FaceGeometry class definition: This policy class is used to allow construction of face elements that ...
A general Finite Element class.
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
void interpolated_zeta(const Vector< double > &s, Vector< double > &zeta) const
Calculate the interpolated value of zeta, the intrinsic coordinate of the element when viewed as a co...
virtual void shape(const Vector< double > &s, Shape &psi) const =0
Calculate the geometric shape functions at local coordinate s. This function must be overloaded for e...
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
unsigned nnode() const
Return the number of nodes.
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.
virtual void shape_at_knot(const unsigned &ipt, Shape &psi) const
Return the geometric shape function at the ipt-th integration point.
bool has_hanging_nodes() const
Return boolean to indicate if any of the element's nodes are geometrically hanging.
A class for all isoparametric elements that solve the Helmholtz equations.
A class for elements that allow the imposition of an prescribed flux (determined from the normal disp...
void fill_in_generic_residual_contribution_helmholtz_flux_from_displacement(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Add the element's contribution to its residual vector. flag=1(or 0): do (or don't) compute the contri...
FourierDecomposedHelmholtzFluxFromNormalDisplacementBCElement(const FourierDecomposedHelmholtzFluxFromNormalDisplacementBCElement &dummy)=delete
Broken copy constructor.
void output(std::ostream &outfile)
Output function.
std::complex< unsigned > U_index_helmholtz_from_displacement
The index at which the unknown is stored at the nodes.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Add the element's contribution to its residual vector and its Jacobian matrix.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Broken assignment operator.
void output(FILE *file_pt)
C-style output function.
double shape_and_test_at_knot(const unsigned &ipt, Shape &psi, Shape &test) const
Function to compute the shape and test functions and to return the Jacobian of mapping between local ...
unsigned Dim
The spatial dimension of the problem.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: flux etc at Gauss points; n_plot is ignored.
double shape_and_test(const Vector< double > &s, Shape &psi, Shape &test) const
Function to compute the shape and test functions and to return the Jacobian of mapping between local ...
A class for elements that allow the imposition of an applied traction in the equations of time-harmon...
void output(FILE *file_pt)
C_style output function.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: Plot traction etc at Gauss points nplot is ignored.
double *& q_pt()
Return a pointer the ratio of stress scales used to non-dimensionalise the fluid and solid equations.
Vector< std::complex< unsigned > > U_index_fourier_decomposed_time_harmonic_linear_elasticity_helmholtz_traction
Index at which the i-th displacement component is stored.
static double Default_Q_Value
Static default value for the ratio of stress scales used in the fluid and solid equations (default is...
void fill_in_contribution_to_residuals_helmholtz_traction(Vector< double > &residuals)
Helper function that actually calculates the residuals.
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.
const double & q() const
Return the ratio of the stress scales used to non-dimensionalise the fluid and elasticity equations....
double * Q_pt
Pointer to the ratio, , of the stress used to non-dimensionalise the fluid stresses to the stress us...
FourierDecomposedTimeHarmonicLinElastLoadedByHelmholtzPressureBCElement(FiniteElement *const &element_pt, const int &face_index)
Constructor, which takes a "bulk" element and the value of the index and its limit.
void output(std::ostream &outfile)
Output function.
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
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.
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...
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).