OGS
NumLib::ShapeTet10 Class Reference

Detailed Description

Shape function for a 10-nodes tetrahedral element

Definition at line 14 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 constexpr std::array reference_element_centre = {0.25, 0.25, 0.25}
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 38 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 22 of file ShapeTet10-impl.h.

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

◆ 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 7 of file ShapeTet10-impl.h.

8{
9 N[0] = 2. * (1 - r[0] - r[1] - r[2]) * (0.5 - r[0] - r[1] - r[2]);
10 N[1] = r[0] * (2. * r[0] - 1);
11 N[2] = r[1] * (2. * r[1] - 1);
12 N[3] = r[2] * (2. * r[2] - 1);
13 N[4] = 4.0 * r[0] * (1.0 - r[0] - r[1] - r[2]);
14 N[5] = 4.0 * r[0] * r[1];
15 N[6] = 4.0 * r[1] * (1.0 - r[0] - r[1] - r[2]);
16 N[7] = 4.0 * r[2] * (1.0 - r[0] - r[1] - r[2]);
17 N[8] = 4.0 * r[0] * r[2];
18 N[9] = 4.0 * r[1] * r[2];
19}

References NumLib::N.

Member Data Documentation

◆ DIM

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

Definition at line 39 of file ShapeTet10.h.

◆ NPOINTS

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

Definition at line 40 of file ShapeTet10.h.

◆ ORDER

int NumLib::ShapeTet10::ORDER = 2
staticconstexpr

Definition at line 41 of file ShapeTet10.h.

◆ reference_element_centre

std::array NumLib::ShapeTet10::reference_element_centre = {0.25, 0.25, 0.25}
staticconstexpr

Definition at line 36 of file ShapeTet10.h.

36{0.25, 0.25, 0.25};

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