35 template<
unsigned DIM>
46 template<
unsigned DIM>
53 std::ostringstream error_message;
54 error_message <<
"Strain matrix is " <<
strain.ncol() <<
" x "
55 <<
strain.nrow() <<
", but dimension of the equations is "
64 if (n_position_type != 1)
66 throw OomphLibError(
"LinearElasticity is not yet implemented for more "
67 "than one position type",
79 for (
unsigned i = 0;
i <
DIM;
i++)
98 for (
unsigned i = 0;
i <
DIM;
i++)
104 for (
unsigned j = 0;
j <
DIM;
j++)
112 for (
unsigned i = 0;
i <
DIM;
i++)
117 for (
int j = (
DIM - 1);
j >=
static_cast<int>(
i);
j--)
120 if (
static_cast<int>(
i) !=
j)
123 0.5 * (interpolated_dudx(
i,
j) + interpolated_dudx(
j,
i));
132 for (
int j = (
i - 1);
j >= 0;
j--)
144 template<
unsigned DIM>
151 std::ostringstream error_message;
152 error_message <<
"Stress matrix is " <<
stress.ncol() <<
" x "
153 <<
stress.nrow() <<
", but dimension of the equations is "
167 for (
unsigned i = 0;
i <
DIM;
i++)
169 for (
unsigned j = 0;
j <
DIM;
j++)
172 for (
unsigned k = 0;
k <
DIM;
k++)
174 for (
unsigned l = 0;
l <
DIM;
l++)
188 template<
unsigned DIM>
200 if (n_position_type != 1)
202 throw OomphLibError(
"LinearElasticity is not yet implemented for more "
203 "than one position type",
209 if (this->Elasticity_tensor_pt == 0)
219 for (
unsigned i = 0;
i <
DIM;
i++)
244 for (
unsigned i = 0;
i <
DIM; ++
i)
269 for (
unsigned i = 0;
i <
DIM;
i++)
277 accel[
i] += this->d2u_dt2_linear_elasticity(
l,
i) *
psi(
l);
284 for (
unsigned j = 0;
j <
DIM;
j++)
304 for (
unsigned a = 0; a <
DIM; a++)
315 for (
unsigned b = 0; b <
DIM; b++)
317 for (
unsigned c = 0; c <
DIM; c++)
319 for (
unsigned d = 0; d <
DIM; d++)
323 interpolated_dudx(c, d) *
336 for (
unsigned c = 0; c <
DIM; c++)
351 for (
unsigned b = 0; b <
DIM; b++)
353 for (
unsigned d = 0; d <
DIM; d++)
376 template<
unsigned DIM>
379 const unsigned&
nplot,
408 for (
unsigned i = 0;
i <
DIM;
i++)
412 for (
unsigned i = 0;
i <
DIM;
i++)
427 template<
unsigned DIM>
430 const unsigned&
nplot,
460 for (
unsigned i = 0;
i <
DIM;
i++)
464 for (
unsigned i = 0;
i <
DIM;
i++)
479 template<
unsigned DIM>
481 const unsigned&
nplot)
500 this->interpolated_u_linear_elasticity(
s, u);
503 for (
unsigned i = 0;
i <
DIM;
i++)
509 for (
unsigned i = 0;
i <
DIM;
i++)
525 template<
unsigned DIM>
527 const unsigned&
nplot)
543 for (
unsigned i = 0;
i <
DIM;
i++)
549 for (
unsigned i = 0;
i <
DIM;
i++)
567 template<
unsigned DIM>
586 outfile <<
"ZONE" << std::endl;
595 for (
unsigned i = 0;
i <
DIM;
i++)
616 for (
unsigned i = 0;
i <
DIM;
i++)
620 (
exact_soln[
i] - this->interpolated_u_linear_elasticity(
s,
i)) *
621 (
exact_soln[
i] - this->interpolated_u_linear_elasticity(
s,
i)) * W;
625 for (
unsigned i = 0;
i <
DIM;
i++)
631 for (
unsigned i = 0;
i <
DIM;
i++)
646 template<
unsigned DIM>
666 outfile <<
"ZONE" << std::endl;
675 for (
unsigned i = 0;
i <
DIM;
i++)
696 for (
unsigned i = 0;
i <
DIM;
i++)
700 (
exact_soln[
i] - this->interpolated_u_linear_elasticity(
s,
i)) *
701 (
exact_soln[
i] - this->interpolated_u_linear_elasticity(
s,
i)) * W;
705 for (
unsigned i = 0;
i <
DIM;
i++)
711 for (
unsigned i = 0;
i <
DIM;
i++)
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
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...
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
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"...
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as .
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.
static double Default_lambda_sq_value
Static default value for timescale ratio (1.0 – for natural scaling)
void get_strain(const Vector< double > &s, DenseMatrix< double > &strain) const
Return the strain tensor.
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact solution x,y,[z],u,v,[w].
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Validate against exact solution. Solution is provided via function pointer. Plot at a given number of...
void output(std::ostream &outfile)
Output: x,y,[z],u,v,[w].
void get_stress(const Vector< double > &s, DenseMatrix< double > &sigma) const
Return the Cauchy stress tensor, as calculated from the elasticity tensor at specified local coordina...
virtual void fill_in_generic_contribution_to_residuals_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.
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.
virtual double weight(const unsigned &i, const unsigned &j) const
Access function for j-th weight for the i-th derivative.
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).