refineable_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-2022 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 refineable linear elasticity elements
27 
28 // Include guard to prevent multiple inclusions of this header
29 #ifndef OOMPH_REFINEABLE_LINEAR_ELASTICITY_ELEMENTS_HEADER
30 #define OOMPH_REFINEABLE_LINEAR_ELASTICITY_ELEMENTS_HEADER
31 
32 // oomph-lib headers
34 #include "../generic/refineable_quad_element.h"
35 #include "../generic/refineable_brick_element.h"
36 #include "../generic/error_estimator.h"
37 
38 namespace oomph
39 {
40  //========================================================================
42  //========================================================================
43  template<unsigned DIM>
45  : public virtual LinearElasticityEquations<DIM>,
46  public virtual RefineableElement,
47  public virtual ElementWithZ2ErrorEstimator
48  {
49  public:
55  {
56  }
57 
58 
63  void get_interpolated_values(const unsigned& t,
64  const Vector<double>& s,
65  Vector<double>& values)
66  {
67  // Create enough initialised storage
68  values.resize(DIM, 0.0);
69 
70  // Find out how many nodes there are
71  unsigned n_node = this->nnode();
72 
73  // Shape functions
74  Shape psi(n_node);
75  this->shape(s, psi);
76 
77  // Calculate displacements
78  for (unsigned i = 0; i < DIM; i++)
79  {
80  // Get the index at which the i-th velocity is stored
81  unsigned u_nodal_index = this->u_index_linear_elasticity(i);
82  for (unsigned l = 0; l < n_node; l++)
83  {
84  values[i] += this->nodal_value(t, l, u_nodal_index) * psi(l);
85  }
86  }
87  }
88 
94  Vector<double>& values)
95  {
96  values.resize(DIM);
97  this->interpolated_u_linear_elasticity(s, values);
98  }
99 
101  unsigned num_Z2_flux_terms()
102  {
103  // DIM Diagonal strain rates and DIM*(DIM-1)/2 off diagonal terms
104  return DIM + DIM * (DIM - 1) / 2;
105  }
106 
110  {
111 #ifdef PARANOID
112  unsigned num_entries = DIM + ((DIM * DIM) - DIM) / 2;
113  if (flux.size() != num_entries)
114  {
115  std::ostringstream error_message;
116  error_message << "The flux vector has the wrong number of entries, "
117  << flux.size() << ", whereas it should be " << num_entries
118  << std::endl;
119  throw OomphLibError(error_message.str(),
122  }
123 #endif
124 
125  // Get strain matrix
126  DenseMatrix<double> strain(DIM);
127  this->get_strain(s, strain);
128 
129  // Pack into flux Vector
130  unsigned icount = 0;
131 
132  // Start with diagonal terms
133  for (unsigned i = 0; i < DIM; i++)
134  {
135  flux[icount] = strain(i, i);
136  icount++;
137  }
138 
139  // Off diagonals row by row
140  for (unsigned i = 0; i < DIM; i++)
141  {
142  for (unsigned j = i + 1; j < DIM; j++)
143  {
144  flux[icount] = strain(i, j);
145  icount++;
146  }
147  }
148  }
149 
151  unsigned ncont_interpolated_values() const
152  {
153  return DIM;
154  }
155 
158  {
159  RefineableLinearElasticityEquations<DIM>* cast_father_element_pt =
161  this->father_element_pt());
162 
163  // Set pointer to body force function
164  this->Body_force_fct_pt = cast_father_element_pt->body_force_fct_pt();
165 
166  // Set pointer to the contitutive law
167  this->Elasticity_tensor_pt =
168  cast_father_element_pt->elasticity_tensor_pt();
169 
170  // Set the timescale ratio (non-dim. density)
171  this->Lambda_sq_pt = cast_father_element_pt->lambda_sq_pt();
172 
174  this->Unsteady = cast_father_element_pt->is_inertia_enabled();
175  }
176 
177 
178  private:
181  Vector<double>& residuals, DenseMatrix<double>& jacobian, unsigned flag);
182  };
183 
184 
188 
189  //========================================================================
191  //========================================================================
192  template<unsigned DIM, unsigned NNODE_1D>
194  : public virtual QLinearElasticityElement<DIM, NNODE_1D>,
195  public virtual RefineableLinearElasticityEquations<DIM>,
196  public virtual RefineableQElement<DIM>
197  {
198  public:
201  : QLinearElasticityElement<DIM, NNODE_1D>(),
205  {
206  }
207 
209  void rebuild_from_sons(Mesh*& mesh_pt) {}
210 
212  unsigned nvertex_node() const
213  {
215  }
216 
218  Node* vertex_node_pt(const unsigned& j) const
219  {
221  }
222 
225  unsigned nrecovery_order()
226  {
227  return NNODE_1D - 1;
228  }
229 
232  };
233 
234 
235  //======================================================================
237  //======================================================================
238  template<unsigned DIM>
240  : public QLinearElasticityElement<DIM, 2>,
241  public virtual RefineableLinearElasticityEquations<DIM>,
242  public virtual PRefineableQElement<DIM>
243  {
244  public:
247  : RefineableElement(),
251  {
252  // Set integration scheme
253  // (To avoid memory leaks in pre-build and p-refine where new
254  // integration schemes are created)
256  }
257 
260  {
261  delete this->integral_pt();
262  }
263 
264 
267  const PRefineableQLinearElasticityElement<DIM>& dummy) = delete;
268 
270  /* void operator=(const PRefineableQLinearElasticityElement<DIM>&) =
271  * delete;*/
272 
273  void further_build();
274 
276  unsigned ncont_interpolated_values() const
277  {
278  return 1;
279  }
280 
282  unsigned nvertex_node() const
283  {
285  }
286 
288  Node* vertex_node_pt(const unsigned& j) const
289  {
291  }
292 
295  // unsigned nrecovery_order()
296  // {
297  // if(this->nnode_1d() < 4) {return (this->nnode_1d()-1);}
298  // else {return 3;}
299  // }
302  unsigned nrecovery_order()
303  {
304  return 3;
305  }
306 
308  std::ostream& outfile,
310  double& error,
311  double& norm);
312  };
313 
314 
315  //==============================================================
317  //==============================================================
318  template<unsigned NNODE_1D>
320  : public virtual QElement<1, NNODE_1D>
321  {
322  public:
323  // Make sure that we call the constructor of the QElement
324  // Only the Intel compiler seems to need this!
325  FaceGeometry() : QElement<1, NNODE_1D>() {}
326  };
327 
328  //==============================================================
331  //==============================================================
332  template<unsigned NNODE_1D>
335  : public virtual PointElement
336  {
337  public:
338  // Make sure that we call the constructor of the QElement
339  // Only the Intel compiler seems to need this!
341  };
342 
343 
344  //==============================================================
346  //==============================================================
347  template<unsigned NNODE_1D>
349  : public virtual QElement<2, NNODE_1D>
350  {
351  public:
352  // Make sure that we call the constructor of the QElement
353  // Only the Intel compiler seems to need this!
354  FaceGeometry() : QElement<2, NNODE_1D>() {}
355  };
356 
357  //==============================================================
360  //==============================================================
361  template<unsigned NNODE_1D>
364  : public virtual QElement<1, NNODE_1D>
365  {
366  public:
367  // Make sure that we call the constructor of the QElement
368  // Only the Intel compiler seems to need this!
369  FaceGeometry() : QElement<1, NNODE_1D>() {}
370  };
371 
372 
373 } // namespace oomph
374 
375 #endif
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Definition: error_estimator.h:79
FaceGeometry()
Definition: refineable_linear_elasticity_elements.h:340
FaceGeometry()
Definition: refineable_linear_elasticity_elements.h:369
FaceGeometry()
Definition: refineable_linear_elasticity_elements.h:325
FaceGeometry()
Definition: refineable_linear_elasticity_elements.h:354
Definition: elements.h:4998
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
virtual unsigned nvertex_node() const
Definition: elements.h:2491
virtual void shape(const Vector< double > &s, Shape &psi) const =0
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Definition: elements.h:1759
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1963
virtual Node * vertex_node_pt(const unsigned &j) const
Definition: elements.h:2500
virtual void set_integration_scheme(Integral *const &integral_pt)
Set the spatial integration scheme.
Definition: elements.cc:3210
Definition: integral.h:1281
ElasticityTensor * Elasticity_tensor_pt
Pointer to the elasticity tensor.
Definition: linear_elasticity_elements.h:352
double *& lambda_sq_pt()
Access function for pointer to timescale ratio (nondim density)
Definition: linear_elasticity_elements.h:195
ElasticityTensor *& elasticity_tensor_pt()
Return the pointer to the elasticity_tensor.
Definition: linear_elasticity_elements.h:174
double * Lambda_sq_pt
Timescale ratio (non-dim. density)
Definition: linear_elasticity_elements.h:355
void get_strain(const Vector< double > &s, DenseMatrix< double > &strain) const
Return the strain tensor.
Definition: linear_elasticity_elements.cc:47
BodyForceFctPt Body_force_fct_pt
Pointer to body force function.
Definition: linear_elasticity_elements.h:361
bool is_inertia_enabled() const
Access function to flag that switches inertia on/off (const version)
Definition: linear_elasticity_elements.h:226
void interpolated_u_linear_elasticity(const Vector< double > &s, Vector< double > &disp) const
Compute vector of FE interpolated displacement u at local coordinate s.
Definition: linear_elasticity_elements.h:98
bool Unsteady
Flag that switches inertia on/off.
Definition: linear_elasticity_elements.h:358
BodyForceFctPt & body_force_fct_pt()
Access function: Pointer to body force function.
Definition: linear_elasticity_elements.h:201
virtual unsigned u_index_linear_elasticity(const unsigned i) const
Definition: linear_elasticity_elements.h:62
Definition: linear_elasticity_elements.h:380
Definition: mesh.h:67
Definition: nodes.h:906
Definition: oomph_definitions.h:222
Definition: Qelements.h:2274
p-refineable version of 2D QLinearElasticityElement elements
Definition: refineable_linear_elasticity_elements.h:243
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: refineable_linear_elasticity_elements.h:288
PRefineableQLinearElasticityElement()
Constructor, simply call the other constructors.
Definition: refineable_linear_elasticity_elements.h:246
PRefineableQLinearElasticityElement(const PRefineableQLinearElasticityElement< DIM > &dummy)=delete
Broken copy constructor.
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: refineable_linear_elasticity_elements.h:282
void further_build()
Broken assignment operator.
Definition: refineable_linear_elasticity_elements.cc:383
void compute_energy_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_grad_pt, double &error, double &norm)
Get error against and norm of exact solution.
Definition: refineable_linear_elasticity_elements.cc:297
~PRefineableQLinearElasticityElement()
Destructor (to avoid memory leaks)
Definition: refineable_linear_elasticity_elements.h:259
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 1.
Definition: refineable_linear_elasticity_elements.h:276
unsigned nrecovery_order()
Definition: refineable_linear_elasticity_elements.h:302
Definition: elements.h:3439
Definition: Qelements.h:459
Definition: linear_elasticity_elements.h:487
Definition: refineable_elements.h:97
virtual RefineableElement * father_element_pt() const
Return a pointer to the father element.
Definition: refineable_elements.h:539
Class for Refineable LinearElasticity equations.
Definition: refineable_linear_elasticity_elements.h:48
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: DIM.
Definition: refineable_linear_elasticity_elements.h:151
void further_build()
Further build function, pass the pointers down to the sons.
Definition: refineable_linear_elasticity_elements.h:157
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
Definition: refineable_linear_elasticity_elements.h:101
RefineableLinearElasticityEquations()
Constructor.
Definition: refineable_linear_elasticity_elements.h:51
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Definition: refineable_linear_elasticity_elements.h:109
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Definition: refineable_linear_elasticity_elements.h:63
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Definition: refineable_linear_elasticity_elements.h:93
void fill_in_generic_contribution_to_residuals_linear_elasticity(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Overloaded helper function to take hanging nodes into account.
Definition: refineable_linear_elasticity_elements.cc:39
Definition: Qelements.h:2259
Class for refineable QLinearElasticityElement elements.
Definition: refineable_linear_elasticity_elements.h:197
void further_setup_hanging_nodes()
No additional hanging node procedures are required.
Definition: refineable_linear_elasticity_elements.h:231
void rebuild_from_sons(Mesh *&mesh_pt)
Empty rebuild from sons, no need to reconstruct anything here.
Definition: refineable_linear_elasticity_elements.h:209
unsigned nrecovery_order()
Definition: refineable_linear_elasticity_elements.h:225
RefineableQLinearElasticityElement()
Constructor:
Definition: refineable_linear_elasticity_elements.h:200
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: refineable_linear_elasticity_elements.h:212
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: refineable_linear_elasticity_elements.h:218
Definition: shape.h:76
RealScalar s
Definition: level1_cplx_impl.h:130
#define DIM
Definition: linearised_navier_stokes_elements.h:44
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
int error
Definition: calibrate.py:297
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
t
Definition: plotPSD.py:36
#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