35 template<
unsigned DIM>
48 template<
unsigned DIM>
55 std::ostringstream error_message;
56 error_message <<
"Strain matrix is " <<
strain.ncol() <<
" x "
57 <<
strain.nrow() <<
", but dimension of the equations is "
66 if (n_position_type != 1)
68 throw OomphLibError(
"TimeHarmonicLinearElasticity is not yet implemented "
69 "for more than one position type",
81 for (
unsigned i = 0;
i <
DIM;
i++)
95 DIM,
DIM, std::complex<double>(0.0, 0.0));
101 for (
unsigned i = 0;
i <
DIM;
i++)
104 const std::complex<double>
u_value =
109 for (
unsigned j = 0;
j <
DIM;
j++)
117 for (
unsigned i = 0;
i <
DIM;
i++)
122 for (
int j = (
DIM - 1);
j >=
static_cast<int>(
i);
j--)
125 if (
static_cast<int>(
i) !=
j)
128 0.5 * (interpolated_dudx(
i,
j) + interpolated_dudx(
j,
i));
137 for (
int j = (
i - 1);
j >= 0;
j--)
149 template<
unsigned DIM>
156 std::ostringstream error_message;
157 error_message <<
"Stress matrix is " <<
stress.ncol() <<
" x "
158 <<
stress.nrow() <<
", but dimension of the equations is "
172 for (
unsigned i = 0;
i <
DIM;
i++)
174 for (
unsigned j = 0;
j <
DIM;
j++)
177 for (
unsigned k = 0;
k <
DIM;
k++)
179 for (
unsigned l = 0;
l <
DIM;
l++)
193 template<
unsigned DIM>
205 if (n_position_type != 1)
207 throw OomphLibError(
"TimeHarmonicLinearElasticity is not yet implemented "
208 "for more than one position type",
214 if (this->Elasticity_tensor_pt == 0)
224 for (
unsigned i = 0;
i <
DIM;
i++)
250 for (
unsigned i = 0;
i <
DIM; ++
i)
266 DIM, std::complex<double>(0.0, 0.0));
271 DIM,
DIM, std::complex<double>(0.0, 0.0));
277 for (
unsigned i = 0;
i <
DIM;
i++)
284 const std::complex<double>
u_value = std::complex<double>(
291 for (
unsigned j = 0;
j <
DIM;
j++)
311 for (
unsigned a = 0; a <
DIM; a++)
325 for (
unsigned b = 0; b <
DIM; b++)
327 for (
unsigned c = 0; c <
DIM; c++)
329 for (
unsigned d = 0; d <
DIM; d++)
333 interpolated_dudx(c, d).real() *
346 for (
unsigned c = 0; c <
DIM; c++)
361 for (
unsigned b = 0; b <
DIM; b++)
363 for (
unsigned d = 0; d <
DIM; d++)
391 for (
unsigned b = 0; b <
DIM; b++)
393 for (
unsigned c = 0; c <
DIM; c++)
395 for (
unsigned d = 0; d <
DIM; d++)
399 interpolated_dudx(c, d).imag() *
412 for (
unsigned c = 0; c <
DIM; c++)
427 for (
unsigned b = 0; b <
DIM; b++)
429 for (
unsigned d = 0; d <
DIM; d++)
452 template<
unsigned DIM>
455 const unsigned&
nplot,
484 for (
unsigned i = 0;
i <
DIM;
i++)
488 for (
unsigned i = 0;
i < 2 *
DIM;
i++)
503 template<
unsigned DIM>
505 const unsigned&
nplot)
524 this->interpolated_u_time_harmonic_linear_elasticity(
s, u);
527 for (
unsigned i = 0;
i <
DIM;
i++)
533 for (
unsigned i = 0;
i <
DIM;
i++)
539 for (
unsigned i = 0;
i <
DIM;
i++)
555 template<
unsigned DIM>
557 const unsigned&
nplot)
573 for (
unsigned i = 0;
i <
DIM;
i++)
579 for (
unsigned i = 0;
i <
DIM;
i++)
584 this->interpolated_u_time_harmonic_linear_elasticity(
s,
i).
real());
586 for (
unsigned i = 0;
i <
DIM;
i++)
591 this->interpolated_u_time_harmonic_linear_elasticity(
s,
i).
imag());
604 template<
unsigned DIM>
631 for (
unsigned i = 0;
i <
DIM;
i++)
646 this->interpolated_u_time_harmonic_linear_elasticity(
s,
disp);
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
Class for dense matrices, storing all the values of the matrix as a pointer to a pointer with assorte...
virtual double J_eulerian(const Vector< double > &s) const
Return the Jacobian of mapping from local to global coordinates at local position s.
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
double nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n. Produces suitably interpolated values for hanging nodes...
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Return string for tecplot zone header (when plotting nplot points in each "coordinate direction")
virtual double dshape_eulerian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidx) const
Return the geometric shape functions and also first derivatives w.r.t. global coordinates at the ipt-...
unsigned nnodal_position_type() const
Return the number of coordinate types that the element requires to interpolate the geometry between t...
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
unsigned nnode() const
Return the number of nodes.
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Get cector of local coordinates of plot point i (when plotting nplot points in each "coordinate direc...
virtual unsigned nplot_points(const unsigned &nplot) const
Return total number of plot points (when plotting nplot points in each "coordinate direction")
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Compute the geometric shape functions and also first derivatives w.r.t. global coordinates at local c...
double raw_nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n but do NOT take hanging nodes into account.
double raw_nodal_position(const unsigned &n, const unsigned &i) const
Return the i-th coordinate at local node n. Do not use the hanging node representation....
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Add tecplot zone "footer" to output stream (when plotting nplot points in each "coordinate direction"...
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.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
An OomphLibError object which should be thrown when an run-time error is encountered....
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...
TAdvectionDiffusionReactionElement()
Constructor: Call constructors for TElement and AdvectionDiffusionReaction equations.
static double Default_omega_sq_value
Static default value for square of frequency.
void get_strain(const Vector< double > &s, DenseMatrix< std::complex< double > > &strain) const
Return the strain tensor.
virtual void fill_in_generic_contribution_to_residuals_time_harmonic_linear_elasticity(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Private helper function to compute residuals and (if requested via flag) also the Jacobian matrix.
void get_stress(const Vector< double > &s, DenseMatrix< std::complex< double > > &sigma) const
Return the Cauchy stress tensor, as calculated from the elasticity tensor at specified local coordina...
void compute_norm(double &norm)
Compute norm of solution: square of the L2 norm.
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact solution x,y,[z],u_r,v_r,[w_r],u_i,v_i,[w_i].
void output(std::ostream &outfile)
Output: x,y,[z],u_r,v_r,[w_r],u_i,v_i,[w_i].
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).