OGS
ShapeMatrices-impl.h
Go to the documentation of this file.
1 
13 namespace NumLib
14 {
15 namespace detail
16 {
17 /*
18  * zero reset functions for Eigen
19  */
20 template <class T>
21 void setMatrixZero(T& mat)
22 {
23  // mat.setZero();
24  const std::size_t n = mat.rows() * mat.cols();
25  auto* v = mat.data();
26  for (std::size_t i = 0; i < n; i++)
27  {
28  v[i] = .0;
29  }
30 }
31 
32 template <class T>
33 void setVectorZero(T& vec)
34 {
35  // vec.setZero();
36  const std::size_t n = vec.size();
37  auto* v = vec.data();
38  for (std::size_t i = 0; i < n; i++)
39  {
40  v[i] = .0;
41  }
42 }
43 
44 /*
45  * "tag dispatching" technique is used below to emulate explicit specialization
46  * of a template function in a template class
47  */
48 template <ShapeMatrixType FIELD_TYPE>
50 {
51 };
52 
53 template <class T_N, class T_DNDR, class T_J, class T_DNDX>
56 {
57  setVectorZero(shape.N);
58 }
59 
60 template <class T_N, class T_DNDR, class T_J, class T_DNDX>
63 {
64  setMatrixZero(shape.dNdr);
65 }
66 
67 template <class T_N, class T_DNDR, class T_J, class T_DNDX>
70 {
72  setMatrixZero(shape.J);
73  shape.detJ = .0;
74  shape.integralMeasure = 0.0;
75 }
76 
77 template <class T_N, class T_DNDR, class T_J, class T_DNDX>
80 {
83 }
84 
85 template <class T_N, class T_DNDR, class T_J, class T_DNDX>
88 {
90  setMatrixZero(shape.invJ);
91  setMatrixZero(shape.dNdx);
92 }
93 
94 template <class T_N, class T_DNDR, class T_J, class T_DNDX>
97 {
100 }
101 
102 } // namespace detail
103 
104 template <class T_N, class T_DNDR, class T_J, class T_DNDX>
106 {
108 }
109 
110 template <class T_N, class T_DNDR, class T_J, class T_DNDX>
111 template <ShapeMatrixType T_SHAPE_MATRIX_TYPE>
113 {
115 }
116 
117 template <class T_N, class T_DNDR, class T_J, class T_DNDX>
119 {
120  out << "N :\n" << N << "\n";
121  out << "dNdr:\n" << dNdr << "\n";
122  out << "J :\n" << J << "\n";
123  out << "|J| : " << detJ << "\n";
124  out << "invJ:\n" << invJ << "\n";
125  out << "dNdx:\n" << dNdx << "\n";
126 }
127 
128 template <class T_N, class T_DNDR, class T_J, class T_DNDX>
129 std::ostream& operator<<(std::ostream& os,
131 {
132  shape.write(os);
133  return os;
134 }
135 
136 } // namespace NumLib
void setMatrixZero(T &mat)
void setVectorZero(T &vec)
void setZero(ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX > &shape, ShapeDataFieldType< ShapeMatrixType::N >)
std::ostream & operator<<(std::ostream &os, LocalToGlobalIndexMap const &map)
Coordinates mapping matrices at particular location.
Definition: ShapeMatrices.h:43
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