PeriDEM 0.2.0
PeriDEM -- Peridynamics-based high-fidelity model for granular media
|
A class for mapping and quadrature related operations for bi-linear quadrangle element. More...
#include <quadElem.h>
Public Member Functions | |
QuadElem (size_t order) | |
Constructor for quadrangle element. | |
double | elemSize (const std::vector< util::Point > &nodes) override |
Returns the area of element. | |
std::vector< fe::QuadData > | getQuadDatas (const std::vector< util::Point > &nodes) override |
Get vector of quadrature data. | |
std::vector< fe::QuadData > | getQuadPoints (const std::vector< util::Point > &nodes) override |
Get vector of quadrature data. | |
Public Member Functions inherited from fe::BaseElem | |
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 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. | |
Private Member Functions | |
std::vector< double > | getShapes (const util::Point &p) override |
Returns the values of shape function at point p on reference element. | |
std::vector< std::vector< double > > | getDerShapes (const util::Point &p) override |
Returns the values of derivative of shape function at point p on reference element. | |
double | getJacobian (const util::Point &p, const std::vector< util::Point > &nodes, std::vector< std::vector< double > > *J) override |
Computes the Jacobian of map \( \Phi: T^0 \to T \). | |
void | init () override |
Compute the quadrature points for quadrangle element. | |
Additional Inherited Members | |
Protected Member Functions inherited from fe::BaseElem | |
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. | |
Protected Attributes inherited from fe::BaseElem | |
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 class for mapping and quadrature related operations for bi-linear quadrangle element.
The reference quadrangle element \( T^0 \) is given by vertices \( (-1, -1), \, (1,-1), \, (1,1), \, (-1,1) \).
\[N^0_1(\xi, \eta) = \frac{(1- \xi)(1 - \eta)}{4}, \quad N^0_2(\xi, \eta) = \frac{(1+ \xi)(1 - \eta)}{4}, \]
\[ N^0_3(\xi, \eta) = \frac{(1+ \xi)(1 + \eta)}{4},\quad N^0_4(\xi, \eta) = \frac{(1- \xi)(1 + \eta)}{4}. \]
\[\frac{d N^0_1(\xi, \eta)}{d\xi} = \frac{-(1 - \eta)}{4}, \quad \frac{d N^0_1(\xi, \eta)}{d\eta} = \frac{-(1 - \xi)}{4}, \]
\[\frac{d N^0_2(\xi, \eta)}{d\xi} = \frac{(1 - \eta)}{4}, \quad \frac{d N^0_2 (\xi, \eta)}{d\eta} = \frac{-(1 + \xi)}{4}, \]
\[\frac{d N^0_3(\xi, \eta)}{d\xi} = \frac{(1 + \eta)}{4}, \quad \frac{d N^0_3 (\xi, \eta)}{d\eta} = \frac{(1 + \xi)}{4}, \]
\[\frac{d N^0_4(\xi, \eta)}{d\xi} = \frac{-(1 + \eta)}{4}, \quad \frac{d N^0_4 (\xi, \eta)}{d\eta} = \frac{(1 - \xi)}{4}. \]
\[ x(\xi, \eta) = \sum_{i=1}^4 N^0_i(\xi, \eta) v^i_x, \quad y(\xi, \eta) = \sum_{i=1}^4 N^0_i(\xi, \eta) v^i_y \]
where \( v^1, v^2, v^3, v^4\) are vertices of element \( T \).\[ J = \left[ { \begin{array}{cc} \frac{dx}{d\xi} &\frac{dy}{d\xi} \\ \frac{dx}{d\eta} & \frac{dy}{d\eta} \\ \end{array} } \right] \]
and determinant of Jacobian is\[ det(J) = \frac{dx}{d\xi} \times \frac{dy}{d\eta} - \frac{dy}{d\xi}\times \frac{dx}{d\eta}. \]
Definition at line 64 of file quadElem.h.
|
explicit |
Constructor for quadrangle element.
order | Order of quadrature point approximation |
Definition at line 14 of file quadElem.cpp.
References init().
|
overridevirtual |
Returns the area of element.
If quadrangle \( T \) is given by points \( v^1, v^2, v^3, v^4 \) then the area is
\[ area(T) = \frac{(-v^1_1 + v^2_1 + v^3_1 - v^4_1) (-v^1_2 - v^2_2 + v^3_2 + v^4_2) - (-v^1_1 - v^2_1 + v^3_1 + v^4_1) (-v^1_2 + v^2_2 + v^3_2 - v^4_2)}{4}, \]
where \( v^i_1, v^i_2 \) are the x and y component of point \( v^i \).
Note that area and Jacobian of map \( \Phi: T^0 \to T \) are related as
\[ area(T) = area(T^0) \times det(J(\xi = 0, \eta = 0)), \]
where \( area(T^0) = 4 \) and \( J(\xi = 0, \eta = 0) \) is the Jacobian of map at point \( (0,0) \in T^0 \).
nodes | Vertices of element |
Implements fe::BaseElem.
Definition at line 21 of file quadElem.cpp.
|
overrideprivatevirtual |
Returns the values of derivative of shape function at point p on reference element.
p | Location of point |
Implements fe::BaseElem.
Definition at line 105 of file quadElem.cpp.
References util::Point::d_x, and util::Point::d_y.
|
overrideprivatevirtual |
Computes the Jacobian of map \( \Phi: T^0 \to T \).
p | Location of point in reference element |
nodes | Vertices of element |
J | Matrix to store the Jacobian (if not nullptr) |
Implements fe::BaseElem.
Definition at line 127 of file quadElem.cpp.
|
overridevirtual |
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 is set in the constructor. List of data
Let \( T \) is the given element with vertices \( v^1, v^2, v^3, v^4 \) and let \( T^0 \) is the reference element.
\[ w_q = w^0_q * det(J(\xi_q, \eta_q)) \]
where \( det(J) \) is the determinant of the Jacobian \( J \) of map \( \Phi: T^0 \to T \). \( J(\xi_q, \eta_q)\) is the value of \( J \) at point \( (\xi_q, \eta_q) \in T^0 \).\[ N_i(x(\xi_q, \eta_q), y(\xi_q, \eta_q)) = N^0_i(\xi_q, \eta_q). \]
nodes | Vector of vertices of an element |
Implements fe::BaseElem.
Definition at line 32 of file quadElem.cpp.
|
overridevirtual |
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 is set in the constructor. List of data
This function is lite version of QuadElem::getQuadDatas.
nodes | Vector of vertices of an element |
Implements fe::BaseElem.
Definition at line 72 of file quadElem.cpp.
|
overrideprivatevirtual |
Returns the values of shape function at point p on reference element.
p | Location of point |
Implements fe::BaseElem.
Definition at line 93 of file quadElem.cpp.
References util::Point::d_x, and util::Point::d_y.
|
overrideprivatevirtual |
Compute the quadrature points for quadrangle element.
Implements fe::BaseElem.
Definition at line 162 of file quadElem.cpp.
References fe::QuadData::d_derShapes, fe::QuadData::d_detJ, fe::QuadData::d_J, fe::QuadData::d_p, fe::QuadData::d_shapes, and fe::QuadData::d_w.
Referenced by QuadElem().