time_harmonic_linear_elasticity_elements.h
Go to the documentation of this file.
1 // LIC// ====================================================================
2 // LIC// This file forms part of oomph-lib, the object-oriented,
3 // LIC// multi-physics finite-element library, available
4 // LIC// at http://www.oomph-lib.org.
5 // LIC//
6 // LIC// Copyright (C) 2006-2024 Matthias Heil and Andrew Hazel
7 // LIC//
8 // LIC// This library is free software; you can redistribute it and/or
9 // LIC// modify it under the terms of the GNU Lesser General Public
10 // LIC// License as published by the Free Software Foundation; either
11 // LIC// version 2.1 of the License, or (at your option) any later version.
12 // LIC//
13 // LIC// This library is distributed in the hope that it will be useful,
14 // LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // LIC// Lesser General Public License for more details.
17 // LIC//
18 // LIC// You should have received a copy of the GNU Lesser General Public
19 // LIC// License along with this library; if not, write to the Free Software
20 // LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21 // LIC// 02110-1301 USA.
22 // LIC//
23 // LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24 // LIC//
25 // LIC//====================================================================
26 // Header file for general linear elasticity elements
27 
28 // Include guards to prevent multiple inclusion of the header
29 #ifndef OOMPH_TIME_HARMONIC_LINEAR_ELASTICITY_ELEMENTS_HEADER
30 #define OOMPH_TIME_HARMONIC_LINEAR_ELASTICITY_ELEMENTS_HEADER
31 
32 // Config header generated by autoconfig
33 #ifdef HAVE_CONFIG_H
34 #include <oomph-lib-config.h>
35 #endif
36 
37 
38 #ifdef OOMPH_HAS_MPI
39 #include "mpi.h"
40 #endif
41 
42 #include <complex>
43 
44 
45 // OOMPH-LIB headers
46 #include "../generic/Qelements.h"
47 #include "../generic/mesh.h"
48 #include "../generic/hermite_elements.h"
50 #include "../generic/projection.h"
51 
52 namespace oomph
53 {
54  //=======================================================================
55  /// A base class for elements that solve the equations of time-harmonic linear
56  /// elasticity in Cartesian coordinates.
57  /// Combines a few generic functions that are shared by
58  /// TimeHarmonicLinearElasticityEquations
59  /// and TimeHarmonicLinearElasticityEquationsWithPressure (Note: The latter
60  /// don't exist yet but will be written as soon as somebody needs them...)
61  //=======================================================================
62  template<unsigned DIM>
64  {
65  public:
66  /// Return the index at which the i-th real or imag unknown
67  /// displacement component is stored. The default value is appropriate for
68  /// single-physics problems:
69  virtual inline std::complex<unsigned> u_index_time_harmonic_linear_elasticity(
70  const unsigned i) const
71  {
72  return std::complex<unsigned>(i, i + DIM);
73  }
74 
75  /// Compute vector of FE interpolated displacement u at local coordinate s
77  const Vector<double>& s, Vector<std::complex<double>>& disp) const
78  {
79  // Find number of nodes
80  unsigned n_node = nnode();
81 
82  // Local shape function
83  Shape psi(n_node);
84 
85  // Find values of shape function
86  shape(s, psi);
87 
88  for (unsigned i = 0; i < DIM; i++)
89  {
90  // Index at which the nodal value is stored
91  std::complex<unsigned> u_nodal_index =
93 
94  // Initialise value of u
95  disp[i] = std::complex<double>(0.0, 0.0);
96 
97  // Loop over the local nodes and sum
98  for (unsigned l = 0; l < n_node; l++)
99  {
100  const std::complex<double> u_value(
101  nodal_value(l, u_nodal_index.real()),
102  nodal_value(l, u_nodal_index.imag()));
103 
104  disp[i] += u_value * psi[l];
105  }
106  }
107  }
108 
109  /// Return FE interpolated displacement u[i] at local coordinate s
111  const Vector<double>& s, const unsigned& i) const
112  {
113  // Find number of nodes
114  unsigned n_node = nnode();
115 
116  // Local shape function
117  Shape psi(n_node);
118 
119  // Find values of shape function
120  shape(s, psi);
121 
122  // Get nodal index at which i-th velocity is stored
123  std::complex<unsigned> u_nodal_index =
125 
126  // Initialise value of u
127  std::complex<double> interpolated_u(0.0, 0.0);
128 
129  // Loop over the local nodes and sum
130  for (unsigned l = 0; l < n_node; l++)
131  {
132  const std::complex<double> u_value(
133  nodal_value(l, u_nodal_index.real()),
134  nodal_value(l, u_nodal_index.imag()));
135 
136  interpolated_u += u_value * psi[l];
137  }
138 
139  return (interpolated_u);
140  }
141 
142 
143  /// Function pointer to function that specifies the body force
144  /// as a function of the Cartesian coordinates and time FCT(t,x,b) --
145  /// x and b are Vectors!
146  typedef void (*BodyForceFctPt)(const double& t,
147  const Vector<double>& x,
149 
150  /// Constructor: Set null pointers for constitutive law and for
151  /// isotropic growth function. Set physical parameter values to
152  /// default values, and set body force to zero.
157  {
158  }
159 
160  /// Return the pointer to the elasticity_tensor
162  {
163  return Elasticity_tensor_pt;
164  }
165 
166  /// Access function to the entries in the elasticity tensor
167  inline double E(const unsigned& i,
168  const unsigned& j,
169  const unsigned& k,
170  const unsigned& l) const
171  {
172  return (*Elasticity_tensor_pt)(i, j, k, l);
173  }
174 
175  /// Access function for square of non-dim frequency
176  const double& omega_sq() const
177  {
178  return *Omega_sq_pt;
179  }
180 
181  /// Access function for square of non-dim frequency
182  double*& omega_sq_pt()
183  {
184  return Omega_sq_pt;
185  }
186 
187  /// Access function: Pointer to body force function
189  {
190  return Body_force_fct_pt;
191  }
192 
193  /// Access function: Pointer to body force function (const version)
195  {
196  return Body_force_fct_pt;
197  }
198 
199  /// Return the Cauchy stress tensor, as calculated
200  /// from the elasticity tensor at specified local coordinate
201  /// Virtual so separaete versions can (and must!) be provided
202  /// for displacement and pressure-displacement formulations.
203  virtual void get_stress(const Vector<double>& s,
204  DenseMatrix<std::complex<double>>& sigma) const = 0;
205 
206  /// Return the strain tensor
207  void get_strain(const Vector<double>& s,
208  DenseMatrix<std::complex<double>>& strain) const;
209 
210  /// Evaluate body force at Eulerian coordinate x at present time
211  /// (returns zero vector if no body force function pointer has been set)
212  inline void body_force(const Vector<double>& x,
213  Vector<std::complex<double>>& b) const
214  {
215  // If no function has been set, return zero vector
216  if (Body_force_fct_pt == 0)
217  {
218  // Get spatial dimension of element
219  unsigned n = dim();
220  for (unsigned i = 0; i < n; i++)
221  {
222  b[i] = std::complex<double>(0.0, 0.0);
223  }
224  }
225  else
226  {
227  // Get time from timestepper of first node (note that this must
228  // work -- body force only makes sense for elements that can be
229  // deformed and thefore store displacements (at their nodes)
230  double time = node_pt(0)->time_stepper_pt()->time_pt()->time();
231 
232  // Get body force
233  (*Body_force_fct_pt)(time, x, b);
234  }
235  }
236 
237  /// The number of "DOF types" that degrees of freedom in this element
238  /// are sub-divided into: for now lump them all into one DOF.
239  /// Can be adjusted later
240  unsigned ndof_types() const
241  {
242  return 1;
243  }
244 
245  /// Create a list of pairs for all unknowns in this element,
246  /// so that the first entry in each pair contains the global equation
247  /// number of the unknown, while the second one contains the number
248  /// of the "DOF type" that this unknown is associated with.
249  /// (Function can obviously only be called if the equation numbering
250  /// scheme has been set up.)
252  std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list) const
253  {
254  // temporary pair (used to store DOF lookup prior to being added
255  // to list)
256  std::pair<unsigned long, unsigned> dof_lookup;
257 
258  // number of nodes
259  const unsigned n_node = this->nnode();
260 
261  // Integer storage for local unknown
262  int local_unknown = 0;
263 
264  // Loop over the nodes
265  for (unsigned n = 0; n < n_node; n++)
266  {
267  // Loop over dimension (real and imag)
268  for (unsigned i = 0; i < 2 * DIM; i++)
269  {
270  // If the variable is free
271  local_unknown = nodal_local_eqn(n, i);
272 
273  // ignore pinned values
274  if (local_unknown >= 0)
275  {
276  // store DOF lookup in temporary pair: First entry in pair
277  // is global equation number; second entry is DOF type
278  dof_lookup.first = this->eqn_number(local_unknown);
279  dof_lookup.second = 0;
280 
281  // add to list
282  dof_lookup_list.push_front(dof_lookup);
283  }
284  }
285  }
286  }
287 
288 
289  protected:
290  /// Pointer to the elasticity tensor
292 
293  /// Square of nondim frequency
294  double* Omega_sq_pt;
295 
296  /// Pointer to body force function
298 
299  /// Static default value for square of frequency
300  static double Default_omega_sq_value;
301  };
302 
303 
304  /// ////////////////////////////////////////////////////////////////////
305  /// ////////////////////////////////////////////////////////////////////
306  /// ////////////////////////////////////////////////////////////////////
307 
308 
309  //=======================================================================
310  /// A class for elements that solve the equations of linear elasticity
311  /// in cartesian coordinates.
312  //=======================================================================
313  template<unsigned DIM>
316  {
317  public:
318  /// Constructor
320 
321  /// Number of values required at node n.
322  unsigned required_nvalue(const unsigned& n) const
323  {
324  return 2 * DIM;
325  }
326 
327  /// Return the residuals for the solid equations (the discretised
328  /// principle of virtual displacements)
330  {
332  residuals, GeneralisedElement::Dummy_matrix, 0);
333  }
334 
335 
336  /// The jacobian is calculated by finite differences by default,
337  /// We need only to take finite differences w.r.t. positional variables
338  /// For this element
340  DenseMatrix<double>& jacobian)
341  {
342  // Add the contribution to the residuals
343  this
344  ->fill_in_generic_contribution_to_residuals_time_harmonic_linear_elasticity(
345  residuals, jacobian, 1);
346  }
347 
348  /// Return the Cauchy stress tensor, as calculated
349  /// from the elasticity tensor at specified local coordinate
350  void get_stress(const Vector<double>& s,
351  DenseMatrix<std::complex<double>>& sigma) const;
352 
353  /// Output exact solution x,y,[z],u_r,v_r,[w_r],u_i,v_i,[w_i]
354  void output_fct(std::ostream& outfile,
355  const unsigned& nplot,
357 
358  /// Output: x,y,[z],u_r,v_r,[w_r],u_i,v_i,[w_i]
359  void output(std::ostream& outfile)
360  {
361  unsigned n_plot = 5;
362  output(outfile, n_plot);
363  }
364 
365  /// Output: x,y,[z],u_r,v_r,[w_r],u_i,v_i,[w_i]
366  void output(std::ostream& outfile, const unsigned& n_plot);
367 
368 
369  /// C-style output: x,y,[z],u_r,v_r,[w_r],u_i,v_i,[w_i]
370  void output(FILE* file_pt)
371  {
372  unsigned n_plot = 5;
373  output(file_pt, n_plot);
374  }
375 
376  /// Output: x,y,[z],u_r,v_r,[w_r],u_i,v_i,[w_i]
377  void output(FILE* file_pt, const unsigned& n_plot);
378 
379 
380  /// Compute norm of solution: square of the L2 norm
381  void compute_norm(double& norm);
382 
383  private:
384  /// Private helper function to compute residuals and (if requested
385  /// via flag) also the Jacobian matrix.
387  Vector<double>& residuals, DenseMatrix<double>& jacobian, unsigned flag);
388  };
389 
390 
391  /// /////////////////////////////////////////////////////////////////////
392  /// /////////////////////////////////////////////////////////////////////
393  /// /////////////////////////////////////////////////////////////////////
394 
395 
396  //===========================================================================
397  /// An Element that solves the equations of linear elasticity
398  /// in Cartesian coordinates, using QElements for the geometry
399  //============================================================================
400  template<unsigned DIM, unsigned NNODE_1D>
402  : public virtual QElement<DIM, NNODE_1D>,
403  public virtual TimeHarmonicLinearElasticityEquations<DIM>
404  {
405  public:
406  /// Constructor
408  : QElement<DIM, NNODE_1D>(), TimeHarmonicLinearElasticityEquations<DIM>()
409  {
410  }
411 
412  /// Output function
413  void output(std::ostream& outfile)
414  {
416  }
417 
418  /// Output function
419  void output(std::ostream& outfile, const unsigned& n_plot)
420  {
422  }
423 
424 
425  /// C-style output function
426  void output(FILE* file_pt)
427  {
429  }
430 
431  /// C-style output function
432  void output(FILE* file_pt, const unsigned& n_plot)
433  {
435  }
436  };
437 
438 
439  //============================================================================
440  /// FaceGeometry of a linear 2D QTimeHarmonicLinearElasticityElement element
441  //============================================================================
442  template<>
444  : public virtual QElement<1, 2>
445  {
446  public:
447  /// Constructor must call the constructor of the underlying solid element
448  FaceGeometry() : QElement<1, 2>() {}
449  };
450 
451 
452  //============================================================================
453  /// FaceGeometry of a quadratic 2D QTimeHarmonicLinearElasticityElement
454  /// element
455  //============================================================================
456  template<>
458  : public virtual QElement<1, 3>
459  {
460  public:
461  /// Constructor must call the constructor of the underlying element
462  FaceGeometry() : QElement<1, 3>() {}
463  };
464 
465 
466  //============================================================================
467  /// FaceGeometry of a cubic 2D QTimeHarmonicLinearElasticityElement element
468  //============================================================================
469  template<>
471  : public virtual QElement<1, 4>
472  {
473  public:
474  /// Constructor must call the constructor of the underlying element
475  FaceGeometry() : QElement<1, 4>() {}
476  };
477 
478 
479  //============================================================================
480  /// FaceGeometry of a linear 3D QTimeHarmonicLinearElasticityElement element
481  //============================================================================
482  template<>
484  : public virtual QElement<2, 2>
485  {
486  public:
487  /// Constructor must call the constructor of the underlying element
488  FaceGeometry() : QElement<2, 2>() {}
489  };
490 
491  //============================================================================
492  /// FaceGeometry of a quadratic 3D QTimeHarmonicLinearElasticityElement
493  /// element
494  //============================================================================
495  template<>
497  : public virtual QElement<2, 3>
498  {
499  public:
500  /// Constructor must call the constructor of the underlying element
501  FaceGeometry() : QElement<2, 3>() {}
502  };
503 
504 
505  //============================================================================
506  /// FaceGeometry of a cubic 3D QTimeHarmonicLinearElasticityElement element
507  //============================================================================
508  template<>
510  : public virtual QElement<2, 4>
511  {
512  public:
513  /// Constructor must call the constructor of the underlying element
514  FaceGeometry() : QElement<2, 4>() {}
515  };
516 
517  /// /////////////////////////////////////////////////////////////////
518  /// /////////////////////////////////////////////////////////////////
519  /// /////////////////////////////////////////////////////////////////
520 
521 
522  //==========================================================
523  /// Time-harmonic linear elasticity upgraded to become projectable
524  //==========================================================
525  template<class TIME_HARMONIC_LINEAR_ELAST_ELEMENT>
527  : public virtual ProjectableElement<TIME_HARMONIC_LINEAR_ELAST_ELEMENT>
528  {
529  public:
530  /// Constructor [this was only required explicitly
531  /// from gcc 4.5.2 onwards...]
533 
534 
535  /// Specify the values associated with field fld.
536  /// The information is returned in a vector of pairs which comprise
537  /// the Data object and the value within it, that correspond to field fld.
538  /// In the underlying time-harmonic linear elasticity elements the
539  /// real and complex parts of the displacements are stored
540  /// at the nodal values
542  {
543  // Create the vector
545 
546  // Loop over all nodes and extract the fld-th nodal value
547  unsigned nnod = this->nnode();
548  for (unsigned j = 0; j < nnod; j++)
549  {
550  // Add the data value associated with the displacement components
551  data_values.push_back(std::make_pair(this->node_pt(j), fld));
552  }
553 
554  // Return the vector
555  return data_values;
556  }
557 
558  /// Number of fields to be projected: 2*dim, corresponding to
559  /// real and imag parts of the displacement components
561  {
562  return 2 * this->dim();
563  }
564 
565  /// Number of history values to be stored for fld-th field.
566  /// (includes present value!)
567  unsigned nhistory_values_for_projection(const unsigned& fld)
568  {
569 #ifdef PARANOID
570  if (fld > 3)
571  {
572  std::stringstream error_stream;
573  error_stream << "Elements only store four fields so fld can't be"
574  << " " << fld << std::endl;
575  throw OomphLibError(
576  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
577  }
578 #endif
579  return this->node_pt(0)->ntstorage();
580  }
581 
582  /// Number of positional history values: Read out from
583  /// positional timestepper
584  /// (Note: count includes current value!)
586  {
587  return this->node_pt(0)->position_time_stepper_pt()->ntstorage();
588  }
589 
590  /// Return Jacobian of mapping and shape functions of field fld
591  /// at local coordinate s
592  double jacobian_and_shape_of_field(const unsigned& fld,
593  const Vector<double>& s,
594  Shape& psi)
595  {
596  unsigned n_dim = this->dim();
597  unsigned n_node = this->nnode();
598  DShape dpsidx(n_node, n_dim);
599 
600  // Call the derivatives of the shape functions and return
601  // the Jacobian
602  return this->dshape_eulerian(s, psi, dpsidx);
603  }
604 
605 
606  /// Return interpolated field fld at local coordinate s, at time
607  /// level t (t=0: present; t>0: history values)
608  double get_field(const unsigned& t,
609  const unsigned& fld,
610  const Vector<double>& s)
611  {
612  unsigned n_node = this->nnode();
613 
614 #ifdef PARANOID
615  unsigned n_dim = this->node_pt(0)->ndim();
616 #endif
617 
618  // Local shape function
619  Shape psi(n_node);
620 
621  // Find values of shape function
622  this->shape(s, psi);
623 
624  // Initialise value of u
625  double interpolated_u = 0.0;
626 
627  // Sum over the local nodes at that time
628  for (unsigned l = 0; l < n_node; l++)
629  {
630 #ifdef PARANOID
631  unsigned nvalue = this->node_pt(l)->nvalue();
632  if (nvalue != 2 * n_dim)
633  {
634  std::stringstream error_stream;
635  error_stream
636  << "Current implementation only works for non-resized nodes\n"
637  << "but nvalue= " << nvalue << "!= 2 dim = " << 2 * n_dim
638  << std::endl;
639  throw OomphLibError(error_stream.str(),
640  OOMPH_CURRENT_FUNCTION,
641  OOMPH_EXCEPTION_LOCATION);
642  }
643 #endif
644  interpolated_u += this->nodal_value(t, l, fld) * psi[l];
645  }
646  return interpolated_u;
647  }
648 
649 
650  /// Return number of values in field fld
651  unsigned nvalue_of_field(const unsigned& fld)
652  {
653  return this->nnode();
654  }
655 
656 
657  /// Return local equation number of value j in field fld.
658  int local_equation(const unsigned& fld, const unsigned& j)
659  {
660 #ifdef PARANOID
661  unsigned n_dim = this->node_pt(0)->ndim();
662  unsigned nvalue = this->node_pt(j)->nvalue();
663  if (nvalue != 2 * n_dim)
664  {
665  std::stringstream error_stream;
666  error_stream
667  << "Current implementation only works for non-resized nodes\n"
668  << "but nvalue= " << nvalue << "!= 2 dim = " << 2 * n_dim
669  << std::endl;
670  throw OomphLibError(
671  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
672  }
673 #endif
674  return this->nodal_local_eqn(j, fld);
675  }
676  };
677 
678 
679  //=======================================================================
680  /// Face geometry for element is the same as that for the underlying
681  /// wrapped element
682  //=======================================================================
683  template<class ELEMENT>
685  : public virtual FaceGeometry<ELEMENT>
686  {
687  public:
688  FaceGeometry() : FaceGeometry<ELEMENT>() {}
689  };
690 
691 
692  //=======================================================================
693  /// Face geometry of the Face Geometry for element is the same as
694  /// that for the underlying wrapped element
695  //=======================================================================
696  template<class ELEMENT>
699  : public virtual FaceGeometry<FaceGeometry<ELEMENT>>
700  {
701  public:
703  };
704 
705 
706 } // namespace oomph
707 
708 #endif
static char t char * s
Definition: cfortran.h:568
cstr elem_len * i
Definition: cfortran.h:603
char t
Definition: cfortran.h:568
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
Definition: shape.h:278
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
unsigned nvalue() const
Return number of values stored in data object (incl pinned ones).
Definition: nodes.h:483
unsigned ntstorage() const
Return total number of doubles stored per value to record time history of each value (one for steady ...
Definition: nodes.cc:879
//////////////////////////////////////////////////////////////////////////// ////////////////////////...
Definition: matrices.h:386
FaceGeometry()
Constructor must call the constructor of the underlying solid element.
FaceGeometry()
Constructor must call the constructor of the underlying element.
FaceGeometry()
Constructor must call the constructor of the underlying element.
FaceGeometry()
Constructor must call the constructor of the underlying element.
FaceGeometry()
Constructor must call the constructor of the underlying element.
FaceGeometry()
Constructor must call the constructor of the underlying element.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
Definition: elements.h:5002
A general Finite Element class.
Definition: elements.h:1317
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2179
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...
Definition: elements.h:2597
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...
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.
Definition: elements.h:1436
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
Definition: elements.h:2615
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2214
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
Definition: elements.h:1763
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...
Definition: elements.cc:3328
unsigned long eqn_number(const unsigned &ieqn_local) const
Return the global equation number corresponding to the ieqn_local-th local equation number.
Definition: elements.h:708
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
Definition: elements.h:227
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
Definition: nodes.h:1054
TimeStepper *& position_time_stepper_pt()
Return a pointer to the position timestepper.
Definition: nodes.h:1022
An OomphLibError object which should be thrown when an run-time error is encountered....
Wrapper class for projectable elements. Adds "projectability" to the underlying ELEMENT.
Definition: projection.h:183
///////////////////////////////////////////////////////////////// ///////////////////////////////////...
unsigned nfields_for_projection()
Number of fields to be projected: 2*dim, corresponding to real and imag parts of the displacement com...
unsigned nhistory_values_for_projection(const unsigned &fld)
Number of history values to be stored for fld-th field. (includes present value!)
unsigned nhistory_values_for_coordinate_projection()
Number of positional history values: Read out from positional timestepper (Note: count includes curre...
ProjectableTimeHarmonicLinearElasticityElement()
Constructor [this was only required explicitly from gcc 4.5.2 onwards...].
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Return interpolated field fld at local coordinate s, at time level t (t=0: present; t>0: history valu...
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld.
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Specify the values associated with field fld. The information is returned in a vector of pairs which ...
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Return Jacobian of mapping and shape functions of field fld at local coordinate s.
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
Definition: Qelements.h:459
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
A base class that represents the fourth-rank elasticity tensor defined such that.
A base class for elements that solve the equations of time-harmonic linear elasticity in Cartesian co...
TimeHarmonicElasticityTensor *& elasticity_tensor_pt()
Return the pointer to the elasticity_tensor.
virtual void get_stress(const Vector< double > &s, DenseMatrix< std::complex< double >> &sigma) const =0
Return the Cauchy stress tensor, as calculated from the elasticity tensor at specified local coordina...
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contain...
TimeHarmonicElasticityTensor * Elasticity_tensor_pt
Pointer to the elasticity tensor.
BodyForceFctPt & body_force_fct_pt()
Access function: Pointer to body force function.
static double Default_omega_sq_value
Static default value for square of frequency.
unsigned ndof_types() const
The number of "DOF types" that degrees of freedom in this element are sub-divided into: for now lump ...
virtual std::complex< unsigned > u_index_time_harmonic_linear_elasticity(const unsigned i) const
Return the index at which the i-th real or imag unknown displacement component is stored....
double *& omega_sq_pt()
Access function for square of non-dim frequency.
TimeHarmonicLinearElasticityEquationsBase()
Constructor: Set null pointers for constitutive law and for isotropic growth function....
const double & omega_sq() const
Access function for square of non-dim frequency.
std::complex< double > interpolated_u_time_harmonic_linear_elasticity(const Vector< double > &s, const unsigned &i) const
Return FE interpolated displacement u[i] at local coordinate s.
double E(const unsigned &i, const unsigned &j, const unsigned &k, const unsigned &l) const
Access function to the entries in the elasticity tensor.
BodyForceFctPt body_force_fct_pt() const
Access function: Pointer to body force function (const version)
void body_force(const Vector< double > &x, Vector< std::complex< double >> &b) const
Evaluate body force at Eulerian coordinate x at present time (returns zero vector if no body force fu...
void get_strain(const Vector< double > &s, DenseMatrix< std::complex< double >> &strain) const
Return the strain tensor.
void interpolated_u_time_harmonic_linear_elasticity(const Vector< double > &s, Vector< std::complex< double >> &disp) const
Compute vector of FE interpolated displacement u at local coordinate s.
void(* BodyForceFctPt)(const double &t, const Vector< double > &x, Vector< std::complex< double >> &b)
Function pointer to function that specifies the body force as a function of the Cartesian coordinates...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
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 compute_norm(double &norm)
Compute norm of solution: square of the L2 norm.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
The jacobian is calculated by finite differences by default, We need only to take finite differences ...
unsigned required_nvalue(const unsigned &n) const
Number of values required at node n.
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(FILE *file_pt)
C-style output: 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].
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 fill_in_contribution_to_residuals(Vector< double > &residuals)
Return the residuals for the solid equations (the discretised principle of virtual displacements)
unsigned ntstorage() const
Return the number of doubles required to represent history (one for steady)
Definition: timesteppers.h:601
Time *const & time_pt() const
Access function for the pointer to time (const version)
Definition: timesteppers.h:572
double & time()
Return the current value of the continuous time.
Definition: timesteppers.h:123
//////////////////////////////////////////////////////////////////// ////////////////////////////////...