MarmotGeometryElement.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  * Matthias Neuner matthias.neuner@uibk.ac.at
15  * Magdalena Schreter magdalena.schreter@uibk.ac.at
16  *
17  * This file is part of the MAteRialMOdellingToolbox (marmot).
18  *
19  * This library is free software; you can redistribute it and/or
20  * modify it under the terms of the GNU Lesser General Public
21  * License as published by the Free Software Foundation; either
22  * version 2.1 of the License, or (at your option) any later version.
23  *
24  * The full text of the license can be found in the file LICENSE.md at
25  * the top level directory of marmot.
26  * ---------------------------------------------------------------------
27  */
28 #pragma once
30 #include "Marmot/MarmotTypedefs.h"
31 #include "Marmot/MarmotVoigt.h"
32 #include <iostream>
33 #include <map>
34 
35 template < int nDim, int nNodes >
37  /* This is the Geometry Base element, which serves as a base for all MarmotElements.
38  * It corresponds to the GeometryElement in mpFEM,
39  * although this as a static templated version.
40  *
41  * MarmotElements (corresponding do DofElements in mpFEM) can inherit from this element,
42  * and access shape functions, derivatives and B Operator
43  *
44  * The element automatically determines its shape by the given nDimension and number of nodes
45  * */
46 
47 public:
48  /*Typedefs*/
49  /* static constexpr VoigtSize voigtSize = ( ( ( nDim * nDim ) + nDim ) / 2 ); */
52 
53  typedef Eigen::Matrix< double, nDim, 1 > XiSized;
54  typedef Eigen::Matrix< double, nDim * nNodes, 1 > CoordinateVector;
55  typedef Eigen::Matrix< double, nDim, nDim > JacobianSized;
56  typedef Eigen::Matrix< double, 1, nNodes > NSized;
57  typedef Eigen::Matrix< double, nDim, nNodes * nDim > NBSized;
58  typedef Eigen::Matrix< double, nDim, nNodes > dNdXiSized;
59  typedef Eigen::Matrix< double, voigtSize, nNodes * nDim > BSized;
60 
61  /*Properties*/
62  Eigen::Map< const CoordinateVector > coordinates;
64 
65  /*Methods*/
67  : coordinates( nullptr ), shape( Marmot::FiniteElement::getElementShapeByMetric( nDim, nNodes ) ){};
68 
69  std::string getElementShape() const
70  {
71  using namespace Marmot::FiniteElement;
72  static std::map< ElementShapes, std::string > shapes = { { Bar2, "bar2" },
73  { Quad4, "quad4" },
74  { Quad8, "quad8" },
75  { Tetra4, "tetra4" },
76  { Tetra10, "tetra10" },
77  { Hexa8, "hexa8" },
78  { Hexa20, "hexa20" } };
79 
80  return shapes[this->shape];
81  }
82 
83  void assignNodeCoordinates( const double* coords )
84  {
85  new ( &coordinates ) Eigen::Map< const CoordinateVector >( coords );
86  }
87 
88  /*Please specialize these functions for each element individially
89  *.cpp file.
90  *Fully specialized templates are precompiled in marmotMechanics (rather than the unspecialized and
91  *partially specialized templates)
92  * */
93  NSized N( const XiSized& xi ) const;
94  dNdXiSized dNdXi( const XiSized& xi ) const;
95  BSized B( const dNdXiSized& dNdX ) const;
96  BSized BGreen( const dNdXiSized& dNdX, const JacobianSized& F ) const;
97 
98  /*These functions are equal for each element and independent of node number and nDimension*/
99  NBSized NB( const NSized& N ) const { return Marmot::FiniteElement::NB< nDim, nNodes >( N ); }
100 
102  {
103  return Marmot::FiniteElement::Jacobian< nDim, nNodes >( dNdXi, coordinates );
104  }
105 
106  dNdXiSized dNdX( const dNdXiSized& dNdXi, const JacobianSized& JacobianInverse ) const
107  {
108  return ( dNdXi.transpose() * JacobianInverse ).transpose();
109  }
110 
111  JacobianSized F( const dNdXiSized& dNdX, const CoordinateVector& Q ) const
112  {
113  return Marmot::FiniteElement::Jacobian< nDim, nNodes >( dNdX, Q ) + JacobianSized::Identity();
114  }
115 };
MarmotGeometryElement::dNdXi
dNdXiSized dNdXi(const XiSized &xi) const
MarmotGeometryElement::B
BSized B(const dNdXiSized &dNdX) const
Marmot::ContinuumMechanics::VoigtNotation::voigtSizeFromDimension
constexpr VoigtSize voigtSizeFromDimension(int x)
Definition: MarmotVoigt.h:47
MarmotGeometryElement
Definition: MarmotGeometryElement.h:36
MarmotGeometryElement::assignNodeCoordinates
void assignNodeCoordinates(const double *coords)
Definition: MarmotGeometryElement.h:83
MarmotTypedefs.h
Marmot::FiniteElement::Hexa20
@ Hexa20
Definition: MarmotFiniteElement.h:46
MarmotGeometryElement::dNdX
dNdXiSized dNdX(const dNdXiSized &dNdXi, const JacobianSized &JacobianInverse) const
Definition: MarmotGeometryElement.h:106
MarmotGeometryElement::MarmotGeometryElement
MarmotGeometryElement()
Definition: MarmotGeometryElement.h:66
MarmotVoigt.h
Marmot::FiniteElement::Quad4
@ Quad4
Definition: MarmotFiniteElement.h:39
Marmot::FiniteElement::Spatial2D::nDim
constexpr int nDim
Definition: MarmotFiniteElement.h:112
MarmotGeometryElement::BGreen
BSized BGreen(const dNdXiSized &dNdX, const JacobianSized &F) const
MarmotGeometryElement::BSized
Eigen::Matrix< double, voigtSize, nNodes *nDim > BSized
Definition: MarmotGeometryElement.h:59
MarmotFiniteElement.h
MarmotGeometryElement::XiSized
Eigen::Matrix< double, nDim, 1 > XiSized
Definition: MarmotGeometryElement.h:53
Marmot::FiniteElement::Tetra4
@ Tetra4
Definition: MarmotFiniteElement.h:43
Marmot::FiniteElement::Hexa8
@ Hexa8
Definition: MarmotFiniteElement.h:45
MarmotGeometryElement::Jacobian
JacobianSized Jacobian(const dNdXiSized &dNdXi) const
Definition: MarmotGeometryElement.h:101
Marmot::FiniteElement::Bar2
@ Bar2
Definition: MarmotFiniteElement.h:37
Marmot
This file includes functions needed for calculations with stress and strain tensors written in voigt ...
Definition: MarmotTesting.h:30
Marmot::FiniteElement
Definition: MarmotDofLayoutTools.h:34
Marmot::FiniteElement::Tetra10
@ Tetra10
Definition: MarmotFiniteElement.h:44
MarmotGeometryElement::getElementShape
std::string getElementShape() const
Definition: MarmotGeometryElement.h:69
MarmotGeometryElement::dNdXiSized
Eigen::Matrix< double, nDim, nNodes > dNdXiSized
Definition: MarmotGeometryElement.h:58
MarmotGeometryElement::NSized
Eigen::Matrix< double, 1, nNodes > NSized
Definition: MarmotGeometryElement.h:56
MarmotGeometryElement::N
NSized N(const XiSized &xi) const
Marmot::FiniteElement::ElementShapes
ElementShapes
Definition: MarmotFiniteElement.h:36
MarmotGeometryElement::coordinates
Eigen::Map< const CoordinateVector > coordinates
Definition: MarmotGeometryElement.h:62
Marmot::ContinuumMechanics::VoigtNotation::VoigtSize
VoigtSize
Definition: MarmotVoigt.h:45
Marmot::FiniteElement::Spatial1D::Bar2::nNodes
constexpr int nNodes
Definition: MarmotFiniteElement.h:92
MarmotGeometryElement::shape
const Marmot::FiniteElement::ElementShapes shape
Definition: MarmotGeometryElement.h:63
MarmotGeometryElement::NB
NBSized NB(const NSized &N) const
Definition: MarmotGeometryElement.h:99
MarmotGeometryElement::F
JacobianSized F(const dNdXiSized &dNdX, const CoordinateVector &Q) const
Definition: MarmotGeometryElement.h:111
Marmot::FiniteElement::getElementShapeByMetric
ElementShapes getElementShapeByMetric(int nDim, int nNodes)
Definition: MarmotFiniteElementBasic.cpp:7
MarmotGeometryElement::NBSized
Eigen::Matrix< double, nDim, nNodes *nDim > NBSized
Definition: MarmotGeometryElement.h:57
MarmotGeometryElement::CoordinateVector
Eigen::Matrix< double, nDim *nNodes, 1 > CoordinateVector
Definition: MarmotGeometryElement.h:54
MarmotGeometryElement::voigtSize
static constexpr Marmot::ContinuumMechanics::VoigtNotation::VoigtSize voigtSize
Definition: MarmotGeometryElement.h:51
Marmot::FiniteElement::Quad8
@ Quad8
Definition: MarmotFiniteElement.h:40
MarmotGeometryElement::JacobianSized
Eigen::Matrix< double, nDim, nDim > JacobianSized
Definition: MarmotGeometryElement.h:55