MarmotFiniteElement.h
Go to the documentation of this file.
1 /* ---------------------------------------------------------------------
2  * _
3  * _ __ ___ __ _ _ __ _ __ ___ ___ | |_
4  * | '_ ` _ \ / _` | '__| '_ ` _ \ / _ \| __|
5  * | | | | | | (_| | | | | | | | | (_) | |_
6  * |_| |_| |_|\__,_|_| |_| |_| |_|\___/ \__|
7  *
8  * Unit of Strength of Materials and Structural Analysis
9  * University of Innsbruck,
10  * 2020 - today
11  *
12  * festigkeitslehre@uibk.ac.at
13  *
14  * Peter Gamnitzer peter.gamnitzer@uibk.ac.at
15  * Matthias Neuner matthias.neuner@uibk.ac.at
16  * Magdalena Schreter magdalena.schreter@uibk.ac.at
17  *
18  * This file is part of the MAteRialMOdellingToolbox (marmot).
19  *
20  * This library is free software; you can redistribute it and/or
21  * modify it under the terms of the GNU Lesser General Public
22  * License as published by the Free Software Foundation; either
23  * version 2.1 of the License, or (at your option) any later version.
24  *
25  * The full text of the license can be found in the file LICENSE.md at
26  * the top level directory of marmot.
27  * ---------------------------------------------------------------------
28  */
29 #pragma once
30 #include "Marmot/MarmotTypedefs.h"
31 #include <vector>
32 
33 namespace Marmot {
34 
35  namespace FiniteElement {
49  };
50 
52 
53  // 'Expanded' N , aka NBold aka multidimensional Interpolation Operator
54  Eigen::MatrixXd NB( const Eigen::VectorXd& N,
55  const int nDoFPerNode ); // Dynamic version
56 
57  template < int nDim, int nNodes >
58  Eigen::Matrix< double, nDim, nDim * nNodes > NB( const Eigen::Matrix< double, 1, nNodes >& N )
59  {
60  // Alternative Templated version of Interpolation operator NBold;
61  Eigen::Matrix< double, nDim, nDim* nNodes > N_ = Eigen::Matrix< double, nDim, nDim * nNodes >::Zero();
62  for ( int i = 0; i < nNodes; i++ ) {
63  for ( int j = 0; j < nDim; j++ ) {
64  N_( j, nDim * i + j ) = N( i );
65  }
66  }
67  return N_;
68  }
69 
70  Eigen::MatrixXd Jacobian( const Eigen::MatrixXd& dN_dXi,
71  const Eigen::VectorXd& coordinates ); // Dynamic version
72 
73  template < int nDim, int nNodes >
74  Eigen::Matrix< double, nDim, nDim > Jacobian( const Eigen::Matrix< double, nDim, nNodes >& dNdXi,
75  const Eigen::Matrix< double, nDim * nNodes, 1 >& coordinates )
76  {
77  // Alternative Templated version of Jacobian for compile time known
78  // sizes
79  Eigen::Matrix< double, nDim, nDim > J_ = Eigen::Matrix< double, nDim, nDim >::Zero();
80  for ( int i = 0; i < nDim; i++ ) // loop over global dimensions
81  for ( int j = 0; j < nDim; j++ ) // loop over local dimensions
82  for ( int k = 0; k < nNodes; k++ ) // Loop over nodes
83  J_( i, j ) += dNdXi( j, k ) * coordinates( i + k * nDim );
84  return J_;
85  }
86 
87  Eigen::VectorXi expandNodeIndicesToCoordinateIndices( const Eigen::VectorXi& nodeIndices, int nDim );
88 
89  namespace Spatial1D {
90  namespace Bar2 {
91 
92  constexpr int nNodes = 2;
93  using NSized = Eigen::Matrix< double, 1, nNodes >;
94  using dNdXiSized = Eigen::Matrix< double, 1, nNodes >;
95 
96  NSized N( double xi );
97  dNdXiSized dNdXi( double xi );
98  } // namespace Bar2
99 
100  namespace Bar3 {
101 
102  constexpr int nNodes = 3;
103  using NSized = Eigen::Matrix< double, 1, nNodes >;
104  using dNdXiSized = Eigen::Matrix< double, 1, nNodes >;
105 
106  NSized N( double xi );
107  dNdXiSized dNdXi( double xi );
108  } // namespace Bar3
109  } // namespace Spatial1D
110 
111  namespace Spatial2D {
112  constexpr int nDim = 2;
113  constexpr int voigtSize = 3;
114 
115  template < int nNodes >
116  Eigen::Matrix< double, voigtSize, nNodes * nDim > B( const Eigen::Matrix< double, nDim, nNodes >& dNdX )
117  {
118 
119  Eigen::Matrix< double, voigtSize, nNodes* nDim > B_ = Eigen::Matrix< double, voigtSize, nNodes * nDim >::Zero();
120  for ( int i = 0; i < nNodes; i++ ) {
121  B_( 0, 2 * i ) = dNdX( 0, i );
122  B_( 1, 2 * i + 1 ) = dNdX( 1, i );
123  B_( 2, 2 * i ) = dNdX( 1, i );
124  B_( 2, 2 * i + 1 ) = dNdX( 0, i );
125  }
126  return B_;
127  }
128 
129  template < int nNodes >
130  Eigen::Matrix< double, voigtSize, nNodes * nDim > BGreen( const Eigen::Matrix< double, nDim, nNodes >& dNdX,
131  const Eigen::Matrix2d& F )
132  {
133  // Green-Lagrange Strain Operator for given dNdX and Deformationgradient F
134  // Belytschko et. al pp. 213
135  Eigen::Matrix< double, voigtSize, nNodes* nDim > B_ = Eigen::Matrix< double, voigtSize, nNodes * nDim >::Zero();
136  for ( int i = 0; i < nNodes; i++ ) {
137  B_( 0, 2 * i ) = dNdX( 0, i ) * F( 0, 0 );
138  B_( 0, 2 * i + 1 ) = dNdX( 0, i ) * F( 1, 0 );
139  B_( 1, 2 * i ) = dNdX( 1, i ) * F( 0, 1 );
140  B_( 1, 2 * i + 1 ) = dNdX( 1, i ) * F( 1, 1 );
141  B_( 2, 2 * i ) = dNdX( 0, i ) * F( 0, 1 ) + dNdX( 1, i ) * F( 0, 0 );
142  B_( 2, 2 * i + 1 ) = dNdX( 0, i ) * F( 1, 1 ) + dNdX( 1, i ) * F( 1, 0 );
143  }
144  return B_;
145  }
146 
147  namespace Quad4 {
148  constexpr int nNodes = 4;
149 
150  using CoordinateSized = Eigen::Matrix< double, nNodes * nDim, 1 >;
151  using NSized = Eigen::Matrix< double, 1, nNodes >;
152  using dNdXiSized = Eigen::Matrix< double, nDim, nNodes >;
153 
154  NSized N( const Eigen::Vector2d& xi );
155  dNdXiSized dNdXi( const Eigen::Vector2d& xi );
156 
157  Eigen::Vector2i getBoundaryElementIndices( int faceID );
158  } // namespace Quad4
159 
160  namespace Quad8 {
161  constexpr int nNodes = 8;
162 
163  using CoordinateSized = Eigen::Matrix< double, nNodes * nDim, 1 >;
164  using NSized = Eigen::Matrix< double, 1, nNodes >;
165  using dNdXiSized = Eigen::Matrix< double, nDim, nNodes >;
166 
167  NSized N( const Eigen::Vector2d& xi );
168  dNdXiSized dNdXi( const Eigen::Vector2d& xi );
169 
170  Eigen::Vector3i getBoundaryElementIndices( int faceID );
171  } // namespace Quad8
172 
173  } // end of namespace Spatial2D
174 
175  namespace Spatial3D {
176  constexpr int nDim = 3;
177  constexpr int voigtSize = 6;
178 
179  template < int nNodes >
180  Eigen::Matrix< double, voigtSize, nNodes * nDim > B( const Eigen::Matrix< double, nDim, nNodes >& dNdX )
181  {
182  // ABAQUS like notation of strain: ( ε11, ε22, ε33, γ12, γ13, γ23 )
183  // _ _
184  // | dN/dx1 0 0 |
185  // | 0 dN/dx2 0 |
186  // | 0 0 dN/dx3 |
187  // | dN/dx2 dN/dx1 0 |
188  // | dN/dx3 0 dN/dx1 |
189  // |_ 0 dN/dx3 dN/dx2 _|
190 
191  Eigen::Matrix< double, voigtSize, nNodes* nDim > B_ = Eigen::Matrix< double, voigtSize, nNodes * nDim >::Zero();
192 
193  for ( int i = 0; i < nNodes; i++ ) {
194  B_( 0, nDim * i ) = dNdX( 0, i );
195  B_( 1, nDim * i + 1 ) = dNdX( 1, i );
196  B_( 2, nDim * i + 2 ) = dNdX( 2, i );
197  B_( 3, nDim * i + 0 ) = dNdX( 1, i );
198  B_( 3, nDim * i + 1 ) = dNdX( 0, i );
199  B_( 4, nDim * i + 0 ) = dNdX( 2, i );
200  B_( 4, nDim * i + 2 ) = dNdX( 0, i );
201  B_( 5, nDim * i + 1 ) = dNdX( 2, i );
202  B_( 5, nDim * i + 2 ) = dNdX( 1, i );
203  }
204 
205  return B_;
206  }
207 
208  template < int nNodes >
209  Eigen::Matrix< double, voigtSize, nNodes * nDim > BGreen( const Eigen::Matrix< double, nDim, nNodes >& dNdX,
210  const Eigen::Matrix3d& F )
211  {
212  // Green-Lagrange Strain Operator for given dNdX and Deformationgradient F
213  // Belytschko et. al pp. 213
214 
215  Eigen::Matrix< double, voigtSize, nNodes* nDim > B_ = Eigen::Matrix< double, voigtSize, nNodes * nDim >::Zero();
216  for ( int i = 0; i < nNodes; i++ ) {
217  B_( 0, nDim * i ) = dNdX( 0, i ) * F( 0, 0 );
218  B_( 0, nDim * i + 1 ) = dNdX( 0, i ) * F( 1, 0 );
219  B_( 0, nDim * i + 2 ) = dNdX( 0, i ) * F( 2, 0 );
220 
221  B_( 1, nDim * i ) = dNdX( 1, i ) * F( 0, 1 );
222  B_( 1, nDim * i + 1 ) = dNdX( 1, i ) * F( 1, 1 );
223  B_( 1, nDim * i + 2 ) = dNdX( 1, i ) * F( 2, 1 );
224 
225  B_( 2, nDim * i ) = dNdX( 2, i ) * F( 0, 2 );
226  B_( 2, nDim * i + 1 ) = dNdX( 2, i ) * F( 1, 2 );
227  B_( 2, nDim * i + 2 ) = dNdX( 2, i ) * F( 2, 2 );
228 
229  B_( 3, nDim * i ) = dNdX( 0, i ) * F( 0, 1 ) + dNdX( 1, i ) * F( 0, 0 );
230  B_( 3, nDim * i + 1 ) = dNdX( 0, i ) * F( 1, 1 ) + dNdX( 1, i ) * F( 1, 0 );
231  B_( 3, nDim * i + 2 ) = dNdX( 0, i ) * F( 2, 1 ) + dNdX( 1, i ) * F( 2, 0 );
232 
233  B_( 4, nDim * i ) = dNdX( 0, i ) * F( 0, 2 ) + dNdX( 2, i ) * F( 0, 0 );
234  B_( 4, nDim * i + 1 ) = dNdX( 0, i ) * F( 1, 2 ) + dNdX( 2, i ) * F( 1, 0 );
235  B_( 4, nDim * i + 2 ) = dNdX( 0, i ) * F( 2, 2 ) + dNdX( 2, i ) * F( 2, 0 );
236 
237  B_( 5, nDim * i ) = dNdX( 2, i ) * F( 0, 1 ) + dNdX( 1, i ) * F( 0, 2 );
238  B_( 5, nDim * i + 1 ) = dNdX( 2, i ) * F( 1, 1 ) + dNdX( 1, i ) * F( 1, 2 );
239  B_( 5, nDim * i + 2 ) = dNdX( 2, i ) * F( 2, 1 ) + dNdX( 1, i ) * F( 2, 2 );
240  }
241  return B_;
242  }
243 
244  namespace Tetra4 {
245 
246  constexpr int nNodes = 4;
247  using CoordinateSized = Eigen::Matrix< double, nNodes * nDim, 1 >;
248  using NSized = Eigen::Matrix< double, 1, nNodes >;
249  using dNdXiSized = Eigen::Matrix< double, nDim, nNodes >;
250 
251  NSized N( const Eigen::Vector3d& xi );
253 
254  Eigen::Vector3i getBoundaryElementIndices( int faceID );
255 
256  } // namespace Tetra4
257 
258  namespace Tetra10 {
259 
260  constexpr int nNodes = 10;
261  using CoordinateSized = Eigen::Matrix< double, nNodes * nDim, 1 >;
262  using NSized = Eigen::Matrix< double, 1, nNodes >;
263  using dNdXiSized = Eigen::Matrix< double, nDim, nNodes >;
264 
265  NSized N( const Eigen::Vector3d& xi );
267 
268  Eigen::Vector3i getBoundaryElementIndices( int faceID );
269 
270  } // namespace Tetra10
271 
272  namespace Hexa8 {
273  constexpr int nNodes = 8;
274  using CoordinateSized = Eigen::Matrix< double, nNodes * nDim, 1 >;
275  using NSized = Eigen::Matrix< double, 1, nNodes >;
276  using dNdXiSized = Eigen::Matrix< double, nDim, nNodes >;
277 
278  NSized N( const Eigen::Vector3d& xi );
280 
281  Eigen::Vector4i getBoundaryElementIndices( int faceID );
282  } // namespace Hexa8
283 
284  namespace Hexa20 {
285  constexpr int nNodes = 20;
286  using CoordinateSized = Eigen::Matrix< double, nNodes * nDim, 1 >;
287  using NSized = Eigen::Matrix< double, 1, nNodes >;
288  using dNdXiSized = Eigen::Matrix< double, nDim, nNodes >;
289 
290  NSized N( const Eigen::Vector3d& xi );
292 
294  } // namespace Hexa20
295  } // namespace Spatial3D
296 
298 
299  /* Boundary element, for instance for distributed surface loads
300  * */
301 
303  double weight;
304  double JxW;
305  Eigen::VectorXd xi;
306  Eigen::VectorXd N;
307  Eigen::MatrixXd dNdXi;
308  Eigen::MatrixXd dx_dXi;
309  Eigen::VectorXd areaVector;
310  };
311 
312  const int nDim;
313 
315  int nNodes;
317 
318  std::vector< BoundaryElementQuadraturePoint > quadraturePoints;
319 
320  Eigen::VectorXi mapBoundaryToParentScalar;
322  Eigen::VectorXd coordinates;
323 
324  public:
325  BoundaryElement( ElementShapes parentShape,
326  int nDim,
327  int parentFaceNumber,
328  const Eigen::VectorXd& parentCoordinates );
329 
330  Eigen::VectorXd computeScalarLoadVector();
331  Eigen::MatrixXd computeDScalarLoadVector_dCoordinates();
332 
334  Eigen::VectorXd computeSurfaceNormalVectorialLoadVector();
336 
338  Eigen::VectorXd computeVectorialLoadVector( const Eigen::VectorXd& direction );
339  Eigen::MatrixXd computeDVectorialLoadVector_dCoordinates( const Eigen::VectorXd& direction );
340 
341  Eigen::VectorXd condenseParentToBoundaryScalar( const Eigen::VectorXd& parentVector );
342  void assembleIntoParentScalar( const Eigen::VectorXd& boundaryVector,
343  Eigen::Ref< Eigen::VectorXd > ParentVector );
344  void assembleIntoParentStiffnessScalar( const Eigen::MatrixXd& KBoundary, Eigen::Ref< Eigen::MatrixXd > KParent );
345 
346  Eigen::VectorXd condenseParentToBoundaryVectorial( const Eigen::VectorXd& parentVector );
347  void assembleIntoParentVectorial( const Eigen::VectorXd& boundaryVector,
348  Eigen::Ref< Eigen::VectorXd > ParentVector );
349  void assembleIntoParentStiffnessVectorial( const Eigen::MatrixXd& KBoundary,
350  Eigen::Ref< Eigen::MatrixXd > KParent );
351  };
352  } // namespace FiniteElement
353 
354  namespace FiniteElement::Quadrature {
355  constexpr double gp2 = 0.577350269189625764509;
356  constexpr double gp3 = 0.774596669241483;
357 
359 
361  Eigen::VectorXd xi;
362  double weight;
363  };
364 
365  const std::vector< QuadraturePointInfo >& getGaussPointInfo( Marmot::FiniteElement::ElementShapes shape,
366  IntegrationTypes integrationType );
367 
369 
370  namespace Spatial1D {
371  constexpr int nDim = 1;
372 
373  // clang-format off
374  const std::vector< QuadraturePointInfo > gaussPointList1 = {
375  { ( Eigen::VectorXd ( 1 ) << 0 ).finished(), 2.0 }
376  };
377 
378  const std::vector< QuadraturePointInfo > gaussPointList2 = {
379  { ( Eigen::VectorXd ( 1 ) << -gp2 ).finished(), 1.0 },
380  { ( Eigen::VectorXd ( 1 ) << +gp2 ).finished(), 1.0 }
381  };
382 
383  const std::vector< QuadraturePointInfo > gaussPointList3 = {
384  { ( Eigen::VectorXd ( 1 ) << -gp3 ).finished(), 5./9 },
385  { ( Eigen::VectorXd ( 1 ) << 0. ).finished(), 8./9 },
386  { ( Eigen::VectorXd ( 1 ) << +gp3 ).finished(), 5./9 }
387  };
388  // clang-format on
389 
390  } // namespace Spatial1D
391 
392  namespace Spatial2D {
393  constexpr int nDim = 2;
394 
395  // clang-format off
396  const std::vector< QuadraturePointInfo > gaussPointList1x1 = {
397  { Eigen::Vector2d::Zero(), 4. }
398  };
399 
400  const std::vector< QuadraturePointInfo > gaussPointList2x2 = {
401  { ( Eigen::Vector2d () << +gp2, +gp2 ).finished(), 1.0 },
402  { ( Eigen::Vector2d () << -gp2, +gp2 ).finished(), 1.0 },
403  { ( Eigen::Vector2d () << -gp2, -gp2 ).finished(), 1.0 },
404  { ( Eigen::Vector2d () << +gp2, -gp2 ).finished(), 1.0 }
405  };
406 
407  const std::vector< QuadraturePointInfo > gaussPointList3x3 = {
408  { ( Eigen::Vector2d () << 0, 0. ).finished(), 64./81},
409  { ( Eigen::Vector2d () << -gp3, -gp3 ).finished(), 25./81},
410  { ( Eigen::Vector2d () << +gp3, -gp3 ).finished(), 25./81},
411  { ( Eigen::Vector2d () << +gp3, +gp3 ).finished(), 25./81},
412  { ( Eigen::Vector2d () << -gp3, +gp3 ).finished(), 25./81},
413  { ( Eigen::Vector2d () << 0, -gp3 ).finished(), 40./81},
414  { ( Eigen::Vector2d () << gp3, 0. ).finished(), 40./81},
415  { ( Eigen::Vector2d () << 0, +gp3 ).finished(), 40./81},
416  { ( Eigen::Vector2d () << -gp3, 0. ).finished(), 40./81},
417  };
418  // clang-format on
419 
420  void modifyCharElemLengthAbaqusLike( double& charElemLength, int intPoint );
421  } // namespace Spatial2D
422 
423  namespace Spatial3D {
424  constexpr int nDim = 3;
425 
426  // clang-format off
427  const inline std::vector< QuadraturePointInfo > gaussPointList1x1x1 = {
428  { Eigen::Vector3d::Zero(), 8.0 }
429  };
430 
431  const inline std::vector< QuadraturePointInfo > gaussPointListTetra4 = {
432  { (Eigen::Vector3d() << 1./4, 1./4, 1./4).finished(), 1./6}
433  };
434 
435  const inline std::vector< QuadraturePointInfo > gaussPointListTetra10 = {
436  { (Eigen::Vector3d() << (5-std::sqrt(5))/20, (5-std::sqrt(5))/20, (5-std::sqrt(5))/20 ).finished(), 1./24},
437  { (Eigen::Vector3d() << (5-std::sqrt(5))/20, (5-std::sqrt(5))/20, (5+3*std::sqrt(5))/20 ).finished(), 1./24},
438  { (Eigen::Vector3d() << (5-std::sqrt(5))/20, (5+3*std::sqrt(5))/20, (5-std::sqrt(5))/20 ).finished(), 1./24},
439  { (Eigen::Vector3d() << (5+3*std::sqrt(5))/20, (5-std::sqrt(5))/20, (5-std::sqrt(5))/20 ).finished(), 1./24},
440  };
441 
442  const inline std::vector< QuadraturePointInfo > gaussPointList2x2x2 = {
443  { ( Eigen::Vector3d () << -gp2, -gp2, -gp2 ).finished(), 1.0},
444  { ( Eigen::Vector3d () << +gp2, -gp2, -gp2 ).finished(), 1.0},
445  { ( Eigen::Vector3d () << +gp2, +gp2, -gp2 ).finished(), 1.0},
446  { ( Eigen::Vector3d () << -gp2, +gp2, -gp2 ).finished(), 1.0},
447  { ( Eigen::Vector3d () << -gp2, -gp2, +gp2 ).finished(), 1.0},
448  { ( Eigen::Vector3d () << +gp2, -gp2, +gp2 ).finished(), 1.0},
449  { ( Eigen::Vector3d () << +gp2, +gp2, +gp2 ).finished(), 1.0},
450  { ( Eigen::Vector3d () << -gp2, +gp2, +gp2 ).finished(), 1.0},
451  };
452 
453  const inline std::vector< QuadraturePointInfo > gaussPointList3x3x3 = {
454  { ( Eigen::Vector3d () << -gp3, -gp3, -gp3 ).finished(), 0.171467764060357},
455  { ( Eigen::Vector3d () << 0, -gp3, -gp3 ).finished(), 0.274348422496571},
456  { ( Eigen::Vector3d () << +gp3, -gp3, -gp3 ).finished(), 0.171467764060357},
457  { ( Eigen::Vector3d () << -gp3, 0, -gp3 ).finished(), 0.274348422496571},
458  { ( Eigen::Vector3d () << 0, 0, -gp3 ).finished(), 0.438957475994513},
459  { ( Eigen::Vector3d () << gp3, 0, -gp3 ).finished(), 0.274348422496571},
460  { ( Eigen::Vector3d () << -gp3, +gp3, -gp3 ).finished(), 0.171467764060357},
461  { ( Eigen::Vector3d () << 0, +gp3, -gp3 ).finished(), 0.274348422496571},
462  { ( Eigen::Vector3d () << +gp3, +gp3, -gp3 ).finished(), 0.171467764060357},
463 
464  { ( Eigen::Vector3d () << -gp3, -gp3, 0 ).finished(), 0.274348422496571},
465  { ( Eigen::Vector3d () << 0, -gp3, 0 ).finished(), 0.438957475994513},
466  { ( Eigen::Vector3d () << +gp3, -gp3, 0 ).finished(), 0.274348422496571},
467  { ( Eigen::Vector3d () << -gp3, 0, 0 ).finished(), 0.438957475994513},
468  { ( Eigen::Vector3d () << 0, 0, 0 ).finished(), 0.702331961591221},
469  { ( Eigen::Vector3d () << gp3, 0, 0 ).finished(), 0.438957475994513},
470  { ( Eigen::Vector3d () << -gp3, +gp3, 0 ).finished(), 0.274348422496571},
471  { ( Eigen::Vector3d () << 0, +gp3, 0 ).finished(), 0.438957475994513},
472  { ( Eigen::Vector3d () << +gp3, +gp3, 0 ).finished(), 0.274348422496571},
473 
474  { ( Eigen::Vector3d () << -gp3, -gp3, +gp3 ).finished(), 0.171467764060357},
475  { ( Eigen::Vector3d () << 0, -gp3, +gp3 ).finished(), 0.274348422496571},
476  { ( Eigen::Vector3d () << +gp3, -gp3, +gp3 ).finished(), 0.171467764060357},
477  { ( Eigen::Vector3d () << -gp3, 0, +gp3 ).finished(), 0.274348422496571},
478  { ( Eigen::Vector3d () << 0, 0, +gp3 ).finished(), 0.438957475994513},
479  { ( Eigen::Vector3d () << gp3, 0, +gp3 ).finished(), 0.274348422496571},
480  { ( Eigen::Vector3d () << -gp3, +gp3, +gp3 ).finished(), 0.171467764060357},
481  { ( Eigen::Vector3d () << 0, +gp3, +gp3 ).finished(), 0.274348422496571},
482  { ( Eigen::Vector3d () << +gp3, +gp3, +gp3 ).finished(), 0.171467764060357}
483  };
484  // clang-format on
485 
486  } // namespace Spatial3D
487  } // namespace FiniteElement::Quadrature
488 } // namespace Marmot
Marmot::FiniteElement::Spatial3D::Tetra10::N
NSized N(const Eigen::Vector3d &xi)
Marmot::FiniteElement::BoundaryElement::BoundaryElementQuadraturePoint
Definition: MarmotFiniteElement.h:302
Marmot::FiniteElement::Spatial2D::Quad8::NSized
Eigen::Matrix< double, 1, nNodes > NSized
Definition: MarmotFiniteElement.h:164
Marmot::FiniteElement::Spatial3D::Tetra10::NSized
Eigen::Matrix< double, 1, nNodes > NSized
Definition: MarmotFiniteElement.h:262
Marmot::FiniteElement::Spatial1D::Bar3::N
NSized N(double xi)
Definition: MarmotFiniteElement1D.cpp:41
Marmot::FiniteElement::Spatial2D::voigtSize
constexpr int voigtSize
Definition: MarmotFiniteElement.h:113
Marmot::FiniteElement::Quadrature::Spatial2D::gaussPointList3x3
const std::vector< QuadraturePointInfo > gaussPointList3x3
Definition: MarmotFiniteElement.h:407
Marmot::FiniteElement::BoundaryElement::computeScalarLoadVector
Eigen::VectorXd computeScalarLoadVector()
Definition: MarmotFiniteElementBoundary.cpp:124
Marmot::FiniteElement::Spatial3D::Hexa20::N
NSized N(const Eigen::Vector3d &xi)
Marmot::FiniteElement::BoundaryElement::computeDSurfaceNormalVectorialLoadVector_dCoordinates
Eigen::MatrixXd computeDSurfaceNormalVectorialLoadVector_dCoordinates()
Definition: MarmotFiniteElementBoundary.cpp:162
Marmot::FiniteElement::Spatial3D::Tetra4::nNodes
constexpr int nNodes
Definition: MarmotFiniteElement.h:246
Marmot::FiniteElement::BoundaryElement::BoundaryElement
BoundaryElement(ElementShapes parentShape, int nDim, int parentFaceNumber, const Eigen::VectorXd &parentCoordinates)
Definition: MarmotFiniteElementBoundary.cpp:10
Marmot::FiniteElement::Quadrature::FullIntegration
@ FullIntegration
Definition: MarmotFiniteElement.h:358
Marmot::FiniteElement::BoundaryElement::computeDScalarLoadVector_dCoordinates
Eigen::MatrixXd computeDScalarLoadVector_dCoordinates()
Definition: MarmotFiniteElementBoundary.cpp:139
Marmot::FiniteElement::BoundaryElement::condenseParentToBoundaryScalar
Eigen::VectorXd condenseParentToBoundaryScalar(const Eigen::VectorXd &parentVector)
Definition: MarmotFiniteElementBoundary.cpp:229
Marmot::FiniteElement::BoundaryElement::assembleIntoParentScalar
void assembleIntoParentScalar(const Eigen::VectorXd &boundaryVector, Eigen::Ref< Eigen::VectorXd > ParentVector)
Definition: MarmotFiniteElementBoundary.cpp:242
Marmot::FiniteElement::Spatial3D::B
Eigen::Matrix< double, voigtSize, nNodes *nDim > B(const Eigen::Matrix< double, nDim, nNodes > &dNdX)
Definition: MarmotFiniteElement.h:180
Marmot::FiniteElement::BoundaryElement::BoundaryElementQuadraturePoint::weight
double weight
Definition: MarmotFiniteElement.h:303
Marmot::FiniteElement::Spatial3D::Hexa20::dNdXiSized
Eigen::Matrix< double, nDim, nNodes > dNdXiSized
Definition: MarmotFiniteElement.h:288
Marmot::FiniteElement::Quadrature::Spatial3D::gaussPointList1x1x1
const std::vector< QuadraturePointInfo > gaussPointList1x1x1
Definition: MarmotFiniteElement.h:427
Marmot::FiniteElement::BoundaryElement::coordinates
Eigen::VectorXd coordinates
Definition: MarmotFiniteElement.h:322
Marmot::FiniteElement::Spatial1D::Bar2::dNdXiSized
Eigen::Matrix< double, 1, nNodes > dNdXiSized
Definition: MarmotFiniteElement.h:94
Marmot::FiniteElement::BoundaryElement::quadraturePoints
std::vector< BoundaryElementQuadraturePoint > quadraturePoints
Definition: MarmotFiniteElement.h:318
Marmot::FiniteElement::Quadrature::ReducedIntegration
@ ReducedIntegration
Definition: MarmotFiniteElement.h:358
Marmot::FiniteElement::Spatial2D::Quad4::dNdXiSized
Eigen::Matrix< double, nDim, nNodes > dNdXiSized
Definition: MarmotFiniteElement.h:152
Marmot::FiniteElement::Spatial2D::Quad8::dNdXiSized
Eigen::Matrix< double, nDim, nNodes > dNdXiSized
Definition: MarmotFiniteElement.h:165
Marmot::FiniteElement::Quad9
@ Quad9
Definition: MarmotFiniteElement.h:41
MarmotTypedefs.h
Marmot::FiniteElement::Quadrature::gp2
constexpr double gp2
Definition: MarmotFiniteElement.h:355
Marmot::FiniteElement::Spatial3D::Tetra4::CoordinateSized
Eigen::Matrix< double, nNodes *nDim, 1 > CoordinateSized
Definition: MarmotFiniteElement.h:247
Marmot::FiniteElement::Hexa20
@ Hexa20
Definition: MarmotFiniteElement.h:46
Marmot::FiniteElement::Spatial2D::Quad8::N
NSized N(const Eigen::Vector2d &xi)
Marmot::FiniteElement::Spatial3D::Hexa8::NSized
Eigen::Matrix< double, 1, nNodes > NSized
Definition: MarmotFiniteElement.h:275
Marmot::FiniteElement::Spatial1D::Bar3::NSized
Eigen::Matrix< double, 1, nNodes > NSized
Definition: MarmotFiniteElement.h:103
Marmot::FiniteElement::Bar3
@ Bar3
Definition: MarmotFiniteElement.h:38
Marmot::FiniteElement::BoundaryElement::BoundaryElementQuadraturePoint::N
Eigen::VectorXd N
Definition: MarmotFiniteElement.h:306
Marmot::FiniteElement::BoundaryElement::assembleIntoParentVectorial
void assembleIntoParentVectorial(const Eigen::VectorXd &boundaryVector, Eigen::Ref< Eigen::VectorXd > ParentVector)
Definition: MarmotFiniteElementBoundary.cpp:274
Marmot::FiniteElement::BoundaryElement::condenseParentToBoundaryVectorial
Eigen::VectorXd condenseParentToBoundaryVectorial(const Eigen::VectorXd &parentVector)
Definition: MarmotFiniteElementBoundary.cpp:261
Marmot::FiniteElement::Spatial3D::Tetra10::nNodes
constexpr int nNodes
Definition: MarmotFiniteElement.h:260
Marmot::FiniteElement::BoundaryElement
Definition: MarmotFiniteElement.h:297
Marmot::FiniteElement::BoundaryElement::BoundaryElementQuadraturePoint::dNdXi
Eigen::MatrixXd dNdXi
Definition: MarmotFiniteElement.h:307
Marmot::FiniteElement::Quad4
@ Quad4
Definition: MarmotFiniteElement.h:39
Marmot::FiniteElement::Spatial2D::nDim
constexpr int nDim
Definition: MarmotFiniteElement.h:112
Marmot::FiniteElement::Spatial3D::Hexa20::nNodes
constexpr int nNodes
Definition: MarmotFiniteElement.h:285
Marmot::Vector3d
Eigen::Matrix< double, 3, 1 > Vector3d
Definition: MarmotTypedefs.h:42
Marmot::FiniteElement::BoundaryElement::nParentCoordinates
int nParentCoordinates
Definition: MarmotFiniteElement.h:316
Marmot::FiniteElement::Quadrature::Spatial1D::gaussPointList1
const std::vector< QuadraturePointInfo > gaussPointList1
Definition: MarmotFiniteElement.h:374
Marmot::FiniteElement::Spatial3D::Tetra10::CoordinateSized
Eigen::Matrix< double, nNodes *nDim, 1 > CoordinateSized
Definition: MarmotFiniteElement.h:261
Marmot::FiniteElement::Spatial1D::Bar2::NSized
Eigen::Matrix< double, 1, nNodes > NSized
Definition: MarmotFiniteElement.h:93
Marmot::FiniteElement::Spatial2D::Quad4::CoordinateSized
Eigen::Matrix< double, nNodes *nDim, 1 > CoordinateSized
Definition: MarmotFiniteElement.h:150
Marmot::FiniteElement::Tetra4
@ Tetra4
Definition: MarmotFiniteElement.h:43
Marmot::FiniteElement::Spatial3D::BGreen
Eigen::Matrix< double, voigtSize, nNodes *nDim > BGreen(const Eigen::Matrix< double, nDim, nNodes > &dNdX, const Eigen::Matrix3d &F)
Definition: MarmotFiniteElement.h:209
Marmot::FiniteElement::BoundaryElement::computeVectorialLoadVector
Eigen::VectorXd computeVectorialLoadVector(const Eigen::VectorXd &direction)
compute the element load vector for a unit vectorial load in a given direction.
Definition: MarmotFiniteElementBoundary.cpp:206
Marmot::FiniteElement::Quadrature::Spatial2D::gaussPointList2x2
const std::vector< QuadraturePointInfo > gaussPointList2x2
Definition: MarmotFiniteElement.h:400
Marmot::FiniteElement::Spatial3D::Hexa20::NSized
Eigen::Matrix< double, 1, nNodes > NSized
Definition: MarmotFiniteElement.h:287
Marmot::FiniteElement::Quadrature::Spatial1D::nDim
constexpr int nDim
Definition: MarmotFiniteElement.h:371
Marmot::FiniteElement::Spatial3D::Tetra4::N
NSized N(const Eigen::Vector3d &xi)
Marmot::FiniteElement::Spatial2D::Quad4::nNodes
constexpr int nNodes
Definition: MarmotFiniteElement.h:148
Marmot::FiniteElement::Hexa8
@ Hexa8
Definition: MarmotFiniteElement.h:45
Marmot::FiniteElement::Spatial3D::voigtSize
constexpr int voigtSize
Definition: MarmotFiniteElement.h:177
Marmot::FiniteElement::Spatial2D::Quad8::dNdXi
dNdXiSized dNdXi(const Eigen::Vector2d &xi)
Marmot::FiniteElement::Spatial2D::Quad8::nNodes
constexpr int nNodes
Definition: MarmotFiniteElement.h:161
Marmot::FiniteElement::Quadrature::Spatial2D::gaussPointList1x1
const std::vector< QuadraturePointInfo > gaussPointList1x1
Definition: MarmotFiniteElement.h:396
Marmot::FiniteElement::Spatial3D::Tetra4::dNdXiSized
Eigen::Matrix< double, nDim, nNodes > dNdXiSized
Definition: MarmotFiniteElement.h:249
Marmot::FiniteElement::Spatial3D::nDim
constexpr int nDim
Definition: MarmotFiniteElement.h:176
Marmot::Matrix3d
Eigen::Matrix< double, 3, 3 > Matrix3d
Definition: MarmotTypedefs.h:40
Marmot::FiniteElement::Bar2
@ Bar2
Definition: MarmotFiniteElement.h:37
Marmot::FiniteElement::BoundaryElement::computeDVectorialLoadVector_dCoordinates
Eigen::MatrixXd computeDVectorialLoadVector_dCoordinates(const Eigen::VectorXd &direction)
Definition: MarmotFiniteElementBoundary.cpp:221
Marmot
This file includes functions needed for calculations with stress and strain tensors written in voigt ...
Definition: MarmotTesting.h:30
Marmot::FiniteElement::Quadrature::Spatial3D::gaussPointList3x3x3
const std::vector< QuadraturePointInfo > gaussPointList3x3x3
Definition: MarmotFiniteElement.h:453
Marmot::FiniteElement::Spatial2D::Quad8::CoordinateSized
Eigen::Matrix< double, nNodes *nDim, 1 > CoordinateSized
Definition: MarmotFiniteElement.h:163
Marmot::FiniteElement::BoundaryElement::computeSurfaceNormalVectorialLoadVector
Eigen::VectorXd computeSurfaceNormalVectorialLoadVector()
compute the element load vector for a unit vectorial load normal to the surface.
Definition: MarmotFiniteElementBoundary.cpp:147
Marmot::FiniteElement::Spatial2D::B
Eigen::Matrix< double, voigtSize, nNodes *nDim > B(const Eigen::Matrix< double, nDim, nNodes > &dNdX)
Definition: MarmotFiniteElement.h:116
Marmot::FiniteElement::expandNodeIndicesToCoordinateIndices
Eigen::VectorXi expandNodeIndicesToCoordinateIndices(const Eigen::VectorXi &nodeIndices, int nDim)
Marmot::FiniteElement::Tetra10
@ Tetra10
Definition: MarmotFiniteElement.h:44
Marmot::FiniteElement::Quadrature::Spatial1D::gaussPointList2
const std::vector< QuadraturePointInfo > gaussPointList2
Definition: MarmotFiniteElement.h:378
Marmot::FiniteElement::Spatial2D::Quad8::getBoundaryElementIndices
Eigen::Vector3i getBoundaryElementIndices(int faceID)
Definition: MarmotFiniteElement2D.cpp:131
Marmot::FiniteElement::BoundaryElement::BoundaryElementQuadraturePoint::JxW
double JxW
Definition: MarmotFiniteElement.h:304
Marmot::FiniteElement::Spatial3D::Hexa8::N
NSized N(const Eigen::Vector3d &xi)
Marmot::FiniteElement::BoundaryElement::BoundaryElementQuadraturePoint::xi
Eigen::VectorXd xi
Definition: MarmotFiniteElement.h:305
Marmot::FiniteElement::Quadrature::Spatial3D::nDim
constexpr int nDim
Definition: MarmotFiniteElement.h:424
Marmot::FiniteElement::Quadrature::Spatial1D::gaussPointList3
const std::vector< QuadraturePointInfo > gaussPointList3
Definition: MarmotFiniteElement.h:383
Marmot::FiniteElement::BoundaryElement::assembleIntoParentStiffnessScalar
void assembleIntoParentStiffnessScalar(const Eigen::MatrixXd &KBoundary, Eigen::Ref< Eigen::MatrixXd > KParent)
Definition: MarmotFiniteElementBoundary.cpp:252
Marmot::FiniteElement::NB
Eigen::MatrixXd NB(const Eigen::VectorXd &N, const int nDoFPerNode)
Marmot::FiniteElement::BoundaryElement::BoundaryElementQuadraturePoint::areaVector
Eigen::VectorXd areaVector
Definition: MarmotFiniteElement.h:309
Marmot::FiniteElement::BoundaryElement::mapBoundaryToParentVectorial
Eigen::VectorXi mapBoundaryToParentVectorial
Definition: MarmotFiniteElement.h:321
Marmot::FiniteElement::ElementShapes
ElementShapes
Definition: MarmotFiniteElement.h:36
Marmot::FiniteElement::Quad16
@ Quad16
Definition: MarmotFiniteElement.h:42
Marmot::FiniteElement::Hexa64
@ Hexa64
Definition: MarmotFiniteElement.h:48
Marmot::FiniteElement::Spatial3D::Tetra4::getBoundaryElementIndices
Eigen::Vector3i getBoundaryElementIndices(int faceID)
Definition: MarmotFiniteElement3D.cpp:68
Marmot::Vector8i
Eigen::Matrix< int, 8, 1 > Vector8i
Definition: MarmotTypedefs.h:47
Marmot::FiniteElement::Spatial2D::Quad4::dNdXi
dNdXiSized dNdXi(const Eigen::Vector2d &xi)
Marmot::FiniteElement::Spatial3D::Tetra4::NSized
Eigen::Matrix< double, 1, nNodes > NSized
Definition: MarmotFiniteElement.h:248
Marmot::FiniteElement::Spatial3D::Hexa20::dNdXi
dNdXiSized dNdXi(const Eigen::Vector3d &xi)
Marmot::FiniteElement::Spatial1D::Bar2::nNodes
constexpr int nNodes
Definition: MarmotFiniteElement.h:92
Marmot::FiniteElement::Quadrature::getNumGaussPoints
int getNumGaussPoints(Marmot::FiniteElement::ElementShapes shape, IntegrationTypes integrationType)
Definition: MarmotFiniteElementBasic.cpp:160
Marmot::FiniteElement::Quadrature::getGaussPointInfo
const std::vector< QuadraturePointInfo > & getGaussPointInfo(Marmot::FiniteElement::ElementShapes shape, IntegrationTypes integrationType)
Definition: MarmotFiniteElementBasic.cpp:100
Marmot::FiniteElement::Spatial3D::Hexa20::CoordinateSized
Eigen::Matrix< double, nNodes *nDim, 1 > CoordinateSized
Definition: MarmotFiniteElement.h:286
Marmot::FiniteElement::Spatial2D::Quad4::getBoundaryElementIndices
Eigen::Vector2i getBoundaryElementIndices(int faceID)
Definition: MarmotFiniteElement2D.cpp:46
Marmot::FiniteElement::Spatial3D::Tetra10::dNdXiSized
Eigen::Matrix< double, nDim, nNodes > dNdXiSized
Definition: MarmotFiniteElement.h:263
Marmot::FiniteElement::Spatial3D::Tetra10::getBoundaryElementIndices
Eigen::Vector3i getBoundaryElementIndices(int faceID)
Definition: MarmotFiniteElement3D.cpp:177
Marmot::FiniteElement::getElementShapeByMetric
ElementShapes getElementShapeByMetric(int nDim, int nNodes)
Definition: MarmotFiniteElementBasic.cpp:7
Marmot::FiniteElement::Spatial1D::Bar2::dNdXi
dNdXiSized dNdXi(double xi)
Definition: MarmotFiniteElement1D.cpp:22
Marmot::FiniteElement::Quadrature::Spatial3D::gaussPointListTetra10
const std::vector< QuadraturePointInfo > gaussPointListTetra10
Definition: MarmotFiniteElement.h:435
Marmot::FiniteElement::Spatial3D::Hexa20::getBoundaryElementIndices
Marmot::Vector8i getBoundaryElementIndices(int faceID)
Definition: MarmotFiniteElement3D.cpp:419
Marmot::FiniteElement::Hexa27
@ Hexa27
Definition: MarmotFiniteElement.h:47
Marmot::FiniteElement::Spatial1D::Bar3::nNodes
constexpr int nNodes
Definition: MarmotFiniteElement.h:102
Marmot::FiniteElement::BoundaryElement::boundaryShape
ElementShapes boundaryShape
Definition: MarmotFiniteElement.h:314
Marmot::FiniteElement::Spatial3D::Tetra10::dNdXi
dNdXiSized dNdXi(const Eigen::Vector3d &xi)
Marmot::FiniteElement::Spatial3D::Hexa8::dNdXi
dNdXiSized dNdXi(const Eigen::Vector3d &xi)
Marmot::FiniteElement::BoundaryElement::nNodes
int nNodes
Definition: MarmotFiniteElement.h:315
Marmot::FiniteElement::EAS::F
Eigen::MatrixXd F(const Eigen::MatrixXd &J)
Marmot::FiniteElement::Spatial3D::Hexa8::getBoundaryElementIndices
Eigen::Vector4i getBoundaryElementIndices(int faceID)
Definition: MarmotFiniteElement3D.cpp:261
Marmot::FiniteElement::Quadrature::IntegrationTypes
IntegrationTypes
Definition: MarmotFiniteElement.h:358
Marmot::FiniteElement::Quadrature::gp3
constexpr double gp3
Definition: MarmotFiniteElement.h:356
Marmot::FiniteElement::BoundaryElement::mapBoundaryToParentScalar
Eigen::VectorXi mapBoundaryToParentScalar
Definition: MarmotFiniteElement.h:320
Marmot::FiniteElement::Quadrature::Spatial2D::nDim
constexpr int nDim
Definition: MarmotFiniteElement.h:393
Marmot::FiniteElement::Quadrature::QuadraturePointInfo::weight
double weight
Definition: MarmotFiniteElement.h:362
Marmot::FiniteElement::Quadrature::QuadraturePointInfo::xi
Eigen::VectorXd xi
Definition: MarmotFiniteElement.h:361
Marmot::FiniteElement::Spatial1D::Bar3::dNdXi
dNdXiSized dNdXi(double xi)
Definition: MarmotFiniteElement1D.cpp:48
Marmot::FiniteElement::BoundaryElement::assembleIntoParentStiffnessVectorial
void assembleIntoParentStiffnessVectorial(const Eigen::MatrixXd &KBoundary, Eigen::Ref< Eigen::MatrixXd > KParent)
Definition: MarmotFiniteElementBoundary.cpp:284
Marmot::FiniteElement::Spatial1D::Bar3::dNdXiSized
Eigen::Matrix< double, 1, nNodes > dNdXiSized
Definition: MarmotFiniteElement.h:104
Marmot::FiniteElement::Quadrature::Spatial3D::gaussPointList2x2x2
const std::vector< QuadraturePointInfo > gaussPointList2x2x2
Definition: MarmotFiniteElement.h:442
Marmot::FiniteElement::BoundaryElement::BoundaryElementQuadraturePoint::dx_dXi
Eigen::MatrixXd dx_dXi
Definition: MarmotFiniteElement.h:308
Marmot::FiniteElement::Spatial3D::Tetra4::dNdXi
dNdXiSized dNdXi(const Eigen::Vector3d &xi)
Marmot::FiniteElement::Spatial3D::Hexa8::nNodes
constexpr int nNodes
Definition: MarmotFiniteElement.h:273
Marmot::FiniteElement::Quad8
@ Quad8
Definition: MarmotFiniteElement.h:40
Marmot::FiniteElement::Jacobian
Eigen::MatrixXd Jacobian(const Eigen::MatrixXd &dN_dXi, const Eigen::VectorXd &coordinates)
Marmot::FiniteElement::Spatial2D::Quad4::NSized
Eigen::Matrix< double, 1, nNodes > NSized
Definition: MarmotFiniteElement.h:151
Marmot::FiniteElement::Spatial1D::Bar2::N
NSized N(double xi)
Definition: MarmotFiniteElement1D.cpp:15
Marmot::FiniteElement::Quadrature::QuadraturePointInfo
Definition: MarmotFiniteElement.h:360
Marmot::FiniteElement::BoundaryElement::nDim
const int nDim
Definition: MarmotFiniteElement.h:312
Marmot::FiniteElement::Spatial2D::BGreen
Eigen::Matrix< double, voigtSize, nNodes *nDim > BGreen(const Eigen::Matrix< double, nDim, nNodes > &dNdX, const Eigen::Matrix2d &F)
Definition: MarmotFiniteElement.h:130
Marmot::FiniteElement::Spatial2D::Quad4::N
NSized N(const Eigen::Vector2d &xi)
Marmot::FiniteElement::Spatial3D::Hexa8::dNdXiSized
Eigen::Matrix< double, nDim, nNodes > dNdXiSized
Definition: MarmotFiniteElement.h:276
Marmot::FiniteElement::Quadrature::Spatial3D::gaussPointListTetra4
const std::vector< QuadraturePointInfo > gaussPointListTetra4
Definition: MarmotFiniteElement.h:431
Marmot::FiniteElement::Quadrature::Spatial2D::modifyCharElemLengthAbaqusLike
void modifyCharElemLengthAbaqusLike(double &charElemLength, int intPoint)
Marmot::FiniteElement::Spatial3D::Hexa8::CoordinateSized
Eigen::Matrix< double, nNodes *nDim, 1 > CoordinateSized
Definition: MarmotFiniteElement.h:274