|
IBAMR
IBAMR version 0.19.
|
#include <ibtk/FEMapping.h>

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< double > | d_Jacobians |
| std::vector< double > | d_JxW |
| std::vector< libMesh::Point > | d_quadrature_points |
| PointMap< dim, spacedim, n_nodes > | d_point_map |
Base class for all nodal finite element mappings (i.e., mappings corresponding to Lagrange-type finite element spaces).
| n_nodes | Number of nodes of the element: defaults to runtime calculation (-1). |
| using IBTK::FENodalMapping< dim, spacedim, n_nodes >::key_type = quadrature_key_type |
Key type. Completely describes (excepting p-refinement) a libMesh quadrature rule.
| IBTK::FENodalMapping< dim, spacedim, n_nodes >::FENodalMapping | ( | const key_type | quad_key, |
| const libMesh::ElemType | mapping_element_type, | ||
| const FEUpdateFlags | update_flags | ||
| ) |
Constructor.
| [in] | quad_key | The quadrature key (i.e., a complete description of the quadrature rule). |
| [in] | mapping_element_type | The 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_flags | An enum describing which values need to be computed on each element. |
|
overridevirtual |
Recalculate relevant quantities for the provided element.
Implements IBTK::FEMapping< dim, dim >.
Reimplemented in IBTK::Hex27Mapping, IBTK::Tet10Mapping, and IBTK::Tri6Mapping.
|
inlineoverridevirtual |
|
inlineoverridevirtual |
|
inlineoverridevirtual |
|
inlineoverridevirtual |
|
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.
|
overrideprotectedvirtual |
Compute determinants of contravariants (the Jacobians). The default implementation given here is usually the correct one.
Implements IBTK::FEMapping< dim, dim >.
|
overrideprotectedvirtual |
Compute JxW values. The default implementation given here is usually the correct one.
Implements IBTK::FEMapping< dim, dim >.
|
overrideprotectedvirtual |
Compute the positions of quadrature points on the current element.
Implements IBTK::FEMapping< dim, dim >.
|
pure virtualinherited |
Get the current jacobian times quadrature weight (JxW) values.
|
pure virtualinherited |
Get the positions of the quadrature points on the current element.
|
pure virtualinherited |
Get the contravariants.
|
pure virtualinherited |
Get the covariants.
|
staticinherited |
Return a pointer to the correct mapping for a given quadrature key and update flags object.
|
inherited |
|
inherited |
|
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 >.
|
protected |
Information on the relevant quadrature rule.
|
protected |
Computed update flags for the mapping.
|
protected |
Array of contravariants.
|
protected |
Array of covariants (i.e., the transpose of the inverse of the covariants when dim == spacedim)
|
protected |
Array of Jacobians.
|
protected |
Array of JxW values.
|
protected |
Array of mapped quadrature points.
|
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.
1.8.17