OGS
NumLib::ShapeTet10 Class Reference

Detailed Description

Shape function for a 10-nodes tetrahedral element

Definition at line 20 of file ShapeTet10.h.

#include <ShapeTet10.h>

Public Types

using MeshElement = MeshLib::Tet10
 

Static Public Member Functions

template<class T_X , class T_N >
static void computeShapeFunction (const T_X &r, T_N &N)
 
template<class T_X , class T_N >
static void computeGradShapeFunction (const T_X &r, T_N &dN)
 

Static Public Attributes

static const unsigned DIM = MeshElement::dimension
 
static const unsigned NPOINTS = MeshElement::n_all_nodes
 
static constexpr int ORDER = 2
 

Member Typedef Documentation

◆ MeshElement

Definition at line 42 of file ShapeTet10.h.

Member Function Documentation

◆ computeGradShapeFunction()

template<class T_X , class T_N >
void NumLib::ShapeTet10::computeGradShapeFunction ( const T_X & r,
T_N & dN )
static

Evaluate derivatives of the shape function at the given point The point coordinates in r are not used.

Parameters
[in]rpoint coordinates
[out]dNa matrix of the derivatives

Definition at line 29 of file ShapeTet10-impl.h.

30{
31 dN[0] = 4.0 * (r[0] + r[1] + r[2]) - 3.0;
32 dN[1] = 4. * r[0] - 1.;
33 dN[2] = 0.0;
34 dN[3] = 0.0;
35 dN[4] = 4.0 * (1.0 - 2.0 * r[0] - r[1] - r[2]);
36 dN[5] = 4.0 * r[1];
37 dN[6] = -4.0 * r[1];
38 dN[7] = -4.0 * r[2];
39 dN[8] = 4.0 * r[2];
40 dN[9] = 0.0;
41
42 dN[10] = 4. * (r[0] + r[1] + r[2]) - 3.;
43 dN[11] = 0.0;
44 dN[12] = 4. * r[1] - 1.;
45 dN[13] = 0.;
46 dN[14] = -4.0 * r[0];
47 dN[15] = 4.0 * r[0];
48 dN[16] = 4.0 * (1.0 - r[0] - 2.0 * r[1] - r[2]);
49 dN[17] = -4.0 * r[2];
50 dN[18] = 0.0;
51 dN[19] = 4.0 * r[2];
52
53 dN[20] = 4. * (r[0] + r[1] + r[2]) - 3.;
54 dN[21] = 0.;
55 dN[22] = 0.;
56 dN[23] = 4. * r[2] - 1.;
57 dN[24] = -4.0 * r[0];
58 dN[25] = 0.0;
59 dN[26] = -4.0 * r[1];
60 dN[27] = 4.0 * (1.0 - r[0] - r[1] - 2.0 * r[2]);
61 dN[28] = 4.0 * r[0];
62 dN[29] = 4.0 * r[1];
63}

◆ computeShapeFunction()

template<class T_X , class T_N >
void NumLib::ShapeTet10::computeShapeFunction ( const T_X & r,
T_N & N )
static

Evaluate the shape function at the given point

Parameters
[in]rpoint coordinates
[out]Na vector of calculated shape function.

Definition at line 14 of file ShapeTet10-impl.h.

15{
16 N[0] = 2. * (1 - r[0] - r[1] - r[2]) * (0.5 - r[0] - r[1] - r[2]);
17 N[1] = r[0] * (2. * r[0] - 1);
18 N[2] = r[1] * (2. * r[1] - 1);
19 N[3] = r[2] * (2. * r[2] - 1);
20 N[4] = 4.0 * r[0] * (1.0 - r[0] - r[1] - r[2]);
21 N[5] = 4.0 * r[0] * r[1];
22 N[6] = 4.0 * r[1] * (1.0 - r[0] - r[1] - r[2]);
23 N[7] = 4.0 * r[2] * (1.0 - r[0] - r[1] - r[2]);
24 N[8] = 4.0 * r[0] * r[2];
25 N[9] = 4.0 * r[1] * r[2];
26}

References NumLib::N.

Member Data Documentation

◆ DIM

const unsigned NumLib::ShapeTet10::DIM = MeshElement::dimension
static

Definition at line 43 of file ShapeTet10.h.

◆ NPOINTS

const unsigned NumLib::ShapeTet10::NPOINTS = MeshElement::n_all_nodes
static

Definition at line 44 of file ShapeTet10.h.

◆ ORDER

constexpr int NumLib::ShapeTet10::ORDER = 2
staticconstexpr

Definition at line 45 of file ShapeTet10.h.


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