IPPL (Independent Parallel Particle Layer)
IPPL
Loading...
Searching...
No Matches
LagrangeSpace.h
Go to the documentation of this file.
1// Class LagrangeSpace
2// This is the LagrangeSpace class. It is a class representing a Lagrange space
3// for finite element methods on a structured grid.
4
5#ifndef IPPL_LAGRANGESPACE_H
6#define IPPL_LAGRANGESPACE_H
7
8#include <cmath>
9
11
12constexpr unsigned getLagrangeNumElementDOFs(unsigned Dim, unsigned Order) {
13 // needs to be constexpr pow function to work at compile time. Kokkos::pow doesn't work.
14 return static_cast<unsigned>(power(static_cast<int>(Order + 1), static_cast<int>(Dim)));
15}
16
17namespace ippl {
18
30 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
31 typename QuadratureType, typename FieldLHS, typename FieldRHS>
32 // requires IsQuadrature<QuadratureType>
34 : public FiniteElementSpace<T, Dim, getLagrangeNumElementDOFs(Dim, Order), ElementType,
35 QuadratureType, FieldLHS, FieldRHS> {
36 public:
37 // The number of degrees of freedom per element
38 static constexpr unsigned numElementDOFs = getLagrangeNumElementDOFs(Dim, Order);
39
40 // The dimension of the mesh
41 static constexpr unsigned dim = FiniteElementSpace<T, Dim, numElementDOFs, ElementType,
42 QuadratureType, FieldLHS, FieldRHS>::dim;
43
44 // The order of the Lagrange space
45 static constexpr unsigned order = Order;
46
47 // The number of mesh vertices per element
48 static constexpr unsigned numElementVertices =
49 FiniteElementSpace<T, Dim, numElementDOFs, ElementType, QuadratureType, FieldLHS,
50 FieldRHS>::numElementVertices;
51
52 // A vector with the position of the element in the mesh in each dimension
53 typedef typename FiniteElementSpace<T, Dim, numElementDOFs, ElementType, QuadratureType,
54 FieldLHS, FieldRHS>::indices_t indices_t;
55
56 // A point in the global coordinate system
57 typedef typename FiniteElementSpace<T, Dim, numElementDOFs, ElementType, QuadratureType,
58 FieldLHS, FieldRHS>::point_t point_t;
59
60 typedef typename FiniteElementSpace<T, Dim, numElementDOFs, ElementType, QuadratureType,
62
63 // Field layout type for domain decomposition info
65
66 // View types
70
72 // Constructors ///////////////////////////////////////////////////////
74
83 LagrangeSpace(UniformCartesian<T, Dim>& mesh, ElementType& ref_element,
84 const QuadratureType& quadrature, const Layout_t& layout);
85
95 LagrangeSpace(UniformCartesian<T, Dim>& mesh, ElementType& ref_element,
96 const QuadratureType& quadrature);
97
104 void initialize(UniformCartesian<T, Dim>& mesh, const Layout_t& layout);
105
107
112
115
121 KOKKOS_FUNCTION size_t numGlobalDOFs() const override;
122
132 KOKKOS_FUNCTION size_t getLocalDOFIndex(const size_t& elementIndex,
133 const size_t& globalDOFIndex) const override;
134
143 KOKKOS_FUNCTION size_t getGlobalDOFIndex(const size_t& elementIndex,
144 const size_t& localDOFIndex) const override;
145
153 KOKKOS_FUNCTION Vector<size_t, numElementDOFs> getLocalDOFIndices() const override;
154
163 const size_t& element_index) const override;
164
168
178 KOKKOS_FUNCTION T evaluateRefElementShapeFunction(const size_t& localDOF,
179 const point_t& localPoint) const;
180
192 const size_t& localDOF, const point_t& localPoint) const;
193
197
206 template <typename F>
207 FieldLHS evaluateAx(FieldLHS& field, F& evalFunction) const;
208
209 template <typename F>
210 FieldLHS evaluateAx_lower(FieldLHS& field, F& evalFunction) const;
211
212 template <typename F>
213 FieldLHS evaluateAx_upper(FieldLHS& field, F& evalFunction) const;
214
215 template <typename F>
216 FieldLHS evaluateAx_upperlower(FieldLHS& field, F& evalFunction) const;
217
218 template <typename F>
219 FieldLHS evaluateAx_inversediag(FieldLHS& field, F& evalFunction) const;
220
221 template <typename F>
222 FieldLHS evaluateAx_diag(FieldLHS& field, F& evalFunction) const;
223
234 template <typename F>
235 FieldLHS evaluateAx_lift(FieldLHS& field, F& evalFunction) const;
236
244 void evaluateLoadVector(FieldRHS& field) const;
245
249
258 template <typename F>
259 T computeErrorL2(const FieldLHS& u_h, const F& u_sol) const;
260
268 T computeAvg(const FieldLHS& u_h) const;
269
274 // members we need to copy for the following functions:
275 // works since numElementDOFs in LagrangeSpace is static constexpr
278 ElementType ref_element_m;
279
280 // these are the functions needed for interpolation to the space
281 KOKKOS_FUNCTION indices_t getMeshVertexNDIndex(const size_t& vertex_index) const;
282
283 KOKKOS_FUNCTION size_t getLocalDOFIndex(const indices_t& elementNDIndex,
284 const size_t& globalDOFIndex) const;
286 const indices_t& elementNDIndex) const;
287
288 KOKKOS_FUNCTION T evaluateRefElementShapeFunction(const size_t& localDOF,
289 const point_t& localPoint) const;
290 };
291
292 DeviceStruct getDeviceMirror() const;
293
294 private:
303 KOKKOS_FUNCTION bool isDOFOnBoundary(const indices_t& ndindex) const {
304 for (size_t d = 0; d < Dim; ++d) {
305 if (ndindex[d] <= 0 || ndindex[d] >= this->nr_m[d] - 1) {
306 return true;
307 }
308 }
309 return false;
310 }
311
316 Kokkos::View<size_t*> elementIndices;
317 };
318
319} // namespace ippl
320
321#include "FEM/LagrangeSpace.hpp"
322
323#endif
typename ippl::detail::ViewType< ippl::Vector< double, Dim >, 1 >::view_type view_type
constexpr unsigned Dim
constexpr T power(T base, unsigned exponent)
Definition Quadrature.h:17
constexpr unsigned getLagrangeNumElementDOFs(unsigned Dim, unsigned Order)
Definition Archive.h:20
FiniteElementSpace(UniformCartesian< T, Dim > &mesh, ElementType &ref_element, const QuadratureType &quadrature)
FiniteElementSpace< T, Dim, numElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::vertex_points_t vertex_points_t
T computeAvg(const FieldLHS &u_h) const
Given a field, compute the average.
FieldLHS evaluateAx_upperlower(FieldLHS &field, F &evalFunction) const
KOKKOS_FUNCTION size_t numGlobalDOFs() const override
Degree of Freedom operations //////////////////////////////////////.
detail::ViewType< T, Dim, Kokkos::MemoryTraits< Kokkos::Atomic > >::view_type AtomicViewType
LagrangeSpace(UniformCartesian< T, Dim > &mesh, ElementType &ref_element, const QuadratureType &quadrature, const Layout_t &layout)
Construct a new LagrangeSpace object.
FieldLHS evaluateAx_inversediag(FieldLHS &field, F &evalFunction) const
FieldLHS evaluateAx_lift(FieldLHS &field, F &evalFunction) const
Assemble the left stiffness matrix A of the system but only for the boundary values,...
detail::ViewType< T, Dim >::view_type ViewType
FieldLayout< Dim > Layout_t
KOKKOS_FUNCTION Vector< size_t, numElementDOFs > getGlobalDOFIndices(const size_t &element_index) const override
Get the global DOF indices (vector of global DOF indices) of an element.
void evaluateLoadVector(FieldRHS &field) const
Assemble the load vector b of the system Ax = b.
DeviceStruct getDeviceMirror() const
Device struct definitions /////////////////////////////////////////.
KOKKOS_FUNCTION size_t getLocalDOFIndex(const size_t &elementIndex, const size_t &globalDOFIndex) const override
Get the elements local DOF from the element index and global DOF index.
void initialize(UniformCartesian< T, Dim > &mesh, const Layout_t &layout)
Initialize a LagrangeSpace object created with the default constructor.
LagrangeSpace(UniformCartesian< T, Dim > &mesh, ElementType &ref_element, const QuadratureType &quadrature)
Construct a new LagrangeSpace object (without layout) This constructor is made to work with the defau...
KOKKOS_FUNCTION bool isDOFOnBoundary(const indices_t &ndindex) const
Check if a DOF is on the boundary of the mesh.
FieldLHS evaluateAx(FieldLHS &field, F &evalFunction) const
Assembly operations ///////////////////////////////////////////////.
KOKKOS_FUNCTION T evaluateRefElementShapeFunction(const size_t &localDOF, const point_t &localPoint) const
Basis functions and gradients /////////////////////////////////////.
KOKKOS_FUNCTION Vector< size_t, numElementDOFs > getLocalDOFIndices() const override
Get the local DOF indices (vector of local DOF indices) They are independent of the specific element ...
T computeErrorL2(const FieldLHS &u_h, const F &u_sol) const
Error norm computations ///////////////////////////////////////////.
KOKKOS_FUNCTION point_t evaluateRefElementShapeFunctionGradient(const size_t &localDOF, const point_t &localPoint) const
Evaluate the gradient of the shape function of a local degree of freedom at a given point in the refe...
FieldLHS evaluateAx_diag(FieldLHS &field, F &evalFunction) const
FiniteElementSpace< T, Dim, numElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::point_t point_t
FiniteElementSpace< T, Dim, numElementDOFs, ElementType, QuadratureType, FieldLHS, FieldRHS >::indices_t indices_t
FieldLHS evaluateAx_upper(FieldLHS &field, F &evalFunction) const
void initializeElementIndices(const Layout_t &layout)
Initialize a Kokkos view containing the element indices. This distributes the elements among MPI rank...
FieldLHS evaluateAx_lower(FieldLHS &field, F &evalFunction) const
KOKKOS_FUNCTION size_t getGlobalDOFIndex(const size_t &elementIndex, const size_t &localDOFIndex) const override
Get the global DOF index from the element index and local DOF.
Device struct for copies //////////////////////////////////////////.
KOKKOS_FUNCTION size_t getLocalDOFIndex(const indices_t &elementNDIndex, const size_t &globalDOFIndex) const
KOKKOS_FUNCTION Vector< size_t, numElementDOFs > getGlobalDOFIndices(const indices_t &elementNDIndex) const
KOKKOS_FUNCTION indices_t getMeshVertexNDIndex(const size_t &vertex_index) const
KOKKOS_FUNCTION T evaluateRefElementShapeFunction(const size_t &localDOF, const point_t &localPoint) const
static constexpr unsigned numElementDOFs
Kokkos::View< typename NPtr< T, Dim >::type, Properties... > view_type
Definition ViewTypes.h:45