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

#include <ibtk/FEMapping.h>

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

Public Types

using key_type = quadrature_key_type
 

Public Member Functions

 FELagrangeMapping (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< double > & getJxW () const =0
 
virtual const std::vector< libMesh::Point > & getQuadraturePoints () const override
 
virtual const std::vector< libMesh::Point > & getQuadraturePoints () const =0
 
virtual const std::vector< Eigen::Matrix< double, spacedim, dim > > & getContravariants () const override
 
virtual const std::vector< Eigen::Matrix< double, spacedim, dim > > & getContravariants () const =0
 
virtual const std::vector< Eigen::Matrix< double, spacedim, dim > > & getCovariants () const override
 
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 void fillTransforms (const libMesh::Elem *elem) override
 
virtual bool isAffine () const
 
virtual void fillJacobians () override
 
virtual void fillJxW () override
 
virtual void fillQuadraturePoints (const libMesh::Elem *elem) override
 

Protected Attributes

const int d_n_nodes
 
boost::multi_array< std::array< double, dim >, 2 > d_dphi
 
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
 

Friends

class Hex27Mapping
 

Detailed Description

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

A generic implementation for Lagrange-type elements: works for all elements in that family but is less efficient than the specialized classes for lower-order or tensor-product elements. Supports nonzero codimension.

Member Typedef Documentation

◆ key_type

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

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

Constructor & Destructor Documentation

◆ FELagrangeMapping()

template<int dim, int spacedim = dim, int n_nodes = -1>
IBTK::FELagrangeMapping< dim, spacedim, n_nodes >::FELagrangeMapping ( 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

◆ fillTransforms()

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

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

Implements IBTK::FEMapping< dim, spacedim >.

◆ reinit()

virtual void IBTK::FENodalMapping< dim, spacedim, n_nodes >::reinit ( const libMesh::Elem *  elem)
overridevirtualinherited

Recalculate relevant quantities for the provided element.

Implements IBTK::FEMapping< dim, spacedim >.

◆ getJxW() [1/2]

virtual const std::vector<double>& IBTK::FENodalMapping< dim, spacedim, n_nodes >::getJxW
inlineoverridevirtualinherited

◆ getJxW() [2/2]

template<int dim, int spacedim = dim>
virtual const std::vector<double>& IBTK::FEMapping< dim, spacedim >::getJxW ( ) const
pure virtualinherited

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

◆ getQuadraturePoints() [1/2]

virtual const std::vector<libMesh::Point>& IBTK::FENodalMapping< dim, spacedim, n_nodes >::getQuadraturePoints
inlineoverridevirtualinherited

◆ getQuadraturePoints() [2/2]

template<int dim, int spacedim = dim>
virtual const std::vector<libMesh::Point>& IBTK::FEMapping< dim, spacedim >::getQuadraturePoints ( ) const
pure virtualinherited

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

◆ getContravariants() [1/2]

virtual const std::vector<Eigen::Matrix<double, spacedim, dim> >& IBTK::FENodalMapping< dim, spacedim, n_nodes >::getContravariants
inlineoverridevirtualinherited

◆ getContravariants() [2/2]

template<int dim, int spacedim = dim>
virtual const std::vector<Eigen::Matrix<double, spacedim, dim> >& IBTK::FEMapping< dim, spacedim >::getContravariants ( ) const
pure virtualinherited

Get the contravariants.

◆ getCovariants() [1/2]

virtual const std::vector<Eigen::Matrix<double, spacedim, dim> >& IBTK::FENodalMapping< dim, spacedim, n_nodes >::getCovariants
inlineoverridevirtualinherited

◆ getCovariants() [2/2]

template<int dim, int spacedim = dim>
virtual const std::vector<Eigen::Matrix<double, spacedim, dim> >& IBTK::FEMapping< dim, spacedim >::getCovariants ( ) const
pure virtualinherited

Get the covariants.

◆ isAffine()

virtual bool IBTK::FENodalMapping< dim, spacedim, n_nodes >::isAffine
protectedvirtualinherited

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()

virtual void IBTK::FENodalMapping< dim, spacedim, n_nodes >::fillJacobians
overrideprotectedvirtualinherited

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

Implements IBTK::FEMapping< dim, spacedim >.

◆ fillJxW()

virtual void IBTK::FENodalMapping< dim, spacedim, n_nodes >::fillJxW
overrideprotectedvirtualinherited

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

Implements IBTK::FEMapping< dim, spacedim >.

◆ fillQuadraturePoints()

virtual void IBTK::FENodalMapping< dim, spacedim, n_nodes >::fillQuadraturePoints ( const libMesh::Elem *  elem)
overrideprotectedvirtualinherited

Compute the positions of quadrature points on the current element.

Implements IBTK::FEMapping< dim, spacedim >.

◆ build() [1/3]

template<int dim, int spacedim = dim>
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

Friends And Related Function Documentation

◆ Hex27Mapping

template<int dim, int spacedim = dim, int n_nodes = -1>
friend class Hex27Mapping
friend

Member Data Documentation

◆ d_n_nodes

template<int dim, int spacedim = dim, int n_nodes = -1>
const int IBTK::FELagrangeMapping< dim, spacedim, n_nodes >::d_n_nodes
protected

Number of nodes for the considered element type.

◆ d_dphi

template<int dim, int spacedim = dim, int n_nodes = -1>
boost::multi_array<std::array<double, dim>, 2> IBTK::FELagrangeMapping< dim, spacedim, n_nodes >::d_dphi
protected

Values of shape function gradients on the reference element at quadrature points.

◆ d_quadrature_data

const QuadratureData IBTK::FENodalMapping< dim, spacedim, n_nodes >::d_quadrature_data
protectedinherited

Information on the relevant quadrature rule.

◆ d_update_flags

FEUpdateFlags IBTK::FENodalMapping< dim, spacedim, n_nodes >::d_update_flags
protectedinherited

Computed update flags for the mapping.

◆ d_contravariants

std::vector<Eigen::Matrix<double, spacedim, dim> > IBTK::FENodalMapping< dim, spacedim, n_nodes >::d_contravariants
protectedinherited

Array of contravariants.

◆ d_covariants

std::vector<Eigen::Matrix<double, spacedim, dim> > IBTK::FENodalMapping< dim, spacedim, n_nodes >::d_covariants
protectedinherited

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

◆ d_Jacobians

std::vector<double> IBTK::FENodalMapping< dim, spacedim, n_nodes >::d_Jacobians
protectedinherited

Array of Jacobians.

◆ d_JxW

std::vector<double> IBTK::FENodalMapping< dim, spacedim, n_nodes >::d_JxW
protectedinherited

Array of JxW values.

◆ d_quadrature_points

std::vector<libMesh::Point> IBTK::FENodalMapping< dim, spacedim, n_nodes >::d_quadrature_points
protectedinherited

Array of mapped quadrature points.

◆ d_point_map

PointMap<dim, spacedim, n_nodes> IBTK::FENodalMapping< dim, spacedim, n_nodes >::d_point_map
protectedinherited

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: