PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
|
A base class which provides methods to map points to/from reference element and to compute quadrature data. More...
#include <baseElem.h>
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::QuadData > | getQuadDatas (const std::vector< util::Point > &nodes)=0 |
Get vector of quadrature data. | |
virtual std::vector< fe::QuadData > | getQuadPoints (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::QuadData > | d_quads |
Quadrature data collection. | |
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
\[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\).\[ 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 \).Definition at line 84 of file baseElem.h.
fe::BaseElem::BaseElem | ( | size_t | order, |
size_t | element_type | ||
) |
Constructor.
order | Order of quadrature point approximation |
element_type | Type of element |
Definition at line 15 of file baseElem.cpp.
|
pure virtual |
Returns the size of element (length in 1-d, area in 2-d, volume in 3-d element)
nodes | Vertices of element |
Implemented in fe::LineElem, fe::QuadElem, fe::TetElem, and fe::TriElem.
Referenced by fe::Mesh::computeVol().
|
protectedpure virtual |
Returns the values of derivative of shape function at point p on reference element.
p | Location of point |
Implemented in fe::LineElem, fe::QuadElem, fe::TetElem, and fe::TriElem.
|
virtual |
Returns the values of derivative of shape function at point p.
p | Location of point |
nodes | Vertices of element |
Reimplemented in fe::LineElem, fe::TetElem, and fe::TriElem.
Definition at line 30 of file baseElem.cpp.
|
inline |
Get element type.
Definition at line 99 of file baseElem.h.
References d_elemType.
|
protectedpure virtual |
Computes Jacobian of map from reference element \( T^0 \) to given element \( T \).
p | Location of point in reference element |
nodes | Vertices of element |
J | Matrix to store the Jacobian |
Implemented in fe::LineElem, fe::QuadElem, fe::TetElem, and fe::TriElem.
|
inline |
Get number of quadrature points in the data.
Definition at line 111 of file baseElem.h.
References d_numQuadPts.
Referenced by fe::getCurrentQuadPoints(), fe::getMaxShearStressAndLoc(), fe::getStrainStress(), twoparticle_demo::Model::run(), and model::DEMModel::setupQuadratureData().
|
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:
nodes | Vector of vertices of an element |
Implemented in fe::LineElem, fe::QuadElem, fe::TetElem, and fe::TriElem.
Referenced by fe::Mesh::computeVol(), fe::getCurrentQuadPoints(), fe::getMaxShearStressAndLoc(), and fe::getStrainStress().
|
inline |
Get order of quadrature approximation.
Definition at line 105 of file baseElem.h.
References d_quadOrder.
|
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:
This function is a lite version of fe::BaseElem::getQuadDatas.
nodes | Vector of vertices of an element |
Implemented in fe::LineElem, fe::QuadElem, fe::TetElem, and fe::TriElem.
|
protectedpure virtual |
Returns the values of shape function at point p on reference element.
p | Location of point |
Implemented in fe::LineElem, fe::QuadElem, fe::TetElem, and fe::TriElem.
|
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.
p | Location of point |
nodes | Vertices of element |
Reimplemented in fe::LineElem, fe::TetElem, and fe::TriElem.
Definition at line 20 of file baseElem.cpp.
|
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.
|
protectedvirtual |
Maps point p in a given element to the reference element and returns the mapped point.
p | Location of point |
nodes | Vertices of element |
Reimplemented in fe::LineElem, fe::TetElem, and fe::TriElem.
Definition at line 40 of file baseElem.cpp.
|
protected |
|
protected |
Number of quadrature points for order d_quadOrder.
Definition at line 246 of file baseElem.h.
Referenced by getNumQuadPoints().
|
protected |
Order of quadrature point integration approximation.
Definition at line 243 of file baseElem.h.
Referenced by getQuadOrder(), and fe::TetElem::TetElem().
|
protected |
Quadrature data collection.
Definition at line 252 of file baseElem.h.