OGS
ShapeMatrices-impl.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4namespace NumLib
5{
6namespace detail
7{
8/*
9 * zero reset functions for Eigen
10 */
11template <class T>
12void setMatrixZero(T& mat)
13{
14 // mat.setZero();
15 const std::size_t n = mat.rows() * mat.cols();
16 auto* v = mat.data();
17 for (std::size_t i = 0; i < n; i++)
18 {
19 v[i] = .0;
20 }
21}
22
23template <class T>
24void setVectorZero(T& vec)
25{
26 // vec.setZero();
27 const std::size_t n = vec.size();
28 auto* v = vec.data();
29 for (std::size_t i = 0; i < n; i++)
30 {
31 v[i] = .0;
32 }
33}
34
35/*
36 * "tag dispatching" technique is used below to emulate explicit specialization
37 * of a template function in a template class
38 */
39template <ShapeMatrixType FIELD_TYPE>
41{
42};
43
44template <class T_N, class T_DNDR, class T_J, class T_DNDX>
50
51template <class T_N, class T_DNDR, class T_J, class T_DNDX>
57
58template <class T_N, class T_DNDR, class T_J, class T_DNDX>
67
68template <class T_N, class T_DNDR, class T_J, class T_DNDX>
75
76template <class T_N, class T_DNDR, class T_J, class T_DNDX>
84
85template <class T_N, class T_DNDR, class T_J, class T_DNDX>
92
93} // namespace detail
94
95template <class T_N, class T_DNDR, class T_J, class T_DNDX>
100
101template <class T_N, class T_DNDR, class T_J, class T_DNDX>
102template <ShapeMatrixType T_SHAPE_MATRIX_TYPE>
107
108template <class T_N, class T_DNDR, class T_J, class T_DNDX>
110{
111 out << "N :\n" << N << "\n";
112 out << "dNdr:\n" << dNdr << "\n";
113 out << "J :\n" << J << "\n";
114 out << "|J| : " << detJ << "\n";
115 out << "invJ:\n" << invJ << "\n";
116 out << "dNdx:\n" << dNdx << "\n";
117}
118
119template <class T_N, class T_DNDR, class T_J, class T_DNDX>
120std::ostream& operator<<(std::ostream& os,
122{
123 shape.write(os);
124 return os;
125}
126
127} // 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.
JacobianType invJ
Inverse matrix of the Jacobian.
ShapeType N
Vector of shape functions, N(r)
double detJ
Determinant of the Jacobian.
JacobianType J
Jacobian matrix, J=dx/dr.
void write(std::ostream &out) const
void setZero()
reset all data with zero