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// 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_TLINEAR_ELASTICITY_ELEMENTS_HEADER
32 #define OOMPH_TLINEAR_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"
46 #include "../generic/error_estimator.h"
47 
48 
49 namespace oomph
50 {
51 
52 /////////////////////////////////////////////////////////////////////////
53 /////////////////////////////////////////////////////////////////////////
54 // TLinearElasticityElement
55 ////////////////////////////////////////////////////////////////////////
56 ////////////////////////////////////////////////////////////////////////
57 
58 
59 
60 //======================================================================
61 /// TLinearElasticityElement<DIM,NNODE_1D> elements are
62 /// isoparametric triangular
63 /// DIM-dimensional LinearElasticity elements with
64 /// NNODE_1D nodal points along each
65 /// element edge. Inherits from TElement and LinearElasticityEquations
66 //======================================================================
67 template <unsigned DIM, unsigned NNODE_1D>
68  class TLinearElasticityElement : public virtual TElement<DIM,NNODE_1D>,
69  public virtual LinearElasticityEquations<DIM>,
70  public virtual ElementWithZ2ErrorEstimator
71  {
72 
73  public:
74 
75  ///\short Constructor: Call constructors for TElement and
76  /// LinearElasticity equations
77  TLinearElasticityElement() : TElement<DIM,NNODE_1D>(),
78  LinearElasticityEquations<DIM>() { }
79 
80 
81  /// Broken copy constructor
83  {
84  BrokenCopy::broken_copy("TLinearElasticityElement");
85  }
86 
87  /// Broken assignment operator
88 //Commented out broken assignment operator because this can lead to a conflict warning
89 //when used in the virtual inheritence hierarchy. Essentially the compiler doesn't
90 //realise that two separate implementations of the broken function are the same and so,
91 //quite rightly, it shouts.
92  /*void operator=(const TLinearElasticityElement<DIM,NNODE_1D>&)
93  {
94  BrokenCopy::broken_assign("TLinearElasticityElement");
95  }*/
96 
97  /// \short Output function:
98  void output(std::ostream &outfile)
99  {
101  }
102 
103  /// \short Output function:
104  void output(std::ostream &outfile, const unsigned &nplot)
105  {
107  }
108 
109 
110  /// \short C-style output function:
111  void output(FILE* file_pt)
112  {
114  }
115 
116  /// \short C-style output function:
117  void output(FILE* file_pt, const unsigned &n_plot)
118  {
120  }
121 
122  /// \short Number of vertex nodes in the element
123  unsigned nvertex_node() const
125 
126  /// \short Pointer to the j-th vertex node in the element
127  Node* vertex_node_pt(const unsigned& j) const
128  {
130  }
131 
132  /// \short Order of recovery shape functions for Z2 error estimation:
133  /// Same order as shape functions.
134  unsigned nrecovery_order() {return NNODE_1D-1;}
135 
136  /// Number of 'flux' terms for Z2 error estimation
137  unsigned num_Z2_flux_terms()
138  {
139  // DIM Diagonal strain rates and DIM*(DIM-1)/2 off diagonal terms
140  return (DIM + DIM*(DIM-1)/2);
141  }
142 
143  /// \short Get 'flux' for Z2 error recovery: Upper triangular entries
144  /// in strain tensor.
146  {
147 #ifdef PARANOID
148  unsigned num_entries=(DIM+((DIM*DIM)-DIM)/2);
149  if (flux.size()!=num_entries)
150  {
151  std::ostringstream error_message;
152  error_message << "The flux vector has the wrong number of entries, "
153  << flux.size() << ", whereas it should be "
154  << num_entries << std::endl;
155  throw OomphLibError(
156  error_message.str(),
157  OOMPH_CURRENT_FUNCTION,
158  OOMPH_EXCEPTION_LOCATION);
159  }
160 #endif
161 
162  // Get strain matrix
163  DenseMatrix<double> strain(DIM);
164  this->get_strain(s,strain);
165 
166  // Pack into flux Vector
167  unsigned icount=0;
168 
169  // Start with diagonal terms
170  for(unsigned i=0;i<DIM;i++)
171  {
172  flux[icount]=strain(i,i);
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);
182  icount++;
183  }
184  }
185  }
186 };
187 
188 //=======================================================================
189 /// Face geometry for the TLinearElasticityElement elements: 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>
195 class FaceGeometry<TLinearElasticityElement<DIM,NNODE_1D> >:
196  public virtual TElement<DIM-1,NNODE_1D>
197 {
198 
199  public:
200 
201  /// \short Constructor: Call the constructor for the
202  /// appropriate lower-dimensional QElement
203  FaceGeometry() : TElement<DIM-1,NNODE_1D>() {}
204 
205 };
206 
207 //=======================================================================
208 /// Face geometry for the 1D TLinearElasticityElement elements: Point elements
209 //=======================================================================
210 template<unsigned NNODE_1D>
212  public virtual PointElement
213 {
214 
215  public:
216 
217  /// \short Constructor: Call the constructor for the
218  /// appropriate lower-dimensional TElement
220 
221 };
222 
223 
224 }
225 
226 #endif
void broken_copy(const std::string &class_name)
Issue error message and terminate execution.
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
cstr elem_len * i
Definition: cfortran.h:607
void output(FILE *file_pt)
C-style output function:
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
Definition: nodes.h:852
TLinearElasticityElement(const TLinearElasticityElement< DIM, NNODE_1D > &dummy)
Broken copy constructor.
void output(std::ostream &outfile)
Broken assignment operator.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional TElement. ...
void output(std::ostream &outfile, const unsigned &nplot)
Output function:
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function:
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement. ...
static char t char * s
Definition: cfortran.h:572
unsigned num_Z2_flux_terms()
Number of &#39;flux&#39; terms for Z2 error estimation.
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 get_strain(const Vector< double > &s, DenseMatrix< double > &strain) const
Return the strain tensor.
unsigned nvertex_node() const
Number of vertex nodes in the element.
void output(std::ostream &outfile)
Output: x,y,[z],u,v,[w].
TLinearElasticityElement()
Constructor: Call constructors for TElement and LinearElasticity equations.