Tpml_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 Tri/Tet linear elasticity elements
27 #ifndef OOMPH_TPML_TIME_HARMONIC_LINEAR_ELASTICITY_ELEMENTS_HEADER
28 #define OOMPH_TPML_TIME_HARMONIC_LINEAR_ELASTICITY_ELEMENTS_HEADER
29 
30 
31 // Config header generated by autoconfig
32 #ifdef HAVE_CONFIG_H
33 #include <oomph-lib-config.h>
34 #endif
35 
36 
37 // OOMPH-LIB headers
38 #include "../generic/nodes.h"
39 #include "../generic/oomph_utilities.h"
40 #include "../generic/Telements.h"
41 #include "../generic/error_estimator.h"
42 
44 
45 namespace oomph
46 {
47  /// //////////////////////////////////////////////////////////////////////
48  /// //////////////////////////////////////////////////////////////////////
49  // TPMLTimeHarmonicLinearElasticityElement
50  /// /////////////////////////////////////////////////////////////////////
51  /// /////////////////////////////////////////////////////////////////////
52 
53 
54  //======================================================================
55  /// TPMLTimeHarmonicLinearElasticityElement<DIM,NNODE_1D> elements are
56  /// isoparametric triangular
57  /// DIM-dimensional PMLTimeHarmonicLinearElasticity elements with
58  /// NNODE_1D nodal points along each
59  /// element edge. Inherits from TElement and
60  /// PMLTimeHarmonicLinearElasticityEquations
61  //======================================================================
62  template<unsigned DIM, unsigned NNODE_1D>
64  : public TElement<DIM, NNODE_1D>,
66  public virtual ElementWithZ2ErrorEstimator
67  {
68  public:
69  /// Constructor: Call constructors for TElement and
70  /// PMLTimeHarmonicLinearElasticity equations
72  : TElement<DIM, NNODE_1D>(),
74  {
75  }
76 
77 
78  /// Broken copy constructor
81  delete;
82 
83  /// Broken assignment operator
84  void operator=(
86 
87  /// Output function:
88  void output(std::ostream& outfile)
89  {
91  }
92 
93  /// Output function:
94  void output(std::ostream& outfile, const unsigned& nplot)
95  {
97  }
98 
99 
100  /// C-style output function:
101  void output(FILE* file_pt)
102  {
104  }
105 
106  /// C-style output function:
107  void output(FILE* file_pt, const unsigned& n_plot)
108  {
110  }
111 
112 
113  /// Number of vertex nodes in the element
114  unsigned nvertex_node() const
115  {
117  }
118 
119  /// Pointer to the j-th vertex node in the element
120  Node* vertex_node_pt(const unsigned& j) const
121  {
123  }
124 
125  /// Order of recovery shape functions for Z2 error estimation:
126  /// Same order as shape functions.
127  unsigned nrecovery_order()
128  {
129  return NNODE_1D - 1;
130  }
131 
132  /// Number of 'flux' terms for Z2 error estimation
133  unsigned num_Z2_flux_terms()
134  {
135  // DIM Diagonal strain rates and DIM*(DIM-1)/2 off diagonal terms
136  return 2 * (DIM + DIM * (DIM - 1) / 2);
137  }
138 
139  /// Get 'flux' for Z2 error recovery: Upper triangular entries
140  /// in strain tensor.
142  {
143 #ifdef PARANOID
144  unsigned num_entries = 2 * (DIM + ((DIM * DIM) - DIM) / 2);
145  if (flux.size() != num_entries)
146  {
147  std::ostringstream error_message;
148  error_message << "The flux vector has the wrong number of entries, "
149  << flux.size() << ", whereas it should be " << num_entries
150  << std::endl;
151  throw OomphLibError(error_message.str(),
152  OOMPH_CURRENT_FUNCTION,
153  OOMPH_EXCEPTION_LOCATION);
154  }
155 #endif
156 
157  // Get strain matrix
159  this->get_strain(s, strain);
160 
161  // Pack into flux Vector
162  unsigned icount = 0;
163 
164  // Start with diagonal terms
165  for (unsigned i = 0; i < DIM; i++)
166  {
167  flux[icount] = strain(i, i).real();
168  icount++;
169  flux[icount] = strain(i, i).imag();
170  icount++;
171  }
172 
173  // Off diagonals row by row
174  for (unsigned i = 0; i < DIM; i++)
175  {
176  for (unsigned j = i + 1; j < DIM; j++)
177  {
178  flux[icount] = strain(i, j).real();
179  icount++;
180  flux[icount] = strain(i, j).imag();
181  icount++;
182  }
183  }
184  }
185  };
186 
187  //=======================================================================
188  /// Face geometry for the TPMLTimeHarmonicLinearElasticityElement elements:
189  /// The spatial
190  /// dimension of the face elements is one lower than that of the
191  /// bulk element but they have the same number of points
192  /// along their 1D edges.
193  //=======================================================================
194  template<unsigned DIM, unsigned NNODE_1D>
196  : public virtual TElement<DIM - 1, NNODE_1D>
197  {
198  public:
199  /// Constructor: Call the constructor for the
200  /// appropriate lower-dimensional QElement
201  FaceGeometry() : TElement<DIM - 1, NNODE_1D>() {}
202  };
203 
204  //=======================================================================
205  /// Face geometry for the 1D TPMLTimeHarmonicLinearElasticityElement
206  /// elements: Point elements
207  //=======================================================================
208  template<unsigned NNODE_1D>
210  : public virtual PointElement
211  {
212  public:
213  /// Constructor: Call the constructor for the
214  /// appropriate lower-dimensional TElement
216  };
217 
218 
219  /// ///////////////////////////////////////////////////////////////////
220  /// ///////////////////////////////////////////////////////////////////
221  /// ///////////////////////////////////////////////////////////////////
222 
223 
224  //=======================================================================
225  /// Policy class defining the elements to be used in the actual
226  /// PML layers. Same spatial dimension and nnode_1d but quads
227  /// rather than triangles.
228  //=======================================================================
229  template<unsigned NNODE_1D>
231  : public virtual QPMLTimeHarmonicLinearElasticityElement<2, NNODE_1D>
232  {
233  public:
234  /// Constructor: Call the constructor for the
235  /// appropriate QElement
237  {
238  }
239  };
240 
241  //=======================================================================
242  /// Face geometry for the TPMLTimeHarmonicLinearElasticityElement elements:
243  /// The spatial
244  /// dimension of the face elements is one lower than that of the
245  /// bulk element but they have the same number of points
246  /// along their 1D edges.
247  //=======================================================================
248  template<unsigned DIM, unsigned NNODE_1D>
251  : public virtual QElement<DIM - 1, NNODE_1D>
252  {
253  public:
254  /// Constructor: Call the constructor for the
255  /// appropriate lower-dimensional QElement
256  FaceGeometry() : QElement<DIM - 1, NNODE_1D>() {}
257  };
258 
259 
260  /// /////////////////////////////////////////////////////////////////////
261  /// /////////////////////////////////////////////////////////////////////
262  /// /////////////////////////////////////////////////////////////////////
263 
264 
265  //=======================================================================
266  /// Policy class defining the elements to be used in the actual
267  /// PML layers. Same spatial dimension and nnode_1d but quads
268  /// rather than triangles.
269  //=======================================================================
270  template<unsigned NNODE_1D>
273  : public virtual QPMLTimeHarmonicLinearElasticityElement<2, NNODE_1D>
274  {
275  public:
276  /// Constructor: Call the constructor for the
277  /// appropriate QElement
279  {
280  }
281  };
282 
283 
284 } // namespace oomph
285 
286 #endif
static char t char * s
Definition: cfortran.h:568
cstr elem_len * i
Definition: cfortran.h:603
//////////////////////////////////////////////////////////////////////////// ////////////////////////...
Definition: matrices.h:386
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional TElement.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
Definition: elements.h:5002
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
Definition: nodes.h:906
An OomphLibError object which should be thrown when an run-time error is encountered....
General definition of policy class defining the elements to be used in the actual PML layers....
Definition: pml_meshes.h:48
void get_strain(const Vector< double > &s, DenseMatrix< std::complex< double >> &strain) const
Return the strain tensor.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
void output(std::ostream &outfile)
Output: x,y,[z],u_r,v_r,[w_r],u_i,v_i,[w_i].
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
Definition: elements.h:3443
///////////////////////////////////////////////////////////////// ///////////////////////////////////...
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
Definition: Qelements.h:459
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
General TElement class.
Definition: Telements.h:1208
//////////////////////////////////////////////////////////////////////
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Upper triangular entries in strain tensor.
TPMLTimeHarmonicLinearElasticityElement(const TPMLTimeHarmonicLinearElasticityElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void output(std::ostream &outfile, const unsigned &nplot)
Output function:
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function:
void operator=(const TPMLTimeHarmonicLinearElasticityElement< DIM, NNODE_1D > &)=delete
Broken assignment operator.
unsigned nvertex_node() const
Number of vertex nodes in the element.
TPMLTimeHarmonicLinearElasticityElement()
Constructor: Call constructors for TElement and PMLTimeHarmonicLinearElasticity equations.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...