OGS
NumLib::ShapeQuad9 Class Reference

Detailed Description

Shape function for a 9-nodes quadrilateral element

Definition at line 20 of file ShapeQuad9.h.

#include <ShapeQuad9.h>

Public Types

using MeshElement = MeshLib::Quad9
 

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 41 of file ShapeQuad9.h.

Member Function Documentation

◆ computeGradShapeFunction()

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

Evaluate derivatives of the shape function at the given point

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

Definition at line 28 of file ShapeQuad9-impl.h.

29{
30 dN[0] = (r[0] + 0.5) * r[1] * (r[1] + 1) / 2;
31 dN[1] = (r[0] - 0.5) * r[1] * (r[1] + 1) / 2;
32 dN[2] = (r[0] - 0.5) * r[1] * (r[1] - 1) / 2;
33 dN[3] = (r[0] + 0.5) * r[1] * (r[1] - 1) / 2;
34 dN[4] = -r[0] * r[1] * (1 + r[1]);
35 dN[5] = (1 - r[1] * r[1]) * (r[0] - 0.5);
36 dN[6] = r[0] * r[1] * (1 - r[1]);
37 dN[7] = (1 - r[1] * r[1]) * (r[0] + 0.5);
38 dN[8] = 2 * r[0] * (r[1] * r[1] - 1);
39
40 dN[10] = (r[1] + 0.5) * r[0] * (r[0] - 1) / 2;
41 dN[11] = (r[1] - 0.5) * r[0] * (r[0] - 1) / 2;
42 dN[12] = (r[1] - 0.5) * r[0] * (r[0] + 1) / 2;
43 dN[13] = (1 - r[0] * r[0]) * (r[1] + 0.5);
44 dN[14] = r[0] * r[1] * (1 - r[0]);
45 dN[15] = (1 - r[0] * r[0]) * (r[1] - 0.5);
46 dN[16] = -r[0] * r[1] * (1 + r[0]);
47 dN[17] = 2 * r[1] * (r[0] * r[0] - 1);
48 dN[9] = (r[1] + 0.5) * r[0] * (r[0] + 1) / 2;
49}

◆ computeShapeFunction()

template<class T_X , class T_N >
void NumLib::ShapeQuad9::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 ShapeQuad9-impl.h.

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

References NumLib::N.

Member Data Documentation

◆ DIM

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

Definition at line 42 of file ShapeQuad9.h.

◆ NPOINTS

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

Definition at line 43 of file ShapeQuad9.h.

◆ ORDER

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

Definition at line 44 of file ShapeQuad9.h.


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