OGS
ShapeMatrices.h
Go to the documentation of this file.
1
12
13#pragma once
14
15#include <Eigen/Core>
16#include <iosfwd>
17
18namespace NumLib
19{
32
43template <class T_N, class T_DNDR, class T_J, class T_DNDX>
45{
46 using ShapeType = T_N;
47 using DrShapeType = T_DNDR;
48 using JacobianType = T_J;
49 using DxShapeType = T_DNDX;
50
55 double detJ;
60
68 ShapeMatrices() = delete;
69
79 ShapeMatrices(std::size_t local_dim, std::size_t global_dim,
80 std::size_t n_nodes)
81 : N(n_nodes),
82 dNdr(local_dim, n_nodes),
83 J(local_dim, local_dim),
84 detJ(.0),
85 invJ(local_dim, local_dim),
86 dNdx(global_dim, n_nodes),
88 {
89 setZero();
90 }
91
93 void setZero();
94
100 template <ShapeMatrixType T_SHAPE_MATRIX_TYPE>
101 void setZero();
102
107 void write(std::ostream& out) const;
108
110}; // ShapeMatrices
111
112} // namespace NumLib
113
114#include "ShapeMatrices-impl.h"
ShapeMatrixType
Shape matrix type to be calculated.
@ DNDR
calculates dNdr
@ N_J
calculates N, dNdr, J, and detJ
@ DNDR_J
calculates dNdr, J, and detJ
@ DNDX
calculates dNdr, J, detJ, invJ, and dNdx
ShapeMatrices(std::size_t local_dim, std::size_t global_dim, std::size_t n_nodes)
void write(std::ostream &out) const
void setZero()
reset all data with zero