IBAMR  IBAMR version 0.19.
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
IBTK::FENodalMapping< dim, spacedim, n_nodes > Class Template Referenceabstract

#include <ibtk/FEMapping.h>

Inheritance diagram for IBTK::FENodalMapping< dim, spacedim, n_nodes >:
Inheritance graph
[legend]

Public Types

using key_type = quadrature_key_type
 

Public Member Functions

 FENodalMapping (const key_type quad_key, const libMesh::ElemType mapping_element_type, const FEUpdateFlags update_flags)
 
virtual void reinit (const libMesh::Elem *elem) override
 
virtual const std::vector< double > & getJxW () const override
 
virtual const std::vector< libMesh::Point > & getQuadraturePoints () const override
 
virtual const std::vector< Eigen::Matrix< double, spacedim, dim > > & getContravariants () const override
 
virtual const std::vector< Eigen::Matrix< double, spacedim, dim > > & getCovariants () const override
 
virtual const std::vector< double > & getJxW () const=0
 
virtual const std::vector< libMesh::Point > & getQuadraturePoints () const=0
 
virtual const std::vector< Eigen::Matrix< double, spacedim, dim > > & getContravariants () const=0
 
virtual const std::vector< Eigen::Matrix< double, spacedim, dim > > & getCovariants () const=0
 
std::unique_ptr< FEMapping< 2, 2 > > build (const key_type key, const FEUpdateFlags update_flags)
 
std::unique_ptr< FEMapping< 3, 3 > > build (const key_type key, const FEUpdateFlags update_flags)
 

Static Public Member Functions

static std::unique_ptr< FEMapping< dim, spacedim > > build (const key_type key, const FEUpdateFlags update_flags)
 

Protected Member Functions

virtual bool isAffine () const
 
virtual void fillJacobians () override
 
virtual void fillJxW () override
 
virtual void fillQuadraturePoints (const libMesh::Elem *elem) override
 
virtual void fillTransforms (const libMesh::Elem *elem)=0
 

Protected Attributes

const QuadratureData d_quadrature_data
 
FEUpdateFlags d_update_flags
 
std::vector< Eigen::Matrix< double, spacedim, dim > > d_contravariants
 
std::vector< Eigen::Matrix< double, spacedim, dim > > d_covariants
 
std::vector< doubled_Jacobians
 
std::vector< doubled_JxW
 
std::vector< libMesh::Point > d_quadrature_points
 
PointMap< dim, spacedim, n_nodes > d_point_map
 

Detailed Description

template<int dim, int spacedim = dim, int n_nodes = -1>
class IBTK::FENodalMapping< dim, spacedim, n_nodes >

Base class for all nodal finite element mappings (i.e., mappings corresponding to Lagrange-type finite element spaces).

Template Parameters
n_nodesNumber of nodes of the element: defaults to runtime calculation (-1).

Member Typedef Documentation

◆ key_type

template<int dim, int spacedim = dim, int n_nodes = -1>
using IBTK::FENodalMapping< dim, spacedim, n_nodes >::key_type = quadrature_key_type

Key type. Completely describes (excepting p-refinement) a libMesh quadrature rule.

Constructor & Destructor Documentation

◆ FENodalMapping()

template<int dim, int spacedim = dim, int n_nodes = -1>
IBTK::FENodalMapping< dim, spacedim, n_nodes >::FENodalMapping ( const key_type  quad_key,
const libMesh::ElemType  mapping_element_type,
const FEUpdateFlags  update_flags 
)

Constructor.

Parameters
[in]quad_keyThe quadrature key (i.e., a complete description of the quadrature rule).
[in]mapping_element_typeThe element type used to compute the mapping from the reference element to the physical element. This may be different from the element type in the quadrature rule - for example, one could provide TRI6 in the quadrature rule and TRI3 here.
[in]update_flagsAn enum describing which values need to be computed on each element.

Member Function Documentation

◆ reinit()

template<int dim, int spacedim = dim, int n_nodes = -1>
virtual void IBTK::FENodalMapping< dim, spacedim, n_nodes >::reinit ( const libMesh::Elem *  elem)
overridevirtual

Recalculate relevant quantities for the provided element.

Implements IBTK::FEMapping< dim, dim >.

Reimplemented in IBTK::Hex27Mapping, IBTK::Tet10Mapping, and IBTK::Tri6Mapping.

◆ getJxW() [1/2]

template<int dim, int spacedim = dim, int n_nodes = -1>
virtual const std::vector<double>& IBTK::FENodalMapping< dim, spacedim, n_nodes >::getJxW ( ) const
inlineoverridevirtual

◆ getQuadraturePoints() [1/2]

template<int dim, int spacedim = dim, int n_nodes = -1>
virtual const std::vector<libMesh::Point>& IBTK::FENodalMapping< dim, spacedim, n_nodes >::getQuadraturePoints ( ) const
inlineoverridevirtual

◆ getContravariants() [1/2]

template<int dim, int spacedim = dim, int n_nodes = -1>
virtual const std::vector<Eigen::Matrix<double, spacedim, dim> >& IBTK::FENodalMapping< dim, spacedim, n_nodes >::getContravariants ( ) const
inlineoverridevirtual

◆ getCovariants() [1/2]

template<int dim, int spacedim = dim, int n_nodes = -1>
virtual const std::vector<Eigen::Matrix<double, spacedim, dim> >& IBTK::FENodalMapping< dim, spacedim, n_nodes >::getCovariants ( ) const
inlineoverridevirtual

◆ isAffine()

template<int dim, int spacedim = dim, int n_nodes = -1>
virtual bool IBTK::FENodalMapping< dim, spacedim, n_nodes >::isAffine ( ) const
protectedvirtual

Boolean indicating that the mapping is affine - if it is we can skip some computations. The default implementation returns false. Inheriting classes should overload this they represent affine mappings.

◆ fillJacobians()

template<int dim, int spacedim = dim, int n_nodes = -1>
virtual void IBTK::FENodalMapping< dim, spacedim, n_nodes >::fillJacobians ( )
overrideprotectedvirtual

Compute determinants of contravariants (the Jacobians). The default implementation given here is usually the correct one.

Implements IBTK::FEMapping< dim, dim >.

◆ fillJxW()

template<int dim, int spacedim = dim, int n_nodes = -1>
virtual void IBTK::FENodalMapping< dim, spacedim, n_nodes >::fillJxW ( )
overrideprotectedvirtual

Compute JxW values. The default implementation given here is usually the correct one.

Implements IBTK::FEMapping< dim, dim >.

◆ fillQuadraturePoints()

template<int dim, int spacedim = dim, int n_nodes = -1>
virtual void IBTK::FENodalMapping< dim, spacedim, n_nodes >::fillQuadraturePoints ( const libMesh::Elem *  elem)
overrideprotectedvirtual

Compute the positions of quadrature points on the current element.

Implements IBTK::FEMapping< dim, dim >.

◆ getJxW() [2/2]

virtual const std::vector<double>& IBTK::FEMapping< dim, spacedim >::getJxW
pure virtualinherited

Get the current jacobian times quadrature weight (JxW) values.

◆ getQuadraturePoints() [2/2]

virtual const std::vector<libMesh::Point>& IBTK::FEMapping< dim, spacedim >::getQuadraturePoints
pure virtualinherited

Get the positions of the quadrature points on the current element.

◆ getContravariants() [2/2]

virtual const std::vector<Eigen::Matrix<double, spacedim, dim> >& IBTK::FEMapping< dim, spacedim >::getContravariants
pure virtualinherited

Get the contravariants.

◆ getCovariants() [2/2]

virtual const std::vector<Eigen::Matrix<double, spacedim, dim> >& IBTK::FEMapping< dim, spacedim >::getCovariants
pure virtualinherited

Get the covariants.

◆ build() [1/3]

static std::unique_ptr<FEMapping<dim, spacedim> > IBTK::FEMapping< dim, spacedim >::build ( const key_type  key,
const FEUpdateFlags  update_flags 
)
staticinherited

Return a pointer to the correct mapping for a given quadrature key and update flags object.

◆ build() [2/3]

std::unique_ptr< FEMapping< 2, 2 > > IBTK::FEMapping< 2, 2 >::build ( const key_type  key,
const FEUpdateFlags  update_flags 
)
inherited

◆ build() [3/3]

std::unique_ptr< FEMapping< 3, 3 > > IBTK::FEMapping< 3, 3 >::build ( const key_type  key,
const FEUpdateFlags  update_flags 
)
inherited

◆ fillTransforms()

virtual void IBTK::FEMapping< dim, spacedim >::fillTransforms ( const libMesh::Elem *  elem)
protectedpure virtualinherited

Compute the contravariants and covariants. In general each mapping will have to overload this function.

Implemented in IBTK::FELagrangeMapping< 2, 2, 6 >, IBTK::FELagrangeMapping< 3, 3, 27 >, IBTK::FELagrangeMapping< 3, 3, 8 >, and IBTK::FELagrangeMapping< 3, 3, 10 >.

Member Data Documentation

◆ d_quadrature_data

template<int dim, int spacedim = dim, int n_nodes = -1>
const QuadratureData IBTK::FENodalMapping< dim, spacedim, n_nodes >::d_quadrature_data
protected

Information on the relevant quadrature rule.

◆ d_update_flags

template<int dim, int spacedim = dim, int n_nodes = -1>
FEUpdateFlags IBTK::FENodalMapping< dim, spacedim, n_nodes >::d_update_flags
protected

Computed update flags for the mapping.

◆ d_contravariants

template<int dim, int spacedim = dim, int n_nodes = -1>
std::vector<Eigen::Matrix<double, spacedim, dim> > IBTK::FENodalMapping< dim, spacedim, n_nodes >::d_contravariants
protected

Array of contravariants.

◆ d_covariants

template<int dim, int spacedim = dim, int n_nodes = -1>
std::vector<Eigen::Matrix<double, spacedim, dim> > IBTK::FENodalMapping< dim, spacedim, n_nodes >::d_covariants
protected

Array of covariants (i.e., the transpose of the inverse of the covariants when dim == spacedim)

◆ d_Jacobians

template<int dim, int spacedim = dim, int n_nodes = -1>
std::vector<double> IBTK::FENodalMapping< dim, spacedim, n_nodes >::d_Jacobians
protected

Array of Jacobians.

◆ d_JxW

template<int dim, int spacedim = dim, int n_nodes = -1>
std::vector<double> IBTK::FENodalMapping< dim, spacedim, n_nodes >::d_JxW
protected

Array of JxW values.

◆ d_quadrature_points

template<int dim, int spacedim = dim, int n_nodes = -1>
std::vector<libMesh::Point> IBTK::FENodalMapping< dim, spacedim, n_nodes >::d_quadrature_points
protected

Array of mapped quadrature points.

◆ d_point_map

template<int dim, int spacedim = dim, int n_nodes = -1>
PointMap<dim, spacedim, n_nodes> IBTK::FENodalMapping< dim, spacedim, n_nodes >::d_point_map
protected

Object that computes quadrature point locations. This is sufficiently different from the rest of the mapping code that it is implemented in another class.


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