NLMech  0.1.0
fe::LineElem Class Reference

A class for mapping and quadrature related operations for linear 2-node line element. More...

#include <lineElem.h>

Inheritance diagram for fe::LineElem:
Inheritance graph
Collaboration diagram for fe::LineElem:
Collaboration graph

Public Member Functions

double elemSize (const std::vector< util::Point3 > &nodes) override
 Returns the length of element. More...
 
std::vector< std::vector< double > > getDerShapes (const util::Point3 &p, const std::vector< util::Point3 > &nodes) override
 Returns the values of derivative of shape function at point p. 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...
 
std::vector< double > getShapes (const util::Point3 &p, const std::vector< util::Point3 > &nodes) override
 Returns the values of shape function at point p. More...
 
 LineElem (size_t order)
 Constructor for line element. 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...
 
- Public Member Functions inherited from fe::BaseElem
 BaseElem (size_t order, size_t element_type)
 Constructor. 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...
 
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 Jacobian of the 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 line element.
 
util::Point3 mapPointToRefElem (const util::Point3 &p, const std::vector< util::Point3 > &nodes) override
 Maps point p in a given element to the reference element. More...
 

Detailed Description

A class for mapping and quadrature related operations for linear 2-node line element.

The reference line element $T^0 $ is given by vertices at $ v^1 = -1, v^2 = 1 $.

  1. The shape functions at point $ \xi \in T^0 $ are

    \[N^0_1(\xi) = \frac{1 - \xi}{2}, \quad N^0_2(\xi) = \frac{1 + \xi}{2}. \]

  2. Derivative of shape functions are constant and are as follows

    \[\frac{d N^0_1(\xi)}{d\xi} = \frac{-1}{2}, \, \frac{d N^0_2(\xi)}{d\xi} = \frac{1}{2}. \]

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

    \[ x(\xi) = \sum_{i=1}^2 N^0_i(\xi) v^i_x, \]

    where $ v^1, v^2$ are vertices of element $ T $. For 1-d points, we simply have $ v^i_x = v^i $.
  4. The Jacobian of the map $ \Phi: T^0 \to T $ is given by

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

    Since it is a 1-d case, Jacobian and its determinant are same. For line element the formula for Jacobian is as follows

    \[ J = \frac{dx}{d\xi} = \frac{v^2_x - v^1_x}{2} = \frac{length(T) }{length(T^0)}. \]

  5. Inverse map $ \Phi^{-1} $ from $ x \in T $ to $ \xi \in T^0$ is given by

    \[ \xi(x) = \frac{2}{v^2_x - v^1_x} (x - \frac{v^2_x + v^1_x}{2}) = \frac{1}{J}(x - \frac{v^2_x + v^1_x}{2}). \]

Constructor & Destructor Documentation

◆ LineElem()

fe::LineElem::LineElem ( size_t  order)
explicit

Constructor for line element.

Parameters
orderOrder of quadrature point approximation

Member Function Documentation

◆ elemSize()

double fe::LineElem::elemSize ( const std::vector< util::Point3 > &  nodes)
overridevirtual

Returns the length of element.

If line $ T $ is given by points $ v^1, v^2$ then the length is simply

\[ length(T) = v^2_x - v^1_x. \]

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

Implements fe::BaseElem.

◆ getDerShapes() [1/2]

std::vector< std::vector< double > > fe::LineElem::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.

◆ getDerShapes() [2/2]

std::vector< std::vector< double > > fe::LineElem::getDerShapes ( const util::Point3 p,
const std::vector< util::Point3 > &  nodes 
)
overridevirtual

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

Let $ x$ is the point on line $ T $ and let $ \xi $ is the point on reference line $ T^0 $. Let shape functions on $ T$ are $ N_1, N_2 $ and on $ T^0 $ are $ N^0_1, N^0_2 $.

We are interested in $ \frac{\partial N_i(x_p)}{\partial x}$. By using the map $ \xi \to x $ we have

\[ N^0_i(\xi) = N_i(x(\xi)) \]

and therefore we can write

\[ \frac{\partial N^0_i(\xi)}{\partial \xi} = \frac{\partial N_i}{\partial x} \frac{\partial x}{\partial \xi}. \]

Since $ \frac{\partial x}{\partial \xi} = J $ is a Jacobian of map which can be computed easily if $ v^1, v^2 $ are known, we can invert and obtain the formula

\[ \frac{\partial N_i(\xi)}{\partial x} = \frac{1}{J} \frac{\partial N^0_i(\xi)}{\partial \xi}. \]

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

Reimplemented from fe::BaseElem.

◆ getJacobian()

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

Computes Jacobian of the 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 Jacobian (same as Jacobian in 1-d)

Implements fe::BaseElem.

◆ getQuadDatas()

std::vector< fe::QuadData > fe::LineElem::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 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

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

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

    \[ w_q = w^0_q * J \]

    where $ J $ is the Jacobian of map $ \Phi $.
  3. We compute shape functions $ N_1, N_2$ associated to $ T $ at quadrature point $ x(\xi_q) $ using formula

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

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

Implements fe::BaseElem.

◆ getQuadPoints()

std::vector< fe::QuadData > fe::LineElem::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 for each quad point:

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

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

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

Implements fe::BaseElem.

◆ getShapes() [1/2]

std::vector< double > fe::LineElem::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.

◆ getShapes() [2/2]

std::vector< double > fe::LineElem::getShapes ( const util::Point3 p,
const std::vector< util::Point3 > &  nodes 
)
overridevirtual

Returns the values of shape function at point p.

Line $ T $ is given by points $ v^1, v^2$. We first map the point p in $ T $ to reference line $ T^0 $ using fe::LineElem::mapPointToRefElem and then compute shape functions at the mapped point using fe::LineElem::getShapes(const util::Point3 &).

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

Reimplemented from fe::BaseElem.

◆ mapPointToRefElem()

util::Point3 fe::LineElem::mapPointToRefElem ( const util::Point3 p,
const std::vector< util::Point3 > &  nodes 
)
overrideprivatevirtual

Maps point p in a given element to the reference element.

Let $ v^1, v^2$ are vertices of element $ T$ and let $T^0 $ is the reference element. Map $\Phi : T^0 \to T$ is given by

\[ \xi(x) = \frac{2}{v^2_x - v^1_x} (x - \frac{v^2_x + v^1_x}{2}) = \frac{1}{J}(x - \frac{v^2_x + v^1_x}{2}). \]

If mapped point $ \xi$ does not satisfy condition

  • \[ -1 \leq \xi \leq 1 \]

    then the point $ \xi $ does not belong to reference line $ T^0$ or equivalently point $x $ does not belong to line $ T $ and the method issues error. Otherwise the method returns point $ \xi$.
Parameters
pLocation of point
nodesVertices of element
Returns
vector Vector of shape functions at point p

Reimplemented from fe::BaseElem.

Here is the call graph for this function:

◆ print()

void fe::LineElem::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::LineElem::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: