Ttime_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_TTIME_HARMONIC_LINEAR_ELASTICITY_ELEMENTS_HEADER
28 #define OOMPH_TTIME_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"
42 #include "../generic/error_estimator.h"
43 
44 namespace oomph
45 {
46  /// //////////////////////////////////////////////////////////////////////
47  /// //////////////////////////////////////////////////////////////////////
48  // TTimeHarmonicLinearElasticityElement
49  /// /////////////////////////////////////////////////////////////////////
50  /// /////////////////////////////////////////////////////////////////////
51 
52 
53  //======================================================================
54  /// TTimeHarmonicLinearElasticityElement<DIM,NNODE_1D> elements are
55  /// isoparametric triangular
56  /// DIM-dimensional TimeHarmonicLinearElasticity elements with
57  /// NNODE_1D nodal points along each
58  /// element edge. Inherits from TElement and
59  /// TimeHarmonicLinearElasticityEquations
60  //======================================================================
61  template<unsigned DIM, unsigned NNODE_1D>
63  : public TElement<DIM, NNODE_1D>,
65  public virtual ElementWithZ2ErrorEstimator
66  {
67  public:
68  /// Constructor: Call constructors for TElement and
69  /// TimeHarmonicLinearElasticity equations
71  : TElement<DIM, NNODE_1D>(), TimeHarmonicLinearElasticityEquations<DIM>()
72  {
73  }
74 
75 
76  /// Broken copy constructor
79  delete;
80 
81  /// Broken assignment operator
82  // Commented out broken assignment operator because this can lead to a
83  // conflict warning when used in the virtual inheritence hierarchy.
84  // Essentially the compiler doesn't realise that two separate
85  // implementations of the broken function are the same and so, quite
86  // rightly, it shouts.
87  /*void operator=(const TTimeHarmonicLinearElasticityElement<DIM,NNODE_1D>&)
88  * = delete;*/
89 
90  /// Output function:
91  void output(std::ostream& outfile)
92  {
94  }
95 
96  /// Output function:
97  void output(std::ostream& outfile, const unsigned& nplot)
98  {
100  }
101 
102 
103  /// C-style output function:
104  void output(FILE* file_pt)
105  {
107  }
108 
109  /// C-style output function:
110  void output(FILE* file_pt, const unsigned& n_plot)
111  {
113  }
114 
115 
116  /// Number of vertex nodes in the element
117  unsigned nvertex_node() const
118  {
120  }
121 
122  /// Pointer to the j-th vertex node in the element
123  Node* vertex_node_pt(const unsigned& j) const
124  {
126  }
127 
128  /// Order of recovery shape functions for Z2 error estimation:
129  /// Same order as shape functions.
130  unsigned nrecovery_order()
131  {
132  return NNODE_1D - 1;
133  }
134 
135  /// Number of 'flux' terms for Z2 error estimation
136  unsigned num_Z2_flux_terms()
137  {
138  // DIM Diagonal strain rates and DIM*(DIM-1)/2 off diagonal terms
139  return 2 * (DIM + DIM * (DIM - 1) / 2);
140  }
141 
142  /// Get 'flux' for Z2 error recovery: Upper triangular entries
143  /// in strain tensor.
145  {
146 #ifdef PARANOID
147  unsigned num_entries = 2 * (DIM + ((DIM * DIM) - DIM) / 2);
148  if (flux.size() != num_entries)
149  {
150  std::ostringstream error_message;
151  error_message << "The flux vector has the wrong number of entries, "
152  << flux.size() << ", whereas it should be " << num_entries
153  << std::endl;
154  throw OomphLibError(error_message.str(),
155  OOMPH_CURRENT_FUNCTION,
156  OOMPH_EXCEPTION_LOCATION);
157  }
158 #endif
159 
160  // Get strain matrix
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  /// Face geometry for the TTimeHarmonicLinearElasticityElement elements:
192  /// The spatial
193  /// dimension of the face elements is one lower than that of the
194  /// bulk element but they have the same number of points
195  /// along their 1D edges.
196  //=======================================================================
197  template<unsigned DIM, unsigned NNODE_1D>
199  : public virtual TElement<DIM - 1, NNODE_1D>
200  {
201  public:
202  /// Constructor: Call the constructor for the
203  /// appropriate lower-dimensional QElement
204  FaceGeometry() : TElement<DIM - 1, NNODE_1D>() {}
205  };
206 
207  //=======================================================================
208  /// Face geometry for the 1D TTimeHarmonicLinearElasticityElement elements:
209  /// Point elements
210  //=======================================================================
211  template<unsigned NNODE_1D>
213  : public virtual PointElement
214  {
215  public:
216  /// Constructor: Call the constructor for the
217  /// appropriate lower-dimensional TElement
219  };
220 
221 
222 } // namespace oomph
223 
224 #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....
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
Definition: elements.h:3443
General TElement class.
Definition: Telements.h:1208
//////////////////////////////////////////////////////////////////////
TTimeHarmonicLinearElasticityElement(const TTimeHarmonicLinearElasticityElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
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.
TTimeHarmonicLinearElasticityElement()
Constructor: Call constructors for TElement and TimeHarmonicLinearElasticity equations.
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
void output(std::ostream &outfile)
Broken assignment operator.
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Upper triangular entries in strain tensor.
void output(std::ostream &outfile, const unsigned &nplot)
Output function:
unsigned nvertex_node() const
Number of vertex nodes in the element.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function:
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].
//////////////////////////////////////////////////////////////////// ////////////////////////////////...