OGS
NumLib::ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX > Struct Template Reference

Detailed Description

template<class T_N, class T_DNDR, class T_J, class T_DNDX>
struct NumLib::ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX >

Coordinates mapping matrices at particular location.

Template Parameters
T_NVector type for shape functions
T_DNDRMatrix type for gradient of shape functions in natural coordinates
T_JJacobian matrix type
T_DNDXMatrix type for gradient of shape functions in physical coordinates

Definition at line 44 of file ShapeMatrices.h.

#include <ShapeMatrices.h>

Public Types

using ShapeType = T_N
 
using DrShapeType = T_DNDR
 
using JacobianType = T_J
 
using DxShapeType = T_DNDX
 

Public Member Functions

 ShapeMatrices ()=delete
 
 ShapeMatrices (std::size_t local_dim, std::size_t global_dim, std::size_t n_nodes)
 
void setZero ()
 reset all data with zero
 
template<ShapeMatrixType T_SHAPE_MATRIX_TYPE>
void setZero ()
 
void write (std::ostream &out) const
 

Public Attributes

ShapeType N
 Vector of shape functions, N(r)
 
DrShapeType dNdr
 
JacobianType J
 Jacobian matrix, J=dx/dr.
 
double detJ
 Determinant of the Jacobian.
 
JacobianType invJ
 Inverse matrix of the Jacobian.
 
DxShapeType dNdx
 
double integralMeasure
 
 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
 

Member Typedef Documentation

◆ DrShapeType

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
using NumLib::ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX >::DrShapeType = T_DNDR

Definition at line 47 of file ShapeMatrices.h.

◆ DxShapeType

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
using NumLib::ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX >::DxShapeType = T_DNDX

Definition at line 49 of file ShapeMatrices.h.

◆ JacobianType

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
using NumLib::ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX >::JacobianType = T_J

Definition at line 48 of file ShapeMatrices.h.

◆ ShapeType

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
using NumLib::ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX >::ShapeType = T_N

Definition at line 46 of file ShapeMatrices.h.

Constructor & Destructor Documentation

◆ ShapeMatrices() [1/2]

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
NumLib::ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX >::ShapeMatrices ( )
delete

Not default constructible, dimensions always must be given.

The default constructor has been deleted explicitly, because with dynamically allocated matrices it is rather easy to forget the required resize() call. Note: the resize() member is also deleted now.

◆ ShapeMatrices() [2/2]

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
NumLib::ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX >::ShapeMatrices ( std::size_t local_dim,
std::size_t global_dim,
std::size_t n_nodes )
inline

Initialize matrices and vectors

Parameters
local_dimSpatial dimension of the element e.g. 1 for line, 2 for quad, and 3 hex etc.
global_dimSpatial dimension of the element's exterior space e.g. 3 for a quad representing a boundary of a hex element.
n_nodesThe number of element nodes

Definition at line 79 of file ShapeMatrices.h.

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 }
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 setZero()
reset all data with zero

References NumLib::ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX >::setZero().

Member Function Documentation

◆ setZero() [1/2]

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
void NumLib::ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX >::setZero ( )
inline

reset all data with zero

Definition at line 105 of file ShapeMatrices-impl.h.

106{
107 detail::setZero(*this, detail::ShapeDataFieldType<ShapeMatrixType::ALL>());
108}
void setZero(ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX > &shape, ShapeDataFieldType< ShapeMatrixType::N >)

References NumLib::detail::setZero().

Referenced by NumLib::ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX >::ShapeMatrices().

◆ setZero() [2/2]

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
template<ShapeMatrixType T_SHAPE_MATRIX_TYPE>
void NumLib::ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX >::setZero ( )
inline

reset specified data with zero

Template Parameters
T_SHAPE_MATRIX_TYPEshape matrix types to be initialized

Definition at line 112 of file ShapeMatrices-impl.h.

113{
114 detail::setZero(*this, detail::ShapeDataFieldType<T_SHAPE_MATRIX_TYPE>());
115}

References NumLib::detail::setZero().

◆ write()

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
void NumLib::ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX >::write ( std::ostream & out) const

writes the matrix entries into the output stream

Parameters
outthe output stream

Definition at line 118 of file ShapeMatrices-impl.h.

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}

References NumLib::N.

Referenced by NumLib::operator<<().

Member Data Documentation

◆ detJ

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
double NumLib::ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX >::detJ

Determinant of the Jacobian.

Definition at line 55 of file ShapeMatrices.h.

Referenced by NumLib::detail::setZero().

◆ dNdr

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
DrShapeType NumLib::ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX >::dNdr

Matrix of gradient of shape functions in natural coordinates, dN(r)/dr

Definition at line 52 of file ShapeMatrices.h.

Referenced by NumLib::detail::setZero().

◆ dNdx

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
DxShapeType NumLib::ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX >::dNdx

Matrix of gradient of shape functions in physical coordinates, dN(r)/dx

Definition at line 57 of file ShapeMatrices.h.

Referenced by NumLib::detail::setZero().

◆ EIGEN_MAKE_ALIGNED_OPERATOR_NEW

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
NumLib::ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX >::EIGEN_MAKE_ALIGNED_OPERATOR_NEW

Definition at line 109 of file ShapeMatrices.h.

◆ integralMeasure

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
double NumLib::ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX >::integralMeasure

Definition at line 59 of file ShapeMatrices.h.

Referenced by NumLib::detail::setZero().

◆ invJ

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
JacobianType NumLib::ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX >::invJ

Inverse matrix of the Jacobian.

Definition at line 56 of file ShapeMatrices.h.

Referenced by NumLib::detail::setZero().

◆ J

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
JacobianType NumLib::ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX >::J

Jacobian matrix, J=dx/dr.

Definition at line 54 of file ShapeMatrices.h.

Referenced by NumLib::detail::setZero().

◆ N

template<class T_N , class T_DNDR , class T_J , class T_DNDX >
ShapeType NumLib::ShapeMatrices< T_N, T_DNDR, T_J, T_DNDX >::N

Vector of shape functions, N(r)

Definition at line 51 of file ShapeMatrices.h.

Referenced by NumLib::detail::setZero().


The documentation for this struct was generated from the following files: