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// Version 1.0; svn revision $LastChangedRevision$
7 //LIC//
8 //LIC// $LastChangedDate$
9 //LIC//
10 //LIC// Copyright (C) 2006-2016 Matthias Heil and Andrew Hazel
11 //LIC//
12 //LIC// This library is free software; you can redistribute it and/or
13 //LIC// modify it under the terms of the GNU Lesser General Public
14 //LIC// License as published by the Free Software Foundation; either
15 //LIC// version 2.1 of the License, or (at your option) any later version.
16 //LIC//
17 //LIC// This library is distributed in the hope that it will be useful,
18 //LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
19 //LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 //LIC// Lesser General Public License for more details.
21 //LIC//
22 //LIC// You should have received a copy of the GNU Lesser General Public
23 //LIC// License along with this library; if not, write to the Free Software
24 //LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
25 //LIC// 02110-1301 USA.
26 //LIC//
27 //LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
28 //LIC//
29 //LIC//====================================================================
30 //Header file for Tri/Tet linear elasticity elements
31 #ifndef OOMPH_TPML_TIME_HARMONIC_LINEAR_ELASTICITY_ELEMENTS_HEADER
32 #define OOMPH_TPML_TIME_HARMONIC_LINEAR_ELASTICITY_ELEMENTS_HEADER
33 
34 
35 // Config header generated by autoconfig
36 #ifdef HAVE_CONFIG_H
37  #include <oomph-lib-config.h>
38 #endif
39 
40 
41 //OOMPH-LIB headers
42 #include "../generic/nodes.h"
43 #include "../generic/oomph_utilities.h"
44 #include "../generic/Telements.h"
45 #include "../generic/error_estimator.h"
46 
48 
49 namespace oomph
50 {
51 
52 /////////////////////////////////////////////////////////////////////////
53 /////////////////////////////////////////////////////////////////////////
54 // TPMLTimeHarmonicLinearElasticityElement
55 ////////////////////////////////////////////////////////////////////////
56 ////////////////////////////////////////////////////////////////////////
57 
58 
59 
60 //======================================================================
61 /// TPMLTimeHarmonicLinearElasticityElement<DIM,NNODE_1D> elements are
62 /// isoparametric triangular
63 /// DIM-dimensional PMLTimeHarmonicLinearElasticity elements with
64 /// NNODE_1D nodal points along each
65 /// element edge. Inherits from TElement and
66 /// PMLTimeHarmonicLinearElasticityEquations
67 //======================================================================
68 template <unsigned DIM, unsigned NNODE_1D, class PML_ELEMENT>
69  class TPMLTimeHarmonicLinearElasticityElement : public TElement<DIM,NNODE_1D>,
70  public PMLTimeHarmonicLinearElasticityEquations<DIM,PML_ELEMENT>,
71  public virtual ElementWithZ2ErrorEstimator
72  {
73 
74  public:
75 
76  ///\short Constructor: Call constructors for TElement and
77  /// PMLTimeHarmonicLinearElasticity equations
79  PMLTimeHarmonicLinearElasticityEquations<DIM, PML_ELEMENT>() { }
80 
81 
82  /// Broken copy constructor
84  {
85  BrokenCopy::broken_copy("TPMLTimeHarmonicLinearElasticityElement");
86  }
87 
88  /// Broken assignment operator
90  {
91  BrokenCopy::broken_assign("TPMLTimeHarmonicLinearElasticityElement");
92  }
93 
94  /// \short Output function:
95  void output(std::ostream &outfile)
96  {
98  }
99 
100  /// \short Output function:
101  void output(std::ostream &outfile, const unsigned &nplot)
102  {
104  }
105 
106 
107  /// \short C-style output function:
108  void output(FILE* file_pt)
109  {
111  }
112 
113  /// \short C-style output function:
114  void output(FILE* file_pt, const unsigned &n_plot)
115  {
117  }
118 
119 
120  /// \short Number of vertex nodes in the element
121  unsigned nvertex_node() const
123 
124  /// \short Pointer to the j-th vertex node in the element
125  Node* vertex_node_pt(const unsigned& j) const
126  {
128  }
129 
130  /// \short Order of recovery shape functions for Z2 error estimation:
131  /// Same order as shape functions.
132  unsigned nrecovery_order() {return NNODE_1D-1;}
133 
134  /// Number of 'flux' terms for Z2 error estimation
135  unsigned num_Z2_flux_terms()
136  {
137  // DIM Diagonal strain rates and DIM*(DIM-1)/2 off diagonal terms
138  return 2*(DIM + DIM*(DIM-1)/2);
139  }
140 
141  /// \short Get 'flux' for Z2 error recovery: Upper triangular entries
142  /// in strain tensor.
144  {
145 #ifdef PARANOID
146  unsigned num_entries=2*(DIM+((DIM*DIM)-DIM)/2);
147  if (flux.size()!=num_entries)
148  {
149  std::ostringstream error_message;
150  error_message << "The flux vector has the wrong number of entries, "
151  << flux.size() << ", whereas it should be "
152  << num_entries << std::endl;
153  throw OomphLibError(
154  error_message.str(),
155  OOMPH_CURRENT_FUNCTION,
156  OOMPH_EXCEPTION_LOCATION);
157  }
158 #endif
159 
160  // Get strain matrix
161  DenseMatrix<std::complex<double> > strain(DIM);
162  this->get_strain(s,strain);
163 
164  // Pack into flux Vector
165  unsigned icount=0;
166 
167  // Start with diagonal terms
168  for(unsigned i=0;i<DIM;i++)
169  {
170  flux[icount]=strain(i,i).real();
171  icount++;
172  flux[icount]=strain(i,i).imag();
173  icount++;
174  }
175 
176  //Off diagonals row by row
177  for(unsigned i=0;i<DIM;i++)
178  {
179  for(unsigned j=i+1;j<DIM;j++)
180  {
181  flux[icount]=strain(i,j).real();
182  icount++;
183  flux[icount]=strain(i,j).imag();
184  icount++;
185  }
186  }
187  }
188 
189 
190 
191 };
192 
193 //=======================================================================
194 /// Face geometry for the TPMLTimeHarmonicLinearElasticityElement elements:
195 /// The spatial
196 /// dimension of the face elements is one lower than that of the
197 /// bulk element but they have the same number of points
198 /// along their 1D edges.
199 //=======================================================================
200 template<unsigned DIM, unsigned NNODE_1D, class PML_ELEMENT>
201 class FaceGeometry<TPMLTimeHarmonicLinearElasticityElement<DIM,NNODE_1D,PML_ELEMENT> >:
202  public virtual TElement<DIM-1,NNODE_1D>
203 {
204 
205  public:
206 
207  /// \short Constructor: Call the constructor for the
208  /// appropriate lower-dimensional QElement
209  FaceGeometry() : TElement<DIM-1,NNODE_1D>() {}
210 
211 };
212 
213 //=======================================================================
214 /// Face geometry for the 1D TPMLTimeHarmonicLinearElasticityElement
215 /// elements: Point elements
216 //=======================================================================
217 template<unsigned NNODE_1D, class PML_ELEMENT>
218 class FaceGeometry<TPMLTimeHarmonicLinearElasticityElement<1,NNODE_1D,PML_ELEMENT> >:
219  public virtual PointElement
220 {
221 
222  public:
223 
224  /// \short Constructor: Call the constructor for the
225  /// appropriate lower-dimensional TElement
227 
228 };
229 
230 
231 
232 //////////////////////////////////////////////////////////////////////
233 //////////////////////////////////////////////////////////////////////
234 //////////////////////////////////////////////////////////////////////
235 
236 
237 //=======================================================================
238 /// Policy class defining the elements to be used in the actual
239 /// PML layers. Same spatial dimension and nnode_1d but quads
240 /// rather than triangles.
241 //=======================================================================
242 template<unsigned NNODE_1D, class PML_ELEMENT>
244  TPMLTimeHarmonicLinearElasticityElement<2,NNODE_1D,PML_ELEMENT> >:
245  public virtual QPMLTimeHarmonicLinearElasticityElement<2,NNODE_1D,PML_ELEMENT>
246 {
247 
248  public:
249 
250  /// \short Constructor: Call the constructor for the
251  /// appropriate QElement
253  QPMLTimeHarmonicLinearElasticityElement<2,NNODE_1D,PML_ELEMENT>()
254  {}
255 
256 };
257 
258 //=======================================================================
259 /// Face geometry for the TPMLTimeHarmonicLinearElasticityElement elements:
260 /// The spatial
261 /// dimension of the face elements is one lower than that of the
262 /// bulk element but they have the same number of points
263 /// along their 1D edges.
264 //=======================================================================
265 template<unsigned DIM, unsigned NNODE_1D, class PML_ELEMENT>
267  TPMLTimeHarmonicLinearElasticityElement<DIM,NNODE_1D,PML_ELEMENT> > >:
268  public virtual QElement<DIM-1,NNODE_1D>
269 {
270 
271  public:
272 
273  /// \short Constructor: Call the constructor for the
274  /// appropriate lower-dimensional QElement
275  FaceGeometry() : QElement<DIM-1,NNODE_1D>() {}
276 
277 };
278 
279 
280 ////////////////////////////////////////////////////////////////////////
281 ////////////////////////////////////////////////////////////////////////
282 ////////////////////////////////////////////////////////////////////////
283 
284 
285 //=======================================================================
286 /// Policy class defining the elements to be used in the actual
287 /// PML layers. Same spatial dimension and nnode_1d but quads
288 /// rather than triangles.
289 //=======================================================================
290 template<unsigned NNODE_1D, class PML_ELEMENT>
293  TPMLTimeHarmonicLinearElasticityElement<2,NNODE_1D,PML_ELEMENT> > >:
294  public virtual QPMLTimeHarmonicLinearElasticityElement<2,NNODE_1D,PML_ELEMENT>
295 {
296 
297  public:
298 
299  /// \short Constructor: Call the constructor for the
300  /// appropriate QElement
302  QPMLTimeHarmonicLinearElasticityElement<2,NNODE_1D,PML_ELEMENT>()
303  {}
304 
305 };
306 
307 
308 
309 
310 }
311 
312 #endif
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional TElement. ...
void broken_copy(const std::string &class_name)
Issue error message and terminate execution.
void operator=(const TPMLTimeHarmonicLinearElasticityElement< DIM, NNODE_1D, PML_ELEMENT > &)
Broken assignment operator.
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
unsigned num_Z2_flux_terms()
Number of &#39;flux&#39; terms for Z2 error estimation.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
cstr elem_len * i
Definition: cfortran.h:607
void output(std::ostream &outfile, const unsigned &nplot)
Output function:
TPMLTimeHarmonicLinearElasticityElement()
Constructor: Call constructors for TElement and PMLTimeHarmonicLinearElasticity equations.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
Definition: nodes.h:852
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get &#39;flux&#39; for Z2 error recovery: Upper triangular entries in strain tensor.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
void output(std::ostream &outfile)
Output: x,y,[z],u_r,v_r,[w_r],u_i,v_i,[w_i].
unsigned nvertex_node() const
Number of vertex nodes in the element.
TPMLTimeHarmonicLinearElasticityElement(const TPMLTimeHarmonicLinearElasticityElement< DIM, NNODE_1D, PML_ELEMENT > &dummy)
Broken copy constructor.
Time-harmonic linear elasticity upgraded to become projectable.
static char t char * s
Definition: cfortran.h:572
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement. ...
void get_strain(const Vector< double > &s, DenseMatrix< std::complex< double > > &strain) const
Return the strain tensor.
void broken_assign(const std::string &class_name)
Issue error message and terminate execution.
Class for dense matrices, storing all the values of the matrix as a pointer to a pointer with assorte...
Definition: communicator.h:50
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function: