Tlinear_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_TLINEAR_ELASTICITY_ELEMENTS_HEADER
28 #define OOMPH_TLINEAR_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 
45 namespace oomph
46 {
49  // TLinearElasticityElement
52 
53 
54  //======================================================================
60  //======================================================================
61  template<unsigned DIM, unsigned NNODE_1D>
63  : public virtual TElement<DIM, NNODE_1D>,
64  public virtual LinearElasticityEquations<DIM>,
65  public virtual ElementWithZ2ErrorEstimator
66  {
67  public:
71  : TElement<DIM, NNODE_1D>(), LinearElasticityEquations<DIM>()
72  {
73  }
74 
75 
78  const TLinearElasticityElement<DIM, NNODE_1D>& dummy) = delete;
79 
81  // Commented out broken assignment operator because this can lead to a
82  // conflict warning when used in the virtual inheritence hierarchy.
83  // Essentially the compiler doesn't realise that two separate
84  // implementations of the broken function are the same and so, quite
85  // rightly, it shouts.
86  /*void operator=(const TLinearElasticityElement<DIM,NNODE_1D>&) = delete;*/
87 
89  void output(std::ostream& outfile)
90  {
92  }
93 
95  void output(std::ostream& outfile, const unsigned& nplot)
96  {
98  }
99 
100 
102  void output(FILE* file_pt)
103  {
105  }
106 
108  void output(FILE* file_pt, const unsigned& n_plot)
109  {
111  }
112 
114  unsigned nvertex_node() const
115  {
117  }
118 
120  Node* vertex_node_pt(const unsigned& j) const
121  {
123  }
124 
127  unsigned nrecovery_order()
128  {
129  return NNODE_1D - 1;
130  }
131 
133  unsigned num_Z2_flux_terms()
134  {
135  // DIM Diagonal strain rates and DIM*(DIM-1)/2 off diagonal terms
136  return (DIM + DIM * (DIM - 1) / 2);
137  }
138 
142  {
143 #ifdef PARANOID
144  unsigned num_entries = (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(),
154  }
155 #endif
156 
157  // Get strain matrix
158  DenseMatrix<double> strain(DIM);
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);
168  icount++;
169  }
170 
171  // Off diagonals row by row
172  for (unsigned i = 0; i < DIM; i++)
173  {
174  for (unsigned j = i + 1; j < DIM; j++)
175  {
176  flux[icount] = strain(i, j);
177  icount++;
178  }
179  }
180  }
181  };
182 
183  //=======================================================================
188  //=======================================================================
189  template<unsigned DIM, unsigned NNODE_1D>
191  : public virtual TElement<DIM - 1, NNODE_1D>
192  {
193  public:
196  FaceGeometry() : TElement<DIM - 1, NNODE_1D>() {}
197  };
198 
199  //=======================================================================
201  //=======================================================================
202  template<unsigned NNODE_1D>
204  : public virtual PointElement
205  {
206  public:
210  };
211 
212 
213 } // namespace oomph
214 
215 #endif
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Definition: error_estimator.h:79
FaceGeometry()
Definition: Tlinear_elasticity_elements.h:209
FaceGeometry()
Definition: Tlinear_elasticity_elements.h:196
Definition: elements.h:4998
void get_strain(const Vector< double > &s, DenseMatrix< double > &strain) const
Return the strain tensor.
Definition: linear_elasticity_elements.cc:47
Definition: linear_elasticity_elements.h:380
void output(std::ostream &outfile)
Output: x,y,[z],u,v,[w].
Definition: linear_elasticity_elements.h:427
Definition: nodes.h:906
Definition: oomph_definitions.h:222
Definition: elements.h:3439
Definition: Telements.h:1208
Definition: Tlinear_elasticity_elements.h:66
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
Definition: Tlinear_elasticity_elements.h:133
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: Tlinear_elasticity_elements.h:114
unsigned nrecovery_order()
Definition: Tlinear_elasticity_elements.h:127
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: Tlinear_elasticity_elements.h:120
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Definition: Tlinear_elasticity_elements.h:141
TLinearElasticityElement(const TLinearElasticityElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void output(std::ostream &outfile)
Broken assignment operator.
Definition: Tlinear_elasticity_elements.h:89
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function:
Definition: Tlinear_elasticity_elements.h:108
void output(std::ostream &outfile, const unsigned &nplot)
Output function:
Definition: Tlinear_elasticity_elements.h:95
void output(FILE *file_pt)
C-style output function:
Definition: Tlinear_elasticity_elements.h:102
TLinearElasticityElement()
Definition: Tlinear_elasticity_elements.h:70
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