|
IPPL (Independent Parallel Particle Layer)
IPPL
|
The FiniteElementSpace class handles the mesh index mapping to vertices and elements and is the base class for other FiniteElementSpace classes (e.g. LagrangeSpace). More...
#include <FiniteElementSpace.h>
Public Types | |
| typedef Vector< size_t, Dim > | indices_t |
| typedef Vector< T, Dim > | point_t |
| typedef Vector< size_t, numElementVertices > | vertex_indices_t |
| typedef Vector< indices_t, numElementVertices > | indices_list_t |
| typedef Vector< point_t, numElementVertices > | vertex_points_t |
Public Member Functions | |
| FiniteElementSpace (UniformCartesian< T, Dim > &mesh, ElementType &ref_element, const QuadratureType &quadrature) | |
| Construct a new FiniteElementSpace object. | |
| void | setMesh (UniformCartesian< T, Dim > &mesh) |
| KOKKOS_FUNCTION size_t | numElements () const |
| Mesh and Element operations ///////////////////////////////////////. | |
| KOKKOS_FUNCTION size_t | numElementsInDim (const size_t &dim) const |
| Get the number of elements in a given dimension. | |
| KOKKOS_FUNCTION indices_t | getMeshVertexNDIndex (const size_t &vertex_index) const |
| Get the NDIndex of a mesh vertex. | |
| KOKKOS_FUNCTION size_t | getMeshVertexIndex (const indices_t &vertex_nd_index) const |
| Get the global index of a mesh vertex given its NDIndex. | |
| KOKKOS_FUNCTION indices_t | getElementNDIndex (const size_t &elementIndex) const |
| Get the NDIndex (vector of indices for each dimension) of a mesh element. | |
| KOKKOS_FUNCTION size_t | getElementIndex (const indices_t &ndindex) const |
| Get the global index of a mesh element given the NDIndex. | |
| KOKKOS_FUNCTION vertex_indices_t | getElementMeshVertexIndices (const indices_t &elementNDIndex) const |
| Get all the global vertex indices of an element (given by its NDIndex). | |
| KOKKOS_FUNCTION indices_list_t | getElementMeshVertexNDIndices (const indices_t &elementNDIndex) const |
| Get all the NDIndices of the vertices of an element (given by its NDIndex). | |
| KOKKOS_FUNCTION vertex_points_t | getElementMeshVertexPoints (const indices_t &elementNDIndex) const |
| Get all the global vertex points of an element (given by its NDIndex). | |
| virtual KOKKOS_FUNCTION size_t | numGlobalDOFs () const =0 |
| Degree of Freedom operations //////////////////////////////////////. | |
| virtual KOKKOS_FUNCTION size_t | getLocalDOFIndex (const size_t &elementIndex, const size_t &globalDOFIndex) const =0 |
| Get the elements local DOF from the element index and global DOF index. | |
| virtual KOKKOS_FUNCTION size_t | getGlobalDOFIndex (const size_t &elementIndex, const size_t &localDOFIndex) const =0 |
| Get the global DOF index from the element index and local DOF. | |
| virtual KOKKOS_FUNCTION Vector< size_t, NumElementDOFs > | getLocalDOFIndices () const =0 |
| Get the local DOF indices (vector of local DOF indices) They are independent of the specific element because it only depends on the reference element type. | |
| virtual KOKKOS_FUNCTION Vector< size_t, NumElementDOFs > | getGlobalDOFIndices (const size_t &elementIndex) const =0 |
| Get the global DOF indices (vector of global DOF indices) of an element. | |
Public Attributes | |
| UniformCartesian< T, Dim > & | mesh_m |
| Member variables //////////////////////////////////////////////////. | |
| ElementType | ref_element_m |
| const QuadratureType & | quadrature_m |
| Vector< size_t, Dim > | nr_m |
| Vector< double, Dim > | hr_m |
| Vector< double, Dim > | origin_m |
Static Public Attributes | |
| static constexpr unsigned | dim = Dim |
| static constexpr unsigned | numElementVertices = calculateNumElementVertices(Dim) |
| static constexpr unsigned | numElementDOFs = NumElementDOFs |
The FiniteElementSpace class handles the mesh index mapping to vertices and elements and is the base class for other FiniteElementSpace classes (e.g. LagrangeSpace).
| T | The floating point type |
| Dim | The dimension of the mesh (same dimension as the space) |
| NumElementDOFs | The number of degrees of freedom per element |
| QuadratureType | The type of the quadrature rule (e.g. MidpointQuadrature, GaussJacobiQuadrature) |
| FieldLHS | The type of the left hand side field |
| FieldRHS | The type of the right hand side field (can be the same as FieldLHS) |
Definition at line 41 of file FiniteElementSpace.h.
| typedef Vector<indices_t, numElementVertices> ippl::FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::indices_list_t |
Definition at line 59 of file FiniteElementSpace.h.
| typedef Vector<size_t, Dim> ippl::FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::indices_t |
Definition at line 51 of file FiniteElementSpace.h.
| typedef Vector<T, Dim> ippl::FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::point_t |
Definition at line 54 of file FiniteElementSpace.h.
| typedef Vector<size_t, numElementVertices> ippl::FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::vertex_indices_t |
Definition at line 57 of file FiniteElementSpace.h.
| typedef Vector<point_t, numElementVertices> ippl::FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::vertex_points_t |
Definition at line 61 of file FiniteElementSpace.h.
| ippl::FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::FiniteElementSpace | ( | UniformCartesian< T, Dim > & | mesh, |
| ElementType & | ref_element, | ||
| const QuadratureType & | quadrature ) |
Construct a new FiniteElementSpace object.
| mesh | The mesh object |
| ref_element | The reference element object |
| quadrature | The quadrature rule object |
Definition at line 6 of file FiniteElementSpace.hpp.
References Dim, ippl::Mesh< T, Dim >::Dimension, FiniteElementSpace(), hr_m, mesh_m, nr_m, origin_m, quadrature_m, and ref_element_m.
Referenced by FiniteElementSpace(), getElementIndex(), getElementMeshVertexIndices(), getElementMeshVertexNDIndices(), getElementMeshVertexPoints(), and getMeshVertexIndex().
| KOKKOS_FUNCTION size_t ippl::FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::getElementIndex | ( | const indices_t & | ndindex | ) | const |
Get the global index of a mesh element given the NDIndex.
| ndindex | indices_t (Vector<size_t, Dim>) - vector of indices for each direction |
Definition at line 154 of file FiniteElementSpace.hpp.
References Dim, and FiniteElementSpace().
| KOKKOS_FUNCTION FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::vertex_indices_t ippl::FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::getElementMeshVertexIndices | ( | const indices_t & | elementNDIndex | ) | const |
Get all the global vertex indices of an element (given by its NDIndex).
| elementNDIndex | The NDIndex of the element |
Definition at line 178 of file FiniteElementSpace.hpp.
References Dim, FiniteElementSpace(), and nr_m.
| KOKKOS_FUNCTION FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::indices_list_t ippl::FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::getElementMeshVertexNDIndices | ( | const indices_t & | elementNDIndex | ) | const |
Get all the NDIndices of the vertices of an element (given by its NDIndex).
| elementNDIndex | The NDIndex of the element |
Definition at line 230 of file FiniteElementSpace.hpp.
References Dim, and FiniteElementSpace().
Referenced by getElementMeshVertexPoints().
| KOKKOS_FUNCTION FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::vertex_points_t ippl::FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::getElementMeshVertexPoints | ( | const indices_t & | elementNDIndex | ) | const |
Get all the global vertex points of an element (given by its NDIndex).
| elementNDIndex | The NDIndex of the element |
Definition at line 270 of file FiniteElementSpace.hpp.
References Dim, ippl::Vector< T, Dim >::dim, FiniteElementSpace(), ippl::NDIndex< Dim >::first(), getElementMeshVertexNDIndices(), hr_m, and origin_m.
| KOKKOS_FUNCTION FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::indices_t ippl::FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::getElementNDIndex | ( | const size_t & | elementIndex | ) | const |
Get the NDIndex (vector of indices for each dimension) of a mesh element.
| elementIndex | indices_t (Vector<size_t, Dim>) - The index of the element |
Definition at line 121 of file FiniteElementSpace.hpp.
References getElementNDIndex().
Referenced by getElementNDIndex().
|
pure virtual |
Get the global DOF index from the element index and local DOF.
| elementIndex | size_t - The index of the element |
| localDOFIndex | size_t - The local DOF index |
Implemented in ippl::LagrangeSpace< T, Dim, Order, ElementType, QuadratureType, FieldLHS, FieldRHS >, ippl::LagrangeSpace< Tlhs, Dim, 1, ElementType, QuadratureType, FieldLHS, FieldLHS >, ippl::LagrangeSpace< Tlhs, Dim, 1, ElementType, QuadratureType, FieldLHS, FieldRHS >, ippl::LagrangeSpace< Tlhs, Dim, Order, ElementType, QuadratureType, FieldLHS, FieldRHS >, ippl::NedelecSpace< T, Dim, Order, ElementType, QuadratureType, FieldType >, and ippl::NedelecSpace< T, Dim, 1, ElementType, QuadratureType, FieldType >.
|
pure virtual |
Get the global DOF indices (vector of global DOF indices) of an element.
| elementIndex | size_t - The index of the element |
Implemented in ippl::LagrangeSpace< T, Dim, Order, ElementType, QuadratureType, FieldLHS, FieldRHS >, ippl::LagrangeSpace< Tlhs, Dim, 1, ElementType, QuadratureType, FieldLHS, FieldLHS >, ippl::LagrangeSpace< Tlhs, Dim, 1, ElementType, QuadratureType, FieldLHS, FieldRHS >, ippl::LagrangeSpace< Tlhs, Dim, Order, ElementType, QuadratureType, FieldLHS, FieldRHS >, ippl::NedelecSpace< T, Dim, Order, ElementType, QuadratureType, FieldType >, and ippl::NedelecSpace< T, Dim, 1, ElementType, QuadratureType, FieldType >.
|
pure virtual |
Get the elements local DOF from the element index and global DOF index.
| elementIndex | size_t - The index of the element |
| globalDOFIndex | size_t - The global DOF index |
Implemented in ippl::LagrangeSpace< T, Dim, Order, ElementType, QuadratureType, FieldLHS, FieldRHS >, ippl::LagrangeSpace< Tlhs, Dim, 1, ElementType, QuadratureType, FieldLHS, FieldLHS >, ippl::LagrangeSpace< Tlhs, Dim, 1, ElementType, QuadratureType, FieldLHS, FieldRHS >, ippl::LagrangeSpace< Tlhs, Dim, Order, ElementType, QuadratureType, FieldLHS, FieldRHS >, ippl::NedelecSpace< T, Dim, Order, ElementType, QuadratureType, FieldType >, and ippl::NedelecSpace< T, Dim, 1, ElementType, QuadratureType, FieldType >.
|
pure virtual |
Get the local DOF indices (vector of local DOF indices) They are independent of the specific element because it only depends on the reference element type.
Implemented in ippl::LagrangeSpace< T, Dim, Order, ElementType, QuadratureType, FieldLHS, FieldRHS >, ippl::LagrangeSpace< Tlhs, Dim, 1, ElementType, QuadratureType, FieldLHS, FieldLHS >, ippl::LagrangeSpace< Tlhs, Dim, 1, ElementType, QuadratureType, FieldLHS, FieldRHS >, ippl::LagrangeSpace< Tlhs, Dim, Order, ElementType, QuadratureType, FieldLHS, FieldRHS >, ippl::NedelecSpace< T, Dim, Order, ElementType, QuadratureType, FieldType >, and ippl::NedelecSpace< T, Dim, 1, ElementType, QuadratureType, FieldType >.
| KOKKOS_FUNCTION size_t ippl::FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::getMeshVertexIndex | ( | const indices_t & | vertex_nd_index | ) | const |
Get the global index of a mesh vertex given its NDIndex.
| vertex_nd_index | indices_t (Vector<size_t, Dim>) - The NDIndex of the vertex (vector of indices for each dimension). |
Definition at line 99 of file FiniteElementSpace.hpp.
References Dim, dim, and FiniteElementSpace().
| KOKKOS_FUNCTION FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::indices_t ippl::FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::getMeshVertexNDIndex | ( | const size_t & | vertex_index | ) | const |
Get the NDIndex of a mesh vertex.
| vertex_index | size_t - The index of the vertex |
Definition at line 68 of file FiniteElementSpace.hpp.
References getMeshVertexNDIndex().
Referenced by getMeshVertexNDIndex().
| KOKKOS_FUNCTION size_t ippl::FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::numElements | ( | ) | const |
Mesh and Element operations ///////////////////////////////////////.
Get the number of elements in the mesh of the space
Definition at line 44 of file FiniteElementSpace.hpp.
References Dim, nr_m, and numElements().
Referenced by numElements().
| KOKKOS_FUNCTION size_t ippl::FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::numElementsInDim | ( | const size_t & | dim | ) | const |
Get the number of elements in a given dimension.
| dim | size_t - representing the dimension |
Definition at line 59 of file FiniteElementSpace.hpp.
References dim, nr_m, and numElementsInDim().
Referenced by numElementsInDim().
|
pure virtual |
Degree of Freedom operations //////////////////////////////////////.
Get the number of global degrees of freedom in the space
Implemented in ippl::LagrangeSpace< T, Dim, Order, ElementType, QuadratureType, FieldLHS, FieldRHS >, ippl::LagrangeSpace< Tlhs, Dim, 1, ElementType, QuadratureType, FieldLHS, FieldLHS >, ippl::LagrangeSpace< Tlhs, Dim, 1, ElementType, QuadratureType, FieldLHS, FieldRHS >, ippl::LagrangeSpace< Tlhs, Dim, Order, ElementType, QuadratureType, FieldLHS, FieldRHS >, ippl::NedelecSpace< T, Dim, Order, ElementType, QuadratureType, FieldType >, and ippl::NedelecSpace< T, Dim, 1, ElementType, QuadratureType, FieldType >.
| void ippl::FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::setMesh | ( | UniformCartesian< T, Dim > & | mesh | ) |
Definition at line 26 of file FiniteElementSpace.hpp.
References Dim, ippl::Mesh< T, Dim >::Dimension, hr_m, mesh_m, nr_m, origin_m, and setMesh().
Referenced by setMesh().
|
staticconstexpr |
Definition at line 43 of file FiniteElementSpace.h.
Referenced by getMeshVertexIndex(), and numElementsInDim().
| Vector<double, Dim> ippl::FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::hr_m |
Definition at line 231 of file FiniteElementSpace.h.
Referenced by FiniteElementSpace(), getElementMeshVertexPoints(), and setMesh().
| UniformCartesian<T, Dim>& ippl::FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::mesh_m |
Member variables //////////////////////////////////////////////////.
Definition at line 227 of file FiniteElementSpace.h.
Referenced by FiniteElementSpace(), and setMesh().
| Vector<size_t, Dim> ippl::FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::nr_m |
Definition at line 230 of file FiniteElementSpace.h.
Referenced by FiniteElementSpace(), getElementMeshVertexIndices(), numElements(), numElementsInDim(), and setMesh().
|
staticconstexpr |
Definition at line 48 of file FiniteElementSpace.h.
|
staticconstexpr |
Definition at line 47 of file FiniteElementSpace.h.
| Vector<double, Dim> ippl::FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::origin_m |
Definition at line 232 of file FiniteElementSpace.h.
Referenced by FiniteElementSpace(), getElementMeshVertexPoints(), and setMesh().
| const QuadratureType& ippl::FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::quadrature_m |
Definition at line 229 of file FiniteElementSpace.h.
Referenced by FiniteElementSpace().
| ElementType ippl::FiniteElementSpace< T, Dim, NumElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::ref_element_m |
Definition at line 228 of file FiniteElementSpace.h.
Referenced by FiniteElementSpace().