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-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 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 {
48  // TTimeHarmonicLinearElasticityElement
51 
52 
53  //======================================================================
60  //======================================================================
61  template<unsigned DIM, unsigned NNODE_1D>
63  : public TElement<DIM, NNODE_1D>,
65  public virtual ElementWithZ2ErrorEstimator
66  {
67  public:
72  {
73  }
74 
75 
79  delete;
80 
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 
91  void output(std::ostream& outfile)
92  {
94  }
95 
97  void output(std::ostream& outfile, const unsigned& nplot)
98  {
100  }
101 
102 
104  void output(FILE* file_pt)
105  {
107  }
108 
110  void output(FILE* file_pt, const unsigned& n_plot)
111  {
113  }
114 
115 
117  unsigned nvertex_node() const
118  {
120  }
121 
123  Node* vertex_node_pt(const unsigned& j) const
124  {
126  }
127 
130  unsigned nrecovery_order()
131  {
132  return NNODE_1D - 1;
133  }
134 
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 
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(),
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  //=======================================================================
196  //=======================================================================
197  template<unsigned DIM, unsigned NNODE_1D>
199  : public virtual TElement<DIM - 1, NNODE_1D>
200  {
201  public:
204  FaceGeometry() : TElement<DIM - 1, NNODE_1D>() {}
205  };
206 
207  //=======================================================================
210  //=======================================================================
211  template<unsigned NNODE_1D>
213  : public virtual PointElement
214  {
215  public:
219  };
220 
221 
222 } // namespace oomph
223 
224 #endif
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Definition: matrices.h:386
Definition: error_estimator.h:79
FaceGeometry()
Definition: Ttime_harmonic_linear_elasticity_elements.h:218
FaceGeometry()
Definition: Ttime_harmonic_linear_elasticity_elements.h:204
Definition: elements.h:4998
Definition: nodes.h:906
Definition: oomph_definitions.h:222
Definition: elements.h:3439
Definition: Telements.h:1208
Definition: Ttime_harmonic_linear_elasticity_elements.h:66
TTimeHarmonicLinearElasticityElement(const TTimeHarmonicLinearElasticityElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
unsigned nrecovery_order()
Definition: Ttime_harmonic_linear_elasticity_elements.h:130
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: Ttime_harmonic_linear_elasticity_elements.h:123
TTimeHarmonicLinearElasticityElement()
Definition: Ttime_harmonic_linear_elasticity_elements.h:70
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
Definition: Ttime_harmonic_linear_elasticity_elements.h:136
void output(std::ostream &outfile)
Broken assignment operator.
Definition: Ttime_harmonic_linear_elasticity_elements.h:91
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Definition: Ttime_harmonic_linear_elasticity_elements.h:144
void output(std::ostream &outfile, const unsigned &nplot)
Output function:
Definition: Ttime_harmonic_linear_elasticity_elements.h:97
void output(FILE *file_pt)
C-style output function:
Definition: Ttime_harmonic_linear_elasticity_elements.h:104
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: Ttime_harmonic_linear_elasticity_elements.h:117
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function:
Definition: Ttime_harmonic_linear_elasticity_elements.h:110
void get_strain(const Vector< double > &s, DenseMatrix< std::complex< double >> &strain) const
Return the strain tensor.
Definition: time_harmonic_linear_elasticity_elements.cc:49
Definition: time_harmonic_linear_elasticity_elements.h:316
void output(std::ostream &outfile)
Output: x,y,[z],u_r,v_r,[w_r],u_i,v_i,[w_i].
Definition: time_harmonic_linear_elasticity_elements.h:359
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
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
#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