OGS
ShapeMatrices.h
Go to the documentation of this file.
1 
13 #pragma once
14 
15 #include <Eigen/Core>
16 #include <ostream>
17 
18 namespace NumLib
19 {
23 enum class ShapeMatrixType
24 {
25  N,
26  DNDR,
27  N_J,
28  DNDR_J,
29  DNDX,
30  ALL
31 };
32 
41 template <class T_N, class T_DNDR, class T_J, class T_DNDX>
43 {
44  using ShapeType = T_N;
45  using DrShapeType = T_DNDR;
46  using JacobianType = T_J;
47  using DxShapeType = T_DNDX;
48 
53  double detJ;
58 
66  ShapeMatrices() = delete;
67 
77  ShapeMatrices(std::size_t local_dim, std::size_t global_dim,
78  std::size_t n_nodes)
79  : N(n_nodes),
80  dNdr(local_dim, n_nodes),
81  J(local_dim, local_dim),
82  detJ(.0),
83  invJ(local_dim, local_dim),
84  dNdx(global_dim, n_nodes),
85  integralMeasure(0.0)
86  {
87  setZero();
88  }
89 
91  void setZero();
92 
98  template <ShapeMatrixType T_SHAPE_MATRIX_TYPE>
99  void setZero();
100 
105  void write(std::ostream& out) const;
106 
108 }; // ShapeMatrices
109 
110 } // namespace NumLib
111 
112 #include "ShapeMatrices-impl.h"
ShapeMatrixType
Shape matrix type to be calculated.
Definition: ShapeMatrices.h:24
@ DNDR
calculates dNdr
@ ALL
calculates all
@ N_J
calculates N, dNdr, J, and detJ
@ DNDR_J
calculates dNdr, J, and detJ
@ DNDX
calculates dNdr, J, detJ, invJ, and dNdx
Coordinates mapping matrices at particular location.
Definition: ShapeMatrices.h:43
ShapeMatrices(std::size_t local_dim, std::size_t global_dim, std::size_t n_nodes)
Definition: ShapeMatrices.h:77
JacobianType invJ
Inverse matrix of the Jacobian.
Definition: ShapeMatrices.h:54
ShapeType N
Vector of shape functions, N(r)
Definition: ShapeMatrices.h:49
double detJ
Determinant of the Jacobian.
Definition: ShapeMatrices.h:53
JacobianType J
Jacobian matrix, J=dx/dr.
Definition: ShapeMatrices.h:52
void write(std::ostream &out) const
void setZero()
reset all data with zero