NLMech  0.1.0
fe::QuadElem Class Reference

A class for mapping and quadrature related operations for bi-linear quadrangle element. More...

#include <quadElem.h>

Inheritance diagram for fe::QuadElem:
Inheritance graph
Collaboration diagram for fe::QuadElem:
Collaboration graph

Public Member Functions

double elemSize (const std::vector< util::Point3 > &nodes) override
 Returns the area of element. More...
 
std::vector< fe::QuadDatagetQuadDatas (const std::vector< util::Point3 > &nodes) override
 Get vector of quadrature data. More...
 
std::vector< fe::QuadDatagetQuadPoints (const std::vector< util::Point3 > &nodes) override
 Get vector of quadrature data. More...
 
void print (int nt=0, int lvl=0) const
 Prints the information about the instance of the object. More...
 
std::string printStr (int nt=0, int lvl=0) const
 Returns the string containing information about the instance of the object. More...
 
 QuadElem (size_t order)
 Constructor for quadrangle element. More...
 
- Public Member Functions inherited from fe::BaseElem
 BaseElem (size_t order, size_t element_type)
 Constructor. More...
 
virtual std::vector< std::vector< double > > getDerShapes (const util::Point3 &p, const std::vector< util::Point3 > &nodes)
 Returns the values of derivative of shape function at point p. More...
 
size_t getElemType ()
 Get element type. More...
 
size_t getNumQuadPoints ()
 Get number of quadrature points in the data. More...
 
size_t getQuadOrder ()
 Get order of quadrature approximation. More...
 
virtual std::vector< double > getShapes (const util::Point3 &p, const std::vector< util::Point3 > &nodes)
 Returns the values of shape function at point p. More...
 
void print (int nt=0, int lvl=0) const
 Prints the information about the instance of the object. More...
 
std::string printStr (int nt=0, int lvl=0) const
 Returns the string containing information about the instance of the object. More...
 

Private Member Functions

std::vector< std::vector< double > > getDerShapes (const util::Point3 &p) override
 Returns the values of derivative of shape function at point p on reference element. More...
 
double getJacobian (const util::Point3 &p, const std::vector< util::Point3 > &nodes, std::vector< std::vector< double >> *J) override
 Computes the Jacobian of map $ \Phi: T^0 \to T $. More...
 
std::vector< double > getShapes (const util::Point3 &p) override
 Returns the values of shape function at point p on reference element. More...
 
void init () override
 Compute the quadrature points for quadrangle element.
 

Detailed Description

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) $.

  1. The shape functions at point $ (\xi, \eta ) \in T^0 $ are

    \[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}. \]

  2. Derivative of shape functions at point $ (\xi, \eta ) \in T^0 $ are as follows

    \[\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}. \]

  3. Map $ \Phi: T^0 \to T$ is given by

    \[ 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 $.
  4. Jacobian of the map $ \Phi: T^0 \to T$ is given by

    \[ 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}. \]

Constructor & Destructor Documentation

◆ QuadElem()

fe::QuadElem::QuadElem ( size_t  order)
explicit

Constructor for quadrangle element.

Parameters
orderOrder of quadrature point approximation

Member Function Documentation

◆ elemSize()

double fe::QuadElem::elemSize ( const std::vector< util::Point3 > &  nodes)
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 $.

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

Implements fe::BaseElem.

◆ getDerShapes()

std::vector< std::vector< double > > fe::QuadElem::getDerShapes ( const util::Point3 p)
overrideprivatevirtual

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

Implements fe::BaseElem.

◆ getJacobian()

double fe::QuadElem::getJacobian ( const util::Point3 p,
const std::vector< util::Point3 > &  nodes,
std::vector< std::vector< double >> *  J 
)
overrideprivatevirtual

Computes the Jacobian of map $ \Phi: T^0 \to T $.

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

Implements fe::BaseElem.

◆ getQuadDatas()

std::vector< fe::QuadData > fe::QuadElem::getQuadDatas ( const std::vector< util::Point3 > &  nodes)
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

  • quad point
  • quad weight
  • shape function evaluated at quad point
  • derivative of shape function evaluated at quad point
  • Jacobian matrix
  • determinant of the Jacobian

Let $ T $ is the given element with vertices $ v^1, v^2, v^3, v^4 $ and let $ T^0 $ is the reference element.

  1. To compute quadrature point, we first compute the quadrature points on $ T^0 $, and then map the point on $ T^0 $ to the point on $ T $ using the map $ \Phi: T^0 \to T$.
  2. To compute the quadrature weight, we compute the quadrature weight associated to the quadrature point in reference triangle $ T^0 $. Suppose $ w^0_q $ is the quadrature weight associated to quadrature point $ (\xi_q, \eta_q) \in T^0 $, then the quadrature point $ w_q $ associated to the mapped point $ (x(\xi_q, \eta_q), y(\xi_q, \eta_q)) \in T $ is given by

    \[ 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 $.
  3. We compute shape functions $ N_1, N_2, N_3, N_4$ associated to $ T $ at the quadrature point $ (x(\xi_q, \eta_q), y(\xi_q,\eta_q)) $ using formula

    \[ N_i(x(\xi_q, \eta_q), y(\xi_q, \eta_q)) = N^0_i(\xi_q, \eta_q). \]

  4. To compute the derivatives of shape functions such as $ \frac{\partial N_i}{\partial x}, \frac{\partial N_i}{\partial y}$ associated to $ T $, we use the relation between derivatives of shape function in $ T $ and $ T^0 $ described in fe::TriElem::getDerShapes.
Parameters
nodesVector of vertices of an element
Returns
vector Vector of QuadData

Implements fe::BaseElem.

◆ getQuadPoints()

std::vector< fe::QuadData > fe::QuadElem::getQuadPoints ( const std::vector< util::Point3 > &  nodes)
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

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

This function is lite version of QuadElem::getQuadDatas.

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

Implements fe::BaseElem.

◆ getShapes()

std::vector< double > fe::QuadElem::getShapes ( const util::Point3 p)
overrideprivatevirtual

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

Implements fe::BaseElem.

◆ print()

void fe::QuadElem::print ( int  nt = 0,
int  lvl = 0 
) const
inline

Prints the information about the instance of the object.

Parameters
ntNumber of tabs to append before each line of string
lvlLevel of information sought (higher level means more information)
Here is the call graph for this function:

◆ printStr()

std::string fe::QuadElem::printStr ( int  nt = 0,
int  lvl = 0 
) const

Returns the string containing information about the instance of the object.

Parameters
ntNumber of tabs to append before each line of string
lvlLevel of information sought (higher level means more information)
Returns
string String containing information about this object

Referenced by print().

Here is the call graph for this function:
Here is the caller graph for this function:

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