Taxisym_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 axisym linear elasticity elements
27 #ifndef OOMPH_TAXISYM_LINEAR_ELASTICITY_ELEMENTS_HEADER
28 #define OOMPH_TAXISYM_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 
45 namespace oomph
46 {
49  // axisym TLinearElasticityElement
52 
53 
54  //======================================================================
60  //======================================================================
61  template<unsigned NNODE_1D>
63  : public virtual TElement<2, NNODE_1D>,
65  public virtual ElementWithZ2ErrorEstimator
66  {
67  public:
72  {
73  }
74 
75 
79 
82  delete;
83 
85  void output(std::ostream& outfile)
86  {
88  }
89 
91  void output(std::ostream& outfile, const unsigned& nplot)
92  {
94  }
95 
96 
98  void output(FILE* file_pt)
99  {
101  }
102 
104  void output(FILE* file_pt, const unsigned& n_plot)
105  {
107  }
108 
110  unsigned nvertex_node() const
111  {
113  }
114 
116  Node* vertex_node_pt(const unsigned& j) const
117  {
119  }
120 
123  unsigned nrecovery_order()
124  {
125  return NNODE_1D - 1;
126  }
127 
129  unsigned num_Z2_flux_terms()
130  {
131  return 6;
132  }
133 
137  {
138 #ifdef PARANOID
139  unsigned num_entries = 6;
140  if (flux.size() != num_entries)
141  {
142  std::ostringstream error_message;
143  error_message << "The flux vector has the wrong number of entries, "
144  << flux.size() << ", whereas it should be " << num_entries
145  << std::endl;
146  throw OomphLibError(error_message.str(),
149  }
150 #endif
151 
152  // Get strain matrix
153  DenseMatrix<double> strain(3);
154  this->get_strain(s, strain);
155 
156  // Pack into flux Vector
157  unsigned icount = 0;
158 
159  // Start with diagonal terms
160  for (unsigned i = 0; i < 3; i++)
161  {
162  flux[icount] = strain(i, i);
163  icount++;
164  }
165 
166  // Off diagonals row by row
167  for (unsigned i = 0; i < 3; i++)
168  {
169  for (unsigned j = i + 1; j < 3; j++)
170  {
171  flux[icount] = strain(i, j);
172  icount++;
173  }
174  }
175  }
176  };
177 
178  //=======================================================================
182  //=======================================================================
183  template<unsigned NNODE_1D>
185  : public virtual TElement<1, NNODE_1D>
186  {
187  public:
190  FaceGeometry() : TElement<1, NNODE_1D>() {}
191  };
192 
193 
194 } // namespace oomph
195 
196 #endif
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Definition: axisym_linear_elasticity_elements.h:434
void get_strain(const Vector< double > &s, DenseMatrix< double > &strain)
Get strain (3x3 entries; r, z, phi)
Definition: axisym_linear_elasticity_elements.cc:55
void output(std::ostream &outfile)
Output: r,z, u_r, u_z, u_theta.
Definition: axisym_linear_elasticity_elements.h:483
Definition: error_estimator.h:79
FaceGeometry()
Definition: Taxisym_linear_elasticity_elements.h:190
Definition: elements.h:4998
Definition: nodes.h:906
Definition: oomph_definitions.h:222
Definition: Taxisym_linear_elasticity_elements.h:66
unsigned nrecovery_order()
Definition: Taxisym_linear_elasticity_elements.h:123
TAxisymmetricLinearElasticityElement(const TAxisymmetricLinearElasticityElement< NNODE_1D > &dummy)=delete
Broken copy constructor.
TAxisymmetricLinearElasticityElement()
Definition: Taxisym_linear_elasticity_elements.h:70
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: Taxisym_linear_elasticity_elements.h:110
void output(std::ostream &outfile)
Output function:
Definition: Taxisym_linear_elasticity_elements.h:85
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function:
Definition: Taxisym_linear_elasticity_elements.h:104
void output(FILE *file_pt)
C-style output function:
Definition: Taxisym_linear_elasticity_elements.h:98
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
Definition: Taxisym_linear_elasticity_elements.h:129
void operator=(const TAxisymmetricLinearElasticityElement< NNODE_1D > &)=delete
Broken assignment operator.
void output(std::ostream &outfile, const unsigned &nplot)
Output function:
Definition: Taxisym_linear_elasticity_elements.h:91
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: Taxisym_linear_elasticity_elements.h:116
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Definition: Taxisym_linear_elasticity_elements.h:136
Definition: Telements.h:1208
RealScalar s
Definition: level1_cplx_impl.h:130
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