PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
Loading...
Searching...
No Matches
fe::BaseElem Class Referenceabstract

A base class which provides methods to map points to/from reference element and to compute quadrature data. More...

#include <baseElem.h>

Inheritance diagram for fe::BaseElem:
Collaboration diagram for fe::BaseElem:

Public Member Functions

 BaseElem (size_t order, size_t element_type)
 Constructor.
 
size_t getElemType ()
 Get element type.
 
size_t getQuadOrder ()
 Get order of quadrature approximation.
 
size_t getNumQuadPoints ()
 Get number of quadrature points in the data.
 
virtual double elemSize (const std::vector< util::Point > &nodes)=0
 Returns the size of element (length in 1-d, area in 2-d, volume in 3-d element)
 
virtual std::vector< double > getShapes (const util::Point &p, const std::vector< util::Point > &nodes)
 Returns the values of shape function at point p.
 
virtual std::vector< std::vector< double > > getDerShapes (const util::Point &p, const std::vector< util::Point > &nodes)
 Returns the values of derivative of shape function at point p.
 
virtual std::vector< fe::QuadDatagetQuadDatas (const std::vector< util::Point > &nodes)=0
 Get vector of quadrature data.
 
virtual std::vector< fe::QuadDatagetQuadPoints (const std::vector< util::Point > &nodes)=0
 Get vector of quadrature data.
 

Protected Member Functions

virtual std::vector< double > getShapes (const util::Point &p)=0
 Returns the values of shape function at point p on reference element.
 
virtual std::vector< std::vector< double > > getDerShapes (const util::Point &p)=0
 Returns the values of derivative of shape function at point p on reference element.
 
virtual util::Point mapPointToRefElem (const util::Point &p, const std::vector< util::Point > &nodes)
 Maps point p in a given element to the reference element and returns the mapped point.
 
virtual double getJacobian (const util::Point &p, const std::vector< util::Point > &nodes, std::vector< std::vector< double > > *J)=0
 Computes Jacobian of map from reference element \( T^0 \) to given element \( T \).
 
virtual void init ()=0
 Compute the quadrature points.
 

Protected Attributes

size_t d_quadOrder
 Order of quadrature point integration approximation.
 
size_t d_numQuadPts
 Number of quadrature points for order d_quadOrder.
 
size_t d_elemType
 Element type.
 
std::vector< fe::QuadDatad_quads
 Quadrature data collection.
 

Detailed Description

A base class which provides methods to map points to/from reference element and to compute quadrature data.

For any type of element, such as fe::LineElem, fe::TriElem, fe::QuadElem, we have following important points

  1. A simple element with predefined vertices is considered as a reference element \( T^0 \). E.g., reference element in fe::TriElem is a triangle with vertices at \((0,0), (1,0), (0,1) \).
  2. Points in reference element \( T^0 \) are described by \( (\xi, \eta, \zeta)\) (in 3-d), \( (\xi, \eta) \) (in 2-d), and \( \xi\) (in 1-d). Points in any other element \( T \) are described by \( (x,y,z)\) (in 3-d), \( (x,y) \) (in 2-d), and \( x\) (in 1-d). Here, we restrict discussion to 3-d element. For 2-d, one should ignore coordinate \(\zeta \) and \( z\), and for 1-d, one should ignore \( \eta, \zeta \) and \( y,z\).
  3. Associated to vertices of the element \( T^0 \) we have shape functions \( N^0_1, N^0_2, ..., N^0_n \), where \( n \) is the number of vertices in the element. Shape functions \( N^0_i \) are functions of point \( (\xi, \eta, \zeta) \in T^0 \). For each type of element, \( n \) is fixed.
  4. For any element \( T \), shape functions are denoted as \( N_1, N_2, ..., N_n \) and are functions of point \( (x,y,z) \in T \).
  5. Element \( T \) is described by vertices \( v^1, v^2, ..., v^n \).
  6. A map \(\Phi : T^0 \to T \) , where \( T \) is a given element formed by vertices \( v^1, v^2, ..., v^n \), is defined as follows:

    \[x(\xi, \eta, \zeta) = \sum_{i=1}^n N^0_i(\xi, \eta, \zeta) v^i_x, \quad y(\xi, \eta, \zeta) = \sum_{i=1}^n N^0_i(\xi, \eta, \zeta) v^i_y, \quad z (\xi, \eta, \zeta) = \sum_{i=1}^n N^0_i(\xi, \eta, \zeta) v^i_z \]

    where \( v^i_x, v^i_y, v^i_z \) are the x, y, and z component of point \( v^i\).
  7. Jacobian of map \( \Phi : T^0 \to T \) is given by

    \[ J = \left[ { \begin{array}{ccc} \frac{dx}{d\xi} &\frac{dy}{d\xi} & \frac{dz}{d\xi} \\ \frac{dx}{d\eta} & \frac{dy}{d\eta} & \frac{dz}{d\eta} \\ \frac{dx}{d\zeta} & \frac{dy}{d\zeta} & \frac{dz}{d\zeta} \\ \end{array} } \right]. \]

    For 2-d element it is

    \[ J = \left[ { \begin{array}{cc} \frac{dx}{d\xi} &\frac{dy}{d\xi} \\ \frac{dx}{d\eta} & \frac{dy}{d\eta} \\ \end{array} } \right]. \]

    For 1-d element it is

    \[ J = \frac{dx}{d\xi}. \]

    Determinant of the Jacobian is an important quantity and is used to compute the quadrature points, and inverse map \( \Phi^{-1} : T \to T^0 \).
See also
fe::LineElem, fe::TriElem, fe::QuadElem

Definition at line 84 of file baseElem.h.

Constructor & Destructor Documentation

◆ BaseElem()

fe::BaseElem::BaseElem ( size_t  order,
size_t  element_type 
)

Constructor.

Parameters
orderOrder of quadrature point approximation
element_typeType of element

Definition at line 15 of file baseElem.cpp.

16 : d_quadOrder(order), d_elemType(element_type),
size_t d_numQuadPts
Number of quadrature points for order d_quadOrder.
Definition baseElem.h:246
size_t d_elemType
Element type.
Definition baseElem.h:249
size_t d_quadOrder
Order of quadrature point integration approximation.
Definition baseElem.h:243
static int vtk_map_element_to_num_nodes[16]
Map from element type to number of nodes (for vtk)

Member Function Documentation

◆ elemSize()

virtual double fe::BaseElem::elemSize ( const std::vector< util::Point > &  nodes)
pure virtual

Returns the size of element (length in 1-d, area in 2-d, volume in 3-d element)

Parameters
nodesVertices of element
Returns
size Size of the element

Implemented in fe::LineElem, fe::QuadElem, fe::TetElem, and fe::TriElem.

Referenced by fe::Mesh::computeVol().

Here is the caller graph for this function:

◆ getDerShapes() [1/2]

virtual std::vector< std::vector< double > > fe::BaseElem::getDerShapes ( const util::Point p)
protectedpure virtual

Returns the values of derivative of shape function at point p on reference element.

Parameters
pLocation of point
Returns
vector Vector of derivative of shape functions

Implemented in fe::LineElem, fe::QuadElem, fe::TetElem, and fe::TriElem.

◆ getDerShapes() [2/2]

std::vector< std::vector< double > > fe::BaseElem::getDerShapes ( const util::Point p,
const std::vector< util::Point > &  nodes 
)
virtual

Returns the values of derivative of shape function at point p.

Parameters
pLocation of point
nodesVertices of element
Returns
vector Vector of derivative of shape functions

Reimplemented in fe::LineElem, fe::TetElem, and fe::TriElem.

Definition at line 30 of file baseElem.cpp.

31 {
32 std::cerr << "Error: For element type = " << d_elemType << " the map from "
33 << "element to reference element is not available.\n"
34 << "Therefore, derivatives of shape function at any "
35 "arbitrary point in the element can not be computed.\n";
36 exit(1);
37}

◆ getElemType()

size_t fe::BaseElem::getElemType ( )
inline

Get element type.

Returns
type Type of element

Definition at line 99 of file baseElem.h.

99{ return d_elemType; }

References d_elemType.

◆ getJacobian()

virtual double fe::BaseElem::getJacobian ( const util::Point p,
const std::vector< util::Point > &  nodes,
std::vector< std::vector< double > > *  J 
)
protectedpure virtual

Computes Jacobian of map from reference element \( T^0 \) to given element \( T \).

Parameters
pLocation of point in reference element
nodesVertices of element
JMatrix to store the Jacobian
Returns
det(J) Determinant of the Jacobain

Implemented in fe::LineElem, fe::QuadElem, fe::TetElem, and fe::TriElem.

◆ getNumQuadPoints()

size_t fe::BaseElem::getNumQuadPoints ( )
inline

Get number of quadrature points in the data.

Returns
N Number of quadrature points

Definition at line 111 of file baseElem.h.

111{ return d_numQuadPts; }

References d_numQuadPts.

Referenced by fe::getCurrentQuadPoints(), fe::getMaxShearStressAndLoc(), fe::getStrainStress(), twoparticle_demo::Model::run(), and model::DEMModel::setupQuadratureData().

Here is the caller graph for this function:

◆ getQuadDatas()

virtual std::vector< fe::QuadData > fe::BaseElem::getQuadDatas ( const std::vector< util::Point > &  nodes)
pure virtual

Get vector of quadrature data.

Given element vertices, this method returns the list of quadrature point and associated quantities. Here, order of quadrature approximation and element type are set in the constructor. List of data for each quad point:

  • quad point
  • quad weight
  • shape function evaluated at quad point
  • derivative of shape function evaluated at quad point
  • Jacobian matrix
  • determinant of the Jacobian
Parameters
nodesVector of vertices of an element
Returns
vector Vector of QuadData

Implemented in fe::LineElem, fe::QuadElem, fe::TetElem, and fe::TriElem.

Referenced by fe::Mesh::computeVol(), fe::getCurrentQuadPoints(), fe::getMaxShearStressAndLoc(), and fe::getStrainStress().

Here is the caller graph for this function:

◆ getQuadOrder()

size_t fe::BaseElem::getQuadOrder ( )
inline

Get order of quadrature approximation.

Returns
order Order of approximation

Definition at line 105 of file baseElem.h.

105{ return d_quadOrder; }

References d_quadOrder.

◆ getQuadPoints()

virtual std::vector< fe::QuadData > fe::BaseElem::getQuadPoints ( const std::vector< util::Point > &  nodes)
pure virtual

Get vector of quadrature data.

Given element vertices, this method returns the list of quadrature point and essential quantities at quadrature points. Here, order of quadrature approximation and element type are set in the constructor. List of data for each quad point:

  • quad point
  • quad weight
  • shape function evaluated at quad point

This function is a lite version of fe::BaseElem::getQuadDatas.

Parameters
nodesVector of vertices of an element
Returns
vector Vector of QuadData

Implemented in fe::LineElem, fe::QuadElem, fe::TetElem, and fe::TriElem.

◆ getShapes() [1/2]

virtual std::vector< double > fe::BaseElem::getShapes ( const util::Point p)
protectedpure virtual

Returns the values of shape function at point p on reference element.

Parameters
pLocation of point
Returns
vector Vector of shape functions at point p

Implemented in fe::LineElem, fe::QuadElem, fe::TetElem, and fe::TriElem.

◆ getShapes() [2/2]

std::vector< double > fe::BaseElem::getShapes ( const util::Point p,
const std::vector< util::Point > &  nodes 
)
virtual

Returns the values of shape function at point p.

The point p is assumed to be inside an element \( T\) given by nodes. The idea is to first map point p to the point \( p^0\) in reference element \( T^0 \) and then compute shape functions at \( p^0 \).

The map from any element \( T \) to reference element is \( T^0 \) gets complex for elements like fe::QuadElem. For elements fe:LineElem and fe::TriElem, this is simple to compute. Therefore, this method is only implemented for fe::TriElem and fe::LineElem.

Parameters
pLocation of point
nodesVertices of element
Returns
vector Vector of shape functions at point p

Reimplemented in fe::LineElem, fe::TetElem, and fe::TriElem.

Definition at line 20 of file baseElem.cpp.

21 {
22 std::cerr << "Error: For element type = " << d_elemType << " the map from "
23 << "element to reference element is not available.\n"
24 << "Therefore, shape function evaluation at any arbitrary point "
25 "in the element is not possible.\n";
26 exit(1);
27}

◆ init()

void fe::BaseElem::init ( )
protectedpure virtual

Compute the quadrature points.

This must be implemented by inheriting classes.

Implemented in fe::LineElem, fe::QuadElem, fe::TetElem, and fe::TriElem.

Definition at line 47 of file baseElem.cpp.

47 {
48
49 std::cerr << "Error: init() of BaseElem must be implemented in inheriting "
50 "class.\n";
51 exit(1);
52}

◆ mapPointToRefElem()

util::Point fe::BaseElem::mapPointToRefElem ( const util::Point p,
const std::vector< util::Point > &  nodes 
)
protectedvirtual

Maps point p in a given element to the reference element and returns the mapped point.

Parameters
pLocation of point
nodesVertices of element
Returns
vector Vector of shape functions at point p

Reimplemented in fe::LineElem, fe::TetElem, and fe::TriElem.

Definition at line 40 of file baseElem.cpp.

41 {
42 std::cerr << "Error: For element type = " << d_elemType << " the map from "
43 << "element to reference element is not available.\n";
44 exit(1);
45}

Field Documentation

◆ d_elemType

size_t fe::BaseElem::d_elemType
protected

Element type.

Definition at line 249 of file baseElem.h.

Referenced by getElemType().

◆ d_numQuadPts

size_t fe::BaseElem::d_numQuadPts
protected

Number of quadrature points for order d_quadOrder.

Definition at line 246 of file baseElem.h.

Referenced by getNumQuadPoints().

◆ d_quadOrder

size_t fe::BaseElem::d_quadOrder
protected

Order of quadrature point integration approximation.

Definition at line 243 of file baseElem.h.

Referenced by getQuadOrder(), and fe::TetElem::TetElem().

◆ d_quads

std::vector<fe::QuadData> fe::BaseElem::d_quads
protected

Quadrature data collection.

Definition at line 252 of file baseElem.h.


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