IPPL (Independent Parallel Particle Layer)
IPPL
Loading...
Searching...
No Matches
FiniteElementSpace.h
Go to the documentation of this file.
1// FiniteElementSpace Class
2// This class is the interface for the solvers with the finite element methods implemented in
3// IPPL. The constructor takes in an IPPL mesh, a reference element and a quadrature rule. The
4// class is templated on the floating point type T, the dimension 'Dim' (1 for 1D, 2 for 2D, 3 for
5// 3D), the number of degrees of freedom per element 'NumElementDOFs', the quadrature rule type
6// 'QuadratureType', the left hand side field type 'FieldLHS' and the right hand side field type
7// 'FieldRHS'.
8
9#ifndef IPPL_FEMSPACE_H
10#define IPPL_FEMSPACE_H
11
12#include "Types/ViewTypes.h"
13
15
19
20namespace ippl {
21
22 constexpr unsigned calculateNumElementVertices(unsigned Dim) {
23 return 1 << Dim; // 2^Dim
24 }
25
38 template <typename T, unsigned Dim, unsigned NumElementDOFs, typename ElementType,
39 typename QuadratureType, typename FieldLHS, typename FieldRHS>
40 // requires IsElement<QuadratureType>
42 public:
43 static constexpr unsigned dim = Dim;
44
45 // the number of mesh vertices per element (not necessarily the same as degrees of freedom,
46 // e.g. a 2D element has 4 vertices)
48 static constexpr unsigned numElementDOFs = NumElementDOFs;
49
50 // A vector with the position of the element in the mesh in each dimension
52
53 // A point in the global coordinate system
55
56 // A vector of vertex indices of the mesh
58
60
62
64 // Constructors ///////////////////////////////////////////////////////
66
74 FiniteElementSpace(UniformCartesian<T, Dim>& mesh, ElementType& ref_element,
75 const QuadratureType& quadrature);
76
78
82
88 KOKKOS_FUNCTION size_t numElements() const;
89
97 KOKKOS_FUNCTION size_t numElementsInDim(const size_t& dim) const;
98
107 KOKKOS_FUNCTION indices_t getMeshVertexNDIndex(const size_t& vertex_index) const;
108
117 KOKKOS_FUNCTION size_t getMeshVertexIndex(const indices_t& vertex_nd_index) const;
118
126 KOKKOS_FUNCTION indices_t getElementNDIndex(const size_t& elementIndex) const;
127
135 KOKKOS_FUNCTION size_t getElementIndex(const indices_t& ndindex) const;
136
145 KOKKOS_FUNCTION vertex_indices_t
146 getElementMeshVertexIndices(const indices_t& elementNDIndex) const;
147
156 KOKKOS_FUNCTION indices_list_t
157 getElementMeshVertexNDIndices(const indices_t& elementNDIndex) const;
158
166 KOKKOS_FUNCTION vertex_points_t
167 getElementMeshVertexPoints(const indices_t& elementNDIndex) const;
168
172
178 KOKKOS_FUNCTION virtual size_t numGlobalDOFs() const = 0;
179
189 KOKKOS_FUNCTION virtual size_t getLocalDOFIndex(const size_t& elementIndex,
190 const size_t& globalDOFIndex) const = 0;
191
200 KOKKOS_FUNCTION virtual size_t getGlobalDOFIndex(const size_t& elementIndex,
201 const size_t& localDOFIndex) const = 0;
202
210 KOKKOS_FUNCTION virtual Vector<size_t, NumElementDOFs> getLocalDOFIndices() const = 0;
211
220 const size_t& elementIndex) const = 0;
221
222
226
228 ElementType ref_element_m;
229 const QuadratureType& quadrature_m;
233 };
234
235} // namespace ippl
236
238
239#endif
constexpr unsigned Dim
Definition Archive.h:20
constexpr unsigned calculateNumElementVertices(unsigned Dim)
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.
virtual KOKKOS_FUNCTION size_t numGlobalDOFs() const =0
Degree of Freedom operations //////////////////////////////////////.
static constexpr unsigned numElementDOFs
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 size_t getElementIndex(const indices_t &ndindex) const
Get the global index of a mesh element given the NDIndex.
FiniteElementSpace(UniformCartesian< T, Dim > &mesh, ElementType &ref_element, const QuadratureType &quadrature)
Construct a new FiniteElementSpace object.
Vector< indices_t, numElementVertices > indices_list_t
static constexpr unsigned dim
Vector< size_t, Dim > indices_t
Vector< size_t, numElementVertices > vertex_indices_t
static constexpr unsigned numElementVertices
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 ...
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).
const QuadratureType & quadrature_m
ElementType ref_element_m
KOKKOS_FUNCTION size_t numElementsInDim(const size_t &dim) const
Get the number of elements in a given dimension.
Vector< size_t, Dim > nr_m
Vector< point_t, numElementVertices > vertex_points_t
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.
KOKKOS_FUNCTION indices_t getMeshVertexNDIndex(const size_t &vertex_index) const
Get the NDIndex of a mesh vertex.
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.
Vector< double, Dim > origin_m
KOKKOS_FUNCTION vertex_points_t getElementMeshVertexPoints(const indices_t &elementNDIndex) const
Get all the global vertex points of an element (given by its NDIndex).
void setMesh(UniformCartesian< T, Dim > &mesh)
UniformCartesian< T, Dim > & mesh_m
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.
Vector< double, Dim > hr_m
KOKKOS_FUNCTION size_t numElements() const
Mesh and Element operations ///////////////////////////////////////.