IPPL (Independent Parallel Particle Layer)
IPPL
Loading...
Searching...
No Matches
LagrangeSpace.hpp
Go to the documentation of this file.
1
2namespace ippl {
3
4 // LagrangeSpace constructor, which calls the FiniteElementSpace constructor,
5 // and decomposes the elements among ranks according to layout.
6 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
7 typename QuadratureType, typename FieldLHS, typename FieldRHS>
9 UniformCartesian<T, Dim>& mesh, ElementType& ref_element, const QuadratureType& quadrature,
10 const Layout_t& layout)
11 : FiniteElementSpace<T, Dim, getLagrangeNumElementDOFs(Dim, Order), ElementType,
12 QuadratureType, FieldLHS, FieldRHS>(mesh, ref_element, quadrature) {
13 // Assert that the dimension is either 1, 2 or 3.
14 static_assert(Dim >= 1 && Dim <= 3,
15 "Finite Element space only supports 1D, 2D and 3D meshes");
16
17 // Initialize the elementIndices view
19 }
20
21 // LagrangeSpace constructor, which calls the FiniteElementSpace constructor.
22 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
23 typename QuadratureType, typename FieldLHS, typename FieldRHS>
25 UniformCartesian<T, Dim>& mesh, ElementType& ref_element, const QuadratureType& quadrature)
26 : FiniteElementSpace<T, Dim, getLagrangeNumElementDOFs(Dim, Order), ElementType,
27 QuadratureType, FieldLHS, FieldRHS>(mesh, ref_element, quadrature) {
28 // Assert that the dimension is either 1, 2 or 3.
29 static_assert(Dim >= 1 && Dim <= 3,
30 "Finite Element space only supports 1D, 2D and 3D meshes");
31 }
32
33 // LagrangeSpace initializer, to be made available to the FEMPoissonSolver
34 // such that we can call it from setRhs.
35 // Sets the correct mesh ad decomposes the elements among ranks according to layout.
36 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
37 typename QuadratureType, typename FieldLHS, typename FieldRHS>
39 UniformCartesian<T, Dim>& mesh, const Layout_t& layout)
40 {
42 QuadratureType, FieldLHS, FieldRHS>::setMesh(mesh);
43
44 // Initialize the elementIndices view
46 }
47
48 // Initialize element indices Kokkos View by distributing elements among MPI ranks.
49 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
50 typename QuadratureType, typename FieldLHS, typename FieldRHS>
51 void LagrangeSpace<T, Dim, Order, ElementType, QuadratureType, FieldLHS,
52 FieldRHS>::initializeElementIndices(const Layout_t& layout) {
53 const auto& ldom = layout.getLocalNDIndex();
54 int npoints = ldom.size();
55 auto first = ldom.first();
56 auto last = ldom.last();
58
59 for (size_t d = 0; d < Dim; ++d) {
60 bounds[d] = this->nr_m[d] - 1;
61 }
62
63 int upperBoundaryPoints = -1;
64
65 // We iterate over the local domain points, getting the corresponding elements,
66 // while tagging upper boundary points such that they can be removed after.
67 Kokkos::View<size_t*> points("npoints", npoints);
68 Kokkos::parallel_reduce(
69 "ComputePoints", npoints,
70 KOKKOS_CLASS_LAMBDA(const int i, int& local) {
71 int idx = i;
72 indices_t val;
73 bool isBoundary = false;
74 for (unsigned int d = 0; d < Dim; ++d) {
75 int range = last[d] - first[d] + 1;
76 val[d] = first[d] + (idx % range);
77 idx /= range;
78 if (val[d] == bounds[d]) {
79 isBoundary = true;
80 }
81 }
82 points(i) = (!isBoundary) * (this->getElementIndex(val));
83 local += isBoundary;
84 },
85 Kokkos::Sum<int>(upperBoundaryPoints));
86 Kokkos::fence();
87
88 // The elementIndices will be the same array as computed above,
89 // with the tagged upper boundary points removed.
90 int elementsPerRank = npoints - upperBoundaryPoints;
91 elementIndices = Kokkos::View<size_t*>("i", elementsPerRank);
92 Kokkos::View<size_t> index("index");
93
94 Kokkos::parallel_for(
95 "RemoveNaNs", npoints, KOKKOS_CLASS_LAMBDA(const int i) {
96 if ((points(i) != 0) || (i == 0)) {
97 const size_t idx = Kokkos::atomic_fetch_add(&index(), 1);
98 elementIndices(idx) = points(i);
99 }
100 });
101 }
102
104 /// Degree of Freedom operations //////////////////////////////////////
105 ///////////////////////////////////////////////////////////////////////
106
107 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
108 typename QuadratureType, typename FieldLHS, typename FieldRHS>
109 KOKKOS_FUNCTION size_t LagrangeSpace<T, Dim, Order, ElementType, QuadratureType, FieldLHS,
110 FieldRHS>::numGlobalDOFs() const {
111 size_t num_global_dofs = 1;
112 for (size_t d = 0; d < Dim; ++d) {
113 num_global_dofs *= this->nr_m[d] * Order;
114 }
115
116 return num_global_dofs;
117 }
118
119 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
120 typename QuadratureType, typename FieldLHS, typename FieldRHS>
121 KOKKOS_FUNCTION
123 (const size_t& elementIndex, const size_t& globalDOFIndex) const {
124 static_assert(Dim == 1 || Dim == 2 || Dim == 3, "Dim must be 1, 2 or 3");
125 // TODO fix not order independent, only works for order 1
126 static_assert(Order == 1, "Only order 1 is supported at the moment");
127
128 // Get all the global DOFs for the element
129 const Vector<size_t, numElementDOFs> global_dofs =
130 this->LagrangeSpace::getGlobalDOFIndices(elementIndex);
131
132 // Find the global DOF in the vector and return the local DOF index
133 // Note: It is important that this only works because the global_dofs
134 // are already arranged in the correct order from getGlobalDOFIndices
135 for (size_t i = 0; i < global_dofs.dim; ++i) {
136 if (global_dofs[i] == globalDOFIndex) {
137 return i;
138 }
139 }
140 // commented this due to this being on device
141 // however, it would be good to throw an error in this case
142 //throw IpplException("LagrangeSpace::getLocalDOFIndex()",
143 // "FEM Lagrange Space: Global DOF not found in specified element");
144 return 0;
145 }
146
147 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
148 typename QuadratureType, typename FieldLHS, typename FieldRHS>
149 KOKKOS_FUNCTION size_t
150 LagrangeSpace<T, Dim, Order, ElementType, QuadratureType, FieldLHS,
151 FieldRHS>::getGlobalDOFIndex(const size_t& elementIndex,
152 const size_t& localDOFIndex) const {
153 const auto global_dofs = this->LagrangeSpace::getGlobalDOFIndices(elementIndex);
154
155 return global_dofs[localDOFIndex];
156 }
157
158 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
159 typename QuadratureType, typename FieldLHS, typename FieldRHS>
160 KOKKOS_FUNCTION Vector<size_t, LagrangeSpace<T, Dim, Order, ElementType, QuadratureType,
161 FieldLHS, FieldRHS>::numElementDOFs>
162 LagrangeSpace<T, Dim, Order, ElementType, QuadratureType, FieldLHS,
163 FieldRHS>::getLocalDOFIndices() const {
165
166 for (size_t dof = 0; dof < numElementDOFs; ++dof) {
167 localDOFs[dof] = dof;
168 }
169
170 return localDOFs;
171 }
172
173 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
174 typename QuadratureType, typename FieldLHS, typename FieldRHS>
175 KOKKOS_FUNCTION Vector<size_t, LagrangeSpace<T, Dim, Order, ElementType, QuadratureType,
176 FieldLHS, FieldRHS>::numElementDOFs>
177 LagrangeSpace<T, Dim, Order, ElementType, QuadratureType, FieldLHS,
178 FieldRHS>::getGlobalDOFIndices(const size_t& elementIndex) const {
179 Vector<size_t, numElementDOFs> globalDOFs(0);
180
181 // get element pos
182 indices_t elementPos = this->getElementNDIndex(elementIndex);
183
184 // Compute the vector to multiply the ndindex with
186 for (size_t d = 1; d < dim; ++d) {
187 for (size_t d2 = d; d2 < Dim; ++d2) {
188 vec[d2] *= this->nr_m[d - 1];
189 }
190 }
191 vec *= Order; // Multiply each dimension by the order
192 size_t smallestGlobalDOF = elementPos.dot(vec);
193
194 // Add the vertex DOFs
195 globalDOFs[0] = smallestGlobalDOF;
196 globalDOFs[1] = smallestGlobalDOF + Order;
197
198 if constexpr (Dim >= 2) {
199 globalDOFs[2] = globalDOFs[1] + this->nr_m[0] * Order;
200 globalDOFs[3] = globalDOFs[0] + this->nr_m[0] * Order;
201 }
202 if constexpr (Dim >= 3) {
203 globalDOFs[4] = globalDOFs[0] + this->nr_m[1] * this->nr_m[0] * Order;
204 globalDOFs[5] = globalDOFs[1] + this->nr_m[1] * this->nr_m[0] * Order;
205 globalDOFs[6] = globalDOFs[2] + this->nr_m[1] * this->nr_m[0] * Order;
206 globalDOFs[7] = globalDOFs[3] + this->nr_m[1] * this->nr_m[0] * Order;
208
209 if constexpr (Order > 1) {
210 // If the order is greater than 1, there are edge and face DOFs, otherwise the work is
211 // done
212
213 // Add the edge DOFs
214 if constexpr (Dim >= 2) {
215 for (size_t i = 0; i < Order - 1; ++i) {
216 globalDOFs[8 + i] = globalDOFs[0] + i + 1;
217 globalDOFs[8 + Order - 1 + i] = globalDOFs[1] + (i + 1) * this->nr_m[1];
218 globalDOFs[8 + 2 * (Order - 1) + i] = globalDOFs[2] - (i + 1);
219 globalDOFs[8 + 3 * (Order - 1) + i] = globalDOFs[3] - (i + 1) * this->nr_m[1];
220 }
221 }
222 if constexpr (Dim >= 3) {
223 // TODO
224 }
225
226 // Add the face DOFs
227 if constexpr (Dim >= 2) {
228 for (size_t i = 0; i < Order - 1; ++i) {
229 for (size_t j = 0; j < Order - 1; ++j) {
230 // TODO CHECK
231 globalDOFs[8 + 4 * (Order - 1) + i * (Order - 1) + j] =
232 globalDOFs[0] + (i + 1) + (j + 1) * this->nr_m[1];
233 globalDOFs[8 + 4 * (Order - 1) + (Order - 1) * (Order - 1) + i * (Order - 1)
234 + j] = globalDOFs[1] + (i + 1) + (j + 1) * this->nr_m[1];
235 globalDOFs[8 + 4 * (Order - 1) + 2 * (Order - 1) * (Order - 1)
236 + i * (Order - 1) + j] =
237 globalDOFs[2] - (i + 1) + (j + 1) * this->nr_m[1];
238 globalDOFs[8 + 4 * (Order - 1) + 3 * (Order - 1) * (Order - 1)
239 + i * (Order - 1) + j] =
240 globalDOFs[3] - (i + 1) + (j + 1) * this->nr_m[1];
241 }
242 }
243 }
245
246 return globalDOFs;
247 }
248
252
253 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
254 typename QuadratureType, typename FieldLHS, typename FieldRHS>
255 KOKKOS_FUNCTION T
258 const size_t& localDOF,
259 const LagrangeSpace<T, Dim, Order, ElementType, QuadratureType, FieldLHS,
260 FieldRHS>::point_t& localPoint) const {
261 static_assert(Order == 1, "Only order 1 is supported at the moment");
262 // Assert that the local vertex index is valid.
263 assert(localDOF < numElementDOFs
264 && "The local vertex index is invalid"); // TODO assumes 1st order Lagrange
265
266 assert(this->ref_element_m.isPointInRefElement(localPoint)
267 && "Point is not in reference element");
269 // Get the local vertex indices for the local vertex index.
270 // TODO fix not order independent, only works for order 1
271 const point_t ref_element_point = this->ref_element_m.getLocalVertices()[localDOF];
272
273 // The variable that accumulates the product of the shape functions.
274 T product = 1;
275
276 for (size_t d = 0; d < Dim; d++) {
277 if (localPoint[d] < ref_element_point[d]) {
278 product *= localPoint[d];
279 } else {
280 product *= 1.0 - localPoint[d];
282 }
284 return product;
286
287 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
288 typename QuadratureType, typename FieldLHS, typename FieldRHS>
289 KOKKOS_FUNCTION typename LagrangeSpace<T, Dim, Order, ElementType, QuadratureType, FieldLHS,
290 FieldRHS>::point_t
293 const size_t& localDOF,
294 const LagrangeSpace<T, Dim, Order, ElementType, QuadratureType, FieldLHS,
295 FieldRHS>::point_t& localPoint) const {
296 // TODO fix not order independent, only works for order 1
297 static_assert(Order == 1 && "Only order 1 is supported at the moment");
298
299 // Assert that the local vertex index is valid.
300 assert(localDOF < numElementDOFs && "The local vertex index is invalid");
301
302 assert(this->ref_element_m.isPointInRefElement(localPoint)
303 && "Point is not in reference element");
304
305 // Get the local dof nd_index
306 const vertex_points_t local_vertex_points = this->ref_element_m.getLocalVertices();
307
308 const point_t& local_vertex_point = local_vertex_points[localDOF];
309
310 point_t gradient(1);
311
312 // To construct the gradient we need to loop over the dimensions and multiply the
313 // shape functions in each dimension except the current one. The one of the current
314 // dimension is replaced by the derivative of the shape function in that dimension,
315 // which is either 1 or -1.
316 for (size_t d = 0; d < Dim; d++) {
317 // The variable that accumulates the product of the shape functions.
318 T product = 1;
319
320 for (size_t d2 = 0; d2 < Dim; d2++) {
321 if (d2 == d) {
322 if (localPoint[d] < local_vertex_point[d]) {
323 product *= 1;
324 } else {
325 product *= -1;
326 }
327 } else {
328 if (localPoint[d2] < local_vertex_point[d2]) {
329 product *= localPoint[d2];
330 } else {
331 product *= 1.0 - localPoint[d2];
332 }
333 }
334 }
335
336 gradient[d] = product;
337 }
338
339 return gradient;
340 }
341
345
346 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
347 typename QuadratureType, typename FieldLHS, typename FieldRHS>
348 template <typename F>
349 FieldLHS LagrangeSpace<T, Dim, Order, ElementType, QuadratureType, FieldLHS,
350 FieldRHS>::evaluateAx(FieldLHS& field, F& evalFunction) const {
351 Inform m("");
352
353 // start a timer
354 static IpplTimings::TimerRef evalAx = IpplTimings::getTimer("evaluateAx");
356
357 // get number of ghost cells in field
358 const int nghost = field.getNghost();
359
360 // create a new field for result with view initialized to zero (views are initialized to
361 // zero by default)
362 FieldLHS resultField(field.get_mesh(), field.getLayout(), nghost);
363
364 // List of quadrature weights
366 this->quadrature_m.getWeightsForRefElement();
367
368 // List of quadrature nodes
370 this->quadrature_m.getIntegrationNodesForRefElement();
371
372 // TODO move outside of evaluateAx (I think it is possible for other problems as well)
373 // Gradients of the basis functions for the DOF at the quadrature nodes
374 Vector<Vector<point_t, numElementDOFs>, QuadratureType::numElementNodes> grad_b_q;
375 for (size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
376 for (size_t i = 0; i < numElementDOFs; ++i) {
377 grad_b_q[k][i] = this->evaluateRefElementShapeFunctionGradient(i, q[k]);
378 }
379 }
380
381 // Make local element matrix -- does not change through the element mesh
382 // Element matrix
384
385 // 1. Compute the Galerkin element matrix A_K
386 for (size_t i = 0; i < numElementDOFs; ++i) {
387 for (size_t j = 0; j < numElementDOFs; ++j) {
388 A_K[i][j] = 0.0;
389 for (size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
390 A_K[i][j] += w[k] * evalFunction(i, j, grad_b_q[k]);
391 }
392 }
393 }
394
395 // Get field data and atomic result data,
396 // since it will be added to during the kokkos loop
397 ViewType view = field.getView();
398 AtomicViewType resultView = resultField.getView();
399
400 // Get boundary conditions from field
401 BConds<FieldLHS, Dim>& bcField = field.getFieldBC();
402 FieldBC bcType = bcField[0]->getBCType();
403
404 // Get domain information
405 auto ldom = (field.getLayout()).getLocalNDIndex();
406
407 using exec_space = typename Kokkos::View<const size_t*>::execution_space;
408 using policy_type = Kokkos::RangePolicy<exec_space>;
409
410 // start a timer
411 static IpplTimings::TimerRef outer_loop = IpplTimings::getTimer("evaluateAx: outer loop");
412 IpplTimings::startTimer(outer_loop);
413
414 // Loop over elements to compute contributions
415 Kokkos::parallel_for(
416 "Loop over elements", policy_type(0, elementIndices.extent(0)),
417 KOKKOS_CLASS_LAMBDA(const size_t index) {
418 const size_t elementIndex = elementIndices(index);
419 const Vector<size_t, numElementDOFs> global_dofs =
420 this->LagrangeSpace::getGlobalDOFIndices(elementIndex);
421 Vector<indices_t, numElementDOFs> global_dof_ndindices;
422
423 for (size_t i = 0; i < numElementDOFs; ++i) {
424 global_dof_ndindices[i] = this->getMeshVertexNDIndex(global_dofs[i]);
425 }
426
427 // local DOF indices (both i and j go from 0 to numDOFs-1 in the element)
428 size_t i, j;
429
430 // global DOF n-dimensional indices (Vector of N indices representing indices in
431 // each dimension)
432 indices_t I_nd, J_nd;
433
434 // 2. Compute the contribution to resultAx = A*x with A_K
435 for (i = 0; i < numElementDOFs; ++i) {
436 I_nd = global_dof_ndindices[i];
437
438 // Handle boundary DOFs
439 // If Zero Dirichlet BCs, skip this DOF
440 // If Constant Dirichlet BCs, identity
441 if ((bcType == CONSTANT_FACE) && (this->isDOFOnBoundary(I_nd))) {
442 for (unsigned d = 0; d < Dim; ++d) {
443 I_nd[d] = I_nd[d] - ldom[d].first() + nghost;
444 }
445 apply(resultView, I_nd) = apply(view, I_nd);
446 continue;
447 } else if ((bcType == ZERO_FACE) && (this->isDOFOnBoundary(I_nd))) {
448 continue;
449 }
450
451 // get the appropriate index for the Kokkos view of the field
452 for (unsigned d = 0; d < Dim; ++d) {
453 I_nd[d] = I_nd[d] - ldom[d].first() + nghost;
454 }
455
456 for (j = 0; j < numElementDOFs; ++j) {
457 J_nd = global_dof_ndindices[j];
458
459 // Skip boundary DOFs (Zero & Constant Dirichlet BCs)
460 if (((bcType == ZERO_FACE) || (bcType == CONSTANT_FACE))
461 && this->isDOFOnBoundary(J_nd)) {
462 continue;
463 }
464
465 // get the appropriate index for the Kokkos view of the field
466 for (unsigned d = 0; d < Dim; ++d) {
467 J_nd[d] = J_nd[d] - ldom[d].first() + nghost;
468 }
469
470 apply(resultView, I_nd) += A_K[i][j] * apply(view, J_nd);
471 }
472 }
473 });
474 IpplTimings::stopTimer(outer_loop);
475
476 if (bcType == PERIODIC_FACE) {
477 resultField.accumulateHalo();
478 bcField.apply(resultField);
479 bcField.assignGhostToPhysical(resultField);
480 } else {
481 resultField.accumulateHalo_noghost();
482 }
483
485
486 return resultField;
487 }
488
489 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
490 typename QuadratureType, typename FieldLHS, typename FieldRHS>
491 template <typename F>
492 FieldLHS LagrangeSpace<T, Dim, Order, ElementType, QuadratureType, FieldLHS,
493 FieldRHS>::evaluateAx_lower(FieldLHS& field, F& evalFunction) const {
494 Inform m("");
495
496 // start a timer
497 static IpplTimings::TimerRef evalAx = IpplTimings::getTimer("evaluateAx");
499
500 // get number of ghost cells in field
501 const int nghost = field.getNghost();
502
503 // create a new field for result with view initialized to zero (views are initialized to
504 // zero by default)
505 FieldLHS resultField(field.get_mesh(), field.getLayout(), nghost);
506
507 // List of quadrature weights
509 this->quadrature_m.getWeightsForRefElement();
510
511 // List of quadrature nodes
513 this->quadrature_m.getIntegrationNodesForRefElement();
514
515 // TODO move outside of evaluateAx (I think it is possible for other problems as well)
516 // Gradients of the basis functions for the DOF at the quadrature nodes
517 Vector<Vector<point_t, numElementDOFs>, QuadratureType::numElementNodes> grad_b_q;
518 for (size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
519 for (size_t i = 0; i < numElementDOFs; ++i) {
520 grad_b_q[k][i] = this->evaluateRefElementShapeFunctionGradient(i, q[k]);
521 }
522 }
523
524 // Make local element matrix -- does not change through the element mesh
525 // Element matrix
527
528 // 1. Compute the Galerkin element matrix A_K
529 for (size_t i = 0; i < numElementDOFs; ++i) {
530 for (size_t j = 0; j < numElementDOFs; ++j) {
531 A_K[i][j] = 0.0;
532 for (size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
533 A_K[i][j] += w[k] * evalFunction(i, j, grad_b_q[k]);
534 }
535 }
536 }
537
538 // Get field data and atomic result data,
539 // since it will be added to during the kokkos loop
540 ViewType view = field.getView();
541 AtomicViewType resultView = resultField.getView();
542
543 // Get boundary conditions from field
544 BConds<FieldLHS, Dim>& bcField = field.getFieldBC();
545 FieldBC bcType = bcField[0]->getBCType();
546
547 // Get domain information
548 auto ldom = (field.getLayout()).getLocalNDIndex();
549
550 using exec_space = typename Kokkos::View<const size_t*>::execution_space;
551 using policy_type = Kokkos::RangePolicy<exec_space>;
552
553 // start a timer
554 static IpplTimings::TimerRef outer_loop = IpplTimings::getTimer("evaluateAx: outer loop");
555 IpplTimings::startTimer(outer_loop);
556
557 // Loop over elements to compute contributions
558 Kokkos::parallel_for(
559 "Loop over elements", policy_type(0, elementIndices.extent(0)),
560 KOKKOS_CLASS_LAMBDA(const size_t index) {
561 const size_t elementIndex = elementIndices(index);
562 const Vector<size_t, numElementDOFs> global_dofs =
563 this->LagrangeSpace::getGlobalDOFIndices(elementIndex);
564 Vector<indices_t, numElementDOFs> global_dof_ndindices;
565
566 for (size_t i = 0; i < numElementDOFs; ++i) {
567 global_dof_ndindices[i] = this->getMeshVertexNDIndex(global_dofs[i]);
568 }
569
570 // local DOF indices
571 size_t i, j;
572
573 // global DOF n-dimensional indices (Vector of N indices representing indices in
574 // each dimension)
575 indices_t I_nd, J_nd;
576
577 // 2. Compute the contribution to resultAx = A*x with A_K
578 for (i = 0; i < numElementDOFs; ++i) {
579 I_nd = global_dof_ndindices[i];
580
581 // Handle boundary DOFs
582 // If Zero Dirichlet BCs, skip this DOF
583 // If Constant Dirichlet BCs, identity
584 if ((bcType == CONSTANT_FACE) && (this->isDOFOnBoundary(I_nd))) {
585 for (unsigned d = 0; d < Dim; ++d) {
586 I_nd[d] = I_nd[d] - ldom[d].first() + nghost;
587 }
588 apply(resultView, I_nd) = apply(view, I_nd);
589 continue;
590 } else if ((bcType == ZERO_FACE) && (this->isDOFOnBoundary(I_nd))) {
591 continue;
592 }
593
594 // get the appropriate index for the Kokkos view of the field
595 for (unsigned d = 0; d < Dim; ++d) {
596 I_nd[d] = I_nd[d] - ldom[d].first() + nghost;
597 }
598
599 for (j = 0; j < numElementDOFs; ++j) {
600 J_nd = global_dof_ndindices[j];
601
602 if (global_dofs[i] >= global_dofs[j]) {
603 continue;
604 }
605
606 // Skip boundary DOFs (Zero & Constant Dirichlet BCs)
607 if (((bcType == ZERO_FACE) || (bcType == CONSTANT_FACE))
608 && this->isDOFOnBoundary(J_nd)) {
609 continue;
610 }
611
612 // get the appropriate index for the Kokkos view of the field
613 for (unsigned d = 0; d < Dim; ++d) {
614 J_nd[d] = J_nd[d] - ldom[d].first() + nghost;
615 }
616
617 apply(resultView, I_nd) += A_K[i][j] * apply(view, J_nd);
618 }
619 }
620 });
621 IpplTimings::stopTimer(outer_loop);
622
623 if (bcType == PERIODIC_FACE) {
624 resultField.accumulateHalo();
625 bcField.apply(resultField);
626 bcField.assignGhostToPhysical(resultField);
627 } else {
628 resultField.accumulateHalo_noghost();
629 }
630
632
633 return resultField;
634 }
635
636 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
637 typename QuadratureType, typename FieldLHS, typename FieldRHS>
638 template <typename F>
639 FieldLHS LagrangeSpace<T, Dim, Order, ElementType, QuadratureType, FieldLHS,
640 FieldRHS>::evaluateAx_upper(FieldLHS& field, F& evalFunction) const {
641 Inform m("");
642
643 // start a timer
644 static IpplTimings::TimerRef evalAx = IpplTimings::getTimer("evaluateAx");
646
647 // get number of ghost cells in field
648 const int nghost = field.getNghost();
649
650 // create a new field for result with view initialized to zero (views are initialized to
651 // zero by default)
652 FieldLHS resultField(field.get_mesh(), field.getLayout(), nghost);
653
654 // List of quadrature weights
656 this->quadrature_m.getWeightsForRefElement();
657
658 // List of quadrature nodes
660 this->quadrature_m.getIntegrationNodesForRefElement();
661
662 // TODO move outside of evaluateAx (I think it is possible for other problems as well)
663 // Gradients of the basis functions for the DOF at the quadrature nodes
664 Vector<Vector<point_t, numElementDOFs>, QuadratureType::numElementNodes> grad_b_q;
665 for (size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
666 for (size_t i = 0; i < numElementDOFs; ++i) {
667 grad_b_q[k][i] = this->evaluateRefElementShapeFunctionGradient(i, q[k]);
668 }
669 }
670
671 // Make local element matrix -- does not change through the element mesh
672 // Element matrix
674
675 // 1. Compute the Galerkin element matrix A_K
676 for (size_t i = 0; i < numElementDOFs; ++i) {
677 for (size_t j = 0; j < numElementDOFs; ++j) {
678 A_K[i][j] = 0.0;
679 for (size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
680 A_K[i][j] += w[k] * evalFunction(i, j, grad_b_q[k]);
681 }
682 }
683 }
684
685 // Get field data and atomic result data,
686 // since it will be added to during the kokkos loop
687 ViewType view = field.getView();
688 AtomicViewType resultView = resultField.getView();
689
690 // Get boundary conditions from field
691 BConds<FieldLHS, Dim>& bcField = field.getFieldBC();
692 FieldBC bcType = bcField[0]->getBCType();
693
694 // Get domain information
695 auto ldom = (field.getLayout()).getLocalNDIndex();
696
697 using exec_space = typename Kokkos::View<const size_t*>::execution_space;
698 using policy_type = Kokkos::RangePolicy<exec_space>;
699
700 // start a timer
701 static IpplTimings::TimerRef outer_loop = IpplTimings::getTimer("evaluateAx: outer loop");
702 IpplTimings::startTimer(outer_loop);
703
704 // Loop over elements to compute contributions
705 Kokkos::parallel_for(
706 "Loop over elements", policy_type(0, elementIndices.extent(0)),
707 KOKKOS_CLASS_LAMBDA(const size_t index) {
708 const size_t elementIndex = elementIndices(index);
709 const Vector<size_t, numElementDOFs> global_dofs =
710 this->LagrangeSpace::getGlobalDOFIndices(elementIndex);
711 Vector<indices_t, numElementDOFs> global_dof_ndindices;
712
713 for (size_t i = 0; i < numElementDOFs; ++i) {
714 global_dof_ndindices[i] = this->getMeshVertexNDIndex(global_dofs[i]);
715 }
716
717 // local DOF indices
718 size_t i, j;
719
720 // global DOF n-dimensional indices (Vector of N indices representing indices in
721 // each dimension)
722 indices_t I_nd, J_nd;
723
724 // 2. Compute the contribution to resultAx = A*x with A_K
725 for (i = 0; i < numElementDOFs; ++i) {
726 I_nd = global_dof_ndindices[i];
727
728 // Handle boundary DOFs
729 // If Zero Dirichlet BCs, skip this DOF
730 // If Constant Dirichlet BCs, identity
731 if ((bcType == CONSTANT_FACE) && (this->isDOFOnBoundary(I_nd))) {
732 for (unsigned d = 0; d < Dim; ++d) {
733 I_nd[d] = I_nd[d] - ldom[d].first() + nghost;
734 }
735 apply(resultView, I_nd) = apply(view, I_nd);
736 continue;
737 } else if ((bcType == ZERO_FACE) && (this->isDOFOnBoundary(I_nd))) {
738 continue;
739 }
740
741 // get the appropriate index for the Kokkos view of the field
742 for (unsigned d = 0; d < Dim; ++d) {
743 I_nd[d] = I_nd[d] - ldom[d].first() + nghost;
744 }
745
746 for (j = 0; j < numElementDOFs; ++j) {
747 J_nd = global_dof_ndindices[j];
748
749 if (global_dofs[i] <= global_dofs[j]) {
750 continue;
751 }
752
753 // Skip boundary DOFs (Zero & Constant Dirichlet BCs)
754 if (((bcType == ZERO_FACE) || (bcType == CONSTANT_FACE))
755 && this->isDOFOnBoundary(J_nd)) {
756 continue;
757 }
758
759 // get the appropriate index for the Kokkos view of the field
760 for (unsigned d = 0; d < Dim; ++d) {
761 J_nd[d] = J_nd[d] - ldom[d].first() + nghost;
762 }
763
764 apply(resultView, I_nd) += A_K[i][j] * apply(view, J_nd);
765 }
766 }
767 });
768 IpplTimings::stopTimer(outer_loop);
769
770 if (bcType == PERIODIC_FACE) {
771 resultField.accumulateHalo();
772 bcField.apply(resultField);
773 bcField.assignGhostToPhysical(resultField);
774 } else {
775 resultField.accumulateHalo_noghost();
776 }
777
779
780 return resultField;
781 }
782
783 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
784 typename QuadratureType, typename FieldLHS, typename FieldRHS>
785 template <typename F>
786 FieldLHS LagrangeSpace<T, Dim, Order, ElementType, QuadratureType, FieldLHS,
787 FieldRHS>::evaluateAx_upperlower(FieldLHS& field, F& evalFunction) const {
788 Inform m("");
789
790 // start a timer
791 static IpplTimings::TimerRef evalAx = IpplTimings::getTimer("evaluateAx");
793
794 // get number of ghost cells in field
795 const int nghost = field.getNghost();
796
797 // create a new field for result with view initialized to zero (views are initialized to
798 // zero by default)
799 FieldLHS resultField(field.get_mesh(), field.getLayout(), nghost);
800
801 // List of quadrature weights
803 this->quadrature_m.getWeightsForRefElement();
804
805 // List of quadrature nodes
807 this->quadrature_m.getIntegrationNodesForRefElement();
808
809 // TODO move outside of evaluateAx (I think it is possible for other problems as well)
810 // Gradients of the basis functions for the DOF at the quadrature nodes
811 Vector<Vector<point_t, numElementDOFs>, QuadratureType::numElementNodes> grad_b_q;
812 for (size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
813 for (size_t i = 0; i < numElementDOFs; ++i) {
814 grad_b_q[k][i] = this->evaluateRefElementShapeFunctionGradient(i, q[k]);
815 }
816 }
817
818 // Make local element matrix -- does not change through the element mesh
819 // Element matrix
821
822 // 1. Compute the Galerkin element matrix A_K
823 for (size_t i = 0; i < numElementDOFs; ++i) {
824 for (size_t j = 0; j < numElementDOFs; ++j) {
825 A_K[i][j] = 0.0;
826 for (size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
827 A_K[i][j] += w[k] * evalFunction(i, j, grad_b_q[k]);
828 }
829 }
830 }
831
832 // Get field data and atomic result data,
833 // since it will be added to during the kokkos loop
834 ViewType view = field.getView();
835 AtomicViewType resultView = resultField.getView();
836
837 // Get boundary conditions from field
838 BConds<FieldLHS, Dim>& bcField = field.getFieldBC();
839 FieldBC bcType = bcField[0]->getBCType();
840
841 // Get domain information
842 auto ldom = (field.getLayout()).getLocalNDIndex();
843
844 using exec_space = typename Kokkos::View<const size_t*>::execution_space;
845 using policy_type = Kokkos::RangePolicy<exec_space>;
846
847 // start a timer
848 static IpplTimings::TimerRef outer_loop = IpplTimings::getTimer("evaluateAx: outer loop");
849 IpplTimings::startTimer(outer_loop);
850
851 // Loop over elements to compute contributions
852 Kokkos::parallel_for(
853 "Loop over elements", policy_type(0, elementIndices.extent(0)),
854 KOKKOS_CLASS_LAMBDA(const size_t index) {
855 const size_t elementIndex = elementIndices(index);
856 const Vector<size_t, numElementDOFs> global_dofs =
857 this->LagrangeSpace::getGlobalDOFIndices(elementIndex);
858 Vector<indices_t, numElementDOFs> global_dof_ndindices;
859
860 for (size_t i = 0; i < numElementDOFs; ++i) {
861 global_dof_ndindices[i] = this->getMeshVertexNDIndex(global_dofs[i]);
862 }
863
864 // local DOF indices
865 size_t i, j;
866
867 // global DOF n-dimensional indices (Vector of N indices representing indices in
868 // each dimension)
869 indices_t I_nd, J_nd;
870
871 // 2. Compute the contribution to resultAx = A*x with A_K
872 for (i = 0; i < numElementDOFs; ++i) {
873 I_nd = global_dof_ndindices[i];
874
875 // Handle boundary DOFs
876 // If Zero Dirichlet BCs, skip this DOF
877 // If Constant Dirichlet BCs, identity
878 if ((bcType == CONSTANT_FACE) && (this->isDOFOnBoundary(I_nd))) {
879 for (unsigned d = 0; d < Dim; ++d) {
880 I_nd[d] = I_nd[d] - ldom[d].first() + nghost;
881 }
882 apply(resultView, I_nd) = apply(view, I_nd);
883 continue;
884 } else if ((bcType == ZERO_FACE) && (this->isDOFOnBoundary(I_nd))) {
885 continue;
886 }
887
888 // get the appropriate index for the Kokkos view of the field
889 for (unsigned d = 0; d < Dim; ++d) {
890 I_nd[d] = I_nd[d] - ldom[d].first() + nghost;
891 }
892
893 for (j = 0; j < i; ++j) {
894 J_nd = global_dof_ndindices[j];
895
896 // Skip boundary DOFs (Zero & Constant Dirichlet BCs)
897 if (((bcType == ZERO_FACE) || (bcType == CONSTANT_FACE))
898 && this->isDOFOnBoundary(J_nd)) {
899 continue;
900 }
901
902 // get the appropriate index for the Kokkos view of the field
903 for (unsigned d = 0; d < Dim; ++d) {
904 J_nd[d] = J_nd[d] - ldom[d].first() + nghost;
905 }
906
907 apply(resultView, I_nd) += A_K[i][j] * apply(view, J_nd);
908 apply(resultView, J_nd) += A_K[j][i] * apply(view, I_nd);
909 }
910 }
911 });
912 IpplTimings::stopTimer(outer_loop);
913
914 if (bcType == PERIODIC_FACE) {
915 resultField.accumulateHalo();
916 bcField.apply(resultField);
917 bcField.assignGhostToPhysical(resultField);
918 } else {
919 resultField.accumulateHalo_noghost();
920 }
921
923
924 return resultField;
925 }
926
927 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
928 typename QuadratureType, typename FieldLHS, typename FieldRHS>
929 template <typename F>
930 FieldLHS LagrangeSpace<T, Dim, Order, ElementType, QuadratureType, FieldLHS,
931 FieldRHS>::evaluateAx_inversediag(FieldLHS& field, F& evalFunction) const {
932 Inform m("");
933
934 // start a timer
935 static IpplTimings::TimerRef evalAx = IpplTimings::getTimer("evaluateAx");
937
938 // get number of ghost cells in field
939 const int nghost = field.getNghost();
940
941 // create a new field for result with view initialized to zero (views are initialized to
942 // zero by default)
943 FieldLHS resultField(field.get_mesh(), field.getLayout(), nghost);
944
945 // List of quadrature weights
947 this->quadrature_m.getWeightsForRefElement();
948
949 // List of quadrature nodes
951 this->quadrature_m.getIntegrationNodesForRefElement();
952
953 // TODO move outside of evaluateAx (I think it is possible for other problems as well)
954 // Gradients of the basis functions for the DOF at the quadrature nodes
955 Vector<Vector<point_t, numElementDOFs>, QuadratureType::numElementNodes> grad_b_q;
956 for (size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
957 for (size_t i = 0; i < numElementDOFs; ++i) {
958 grad_b_q[k][i] = this->evaluateRefElementShapeFunctionGradient(i, q[k]);
959 }
960 }
961
962 // Make local element matrix -- does not change through the element mesh
963 // Element matrix
965
966 // 1. Compute the Galerkin element matrix A_K
967 for (size_t i = 0; i < numElementDOFs; ++i) {
968 A_K_diag[i] = 0.0;
969 for (size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
970 A_K_diag[i] += w[k] * evalFunction(i, i, grad_b_q[k]);
971 }
972 }
973
974 // Get field data and atomic result data,
975 // since it will be added to during the kokkos loop
976 ViewType view = field.getView();
977 AtomicViewType resultView = resultField.getView();
978
979 // Get boundary conditions from field
980 BConds<FieldLHS, Dim>& bcField = field.getFieldBC();
981 FieldBC bcType = bcField[0]->getBCType();
982
983 // Get domain information
984 auto ldom = (field.getLayout()).getLocalNDIndex();
985
986 using exec_space = typename Kokkos::View<const size_t*>::execution_space;
987 using policy_type = Kokkos::RangePolicy<exec_space>;
988
989 // start a timer
990 static IpplTimings::TimerRef outer_loop = IpplTimings::getTimer("evaluateAx: outer loop");
991 IpplTimings::startTimer(outer_loop);
992
993 // Loop over elements to compute contributions
994 Kokkos::parallel_for(
995 "Loop over elements", policy_type(0, elementIndices.extent(0)),
996 KOKKOS_CLASS_LAMBDA(const size_t index) {
997 const size_t elementIndex = elementIndices(index);
998 const Vector<size_t, numElementDOFs> global_dofs =
999 this->LagrangeSpace::getGlobalDOFIndices(elementIndex);
1000 Vector<indices_t, numElementDOFs> global_dof_ndindices;
1001
1002 for (size_t i = 0; i < numElementDOFs; ++i) {
1003 global_dof_ndindices[i] = this->getMeshVertexNDIndex(global_dofs[i]);
1004 }
1005
1006 // local DOF indices
1007 size_t i;
1008
1009 // global DOF n-dimensional indices (Vector of N indices representing indices in
1010 // each dimension)
1011 indices_t I_nd;
1012
1013 // 2. Compute the contribution to resultAx = A*x with A_K
1014 for (i = 0; i < numElementDOFs; ++i) {
1015 I_nd = global_dof_ndindices[i];
1016
1017 // Handle boundary DOFs
1018 // If Zero Dirichlet BCs, skip this DOF
1019 // If Constant Dirichlet BCs, identity
1020 if ((bcType == CONSTANT_FACE) && (this->isDOFOnBoundary(I_nd))) {
1021 for (unsigned d = 0; d < Dim; ++d) {
1022 I_nd[d] = I_nd[d] - ldom[d].first() + nghost;
1023 }
1024 apply(resultView, I_nd) = 1.0;
1025 continue;
1026 } else if ((bcType == ZERO_FACE) && (this->isDOFOnBoundary(I_nd))) {
1027 continue;
1028 }
1029
1030 // get the appropriate index for the Kokkos view of the field
1031 for (unsigned d = 0; d < Dim; ++d) {
1032 I_nd[d] = I_nd[d] - ldom[d].first() + nghost;
1033 }
1034 // sum up all contributions of element matrix
1035 apply(resultView, I_nd) += A_K_diag[i];
1036 }
1037 });
1038
1039 if (bcType == PERIODIC_FACE) {
1040 resultField.accumulateHalo();
1041 bcField.apply(resultField);
1042 bcField.assignGhostToPhysical(resultField);
1043 } else {
1044 resultField.accumulateHalo_noghost();
1045 }
1046
1047 // apply the inverse diagonal after already summed all contributions from element matrices
1048 using index_array_type = typename RangePolicy<Dim, exec_space>::index_array_type;
1049 ippl::parallel_for("Loop over result view to apply inverse", field.getFieldRangePolicy(),
1050 KOKKOS_LAMBDA(const index_array_type& args) {
1051 if (apply(resultView, args) != 0.0) {
1052 apply(resultView, args) = (1.0 / apply(resultView, args)) * apply(view, args);
1053 }
1054 });
1055 IpplTimings::stopTimer(outer_loop);
1056
1057 IpplTimings::stopTimer(evalAx);
1058
1059 return resultField;
1060 }
1061
1062 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
1063 typename QuadratureType, typename FieldLHS, typename FieldRHS>
1064 template <typename F>
1065 FieldLHS LagrangeSpace<T, Dim, Order, ElementType, QuadratureType, FieldLHS,
1066 FieldRHS>::evaluateAx_diag(FieldLHS& field, F& evalFunction) const {
1067 Inform m("");
1068
1069 // start a timer
1070 static IpplTimings::TimerRef evalAx = IpplTimings::getTimer("evaluateAx");
1072
1073 // get number of ghost cells in field
1074 const int nghost = field.getNghost();
1075
1076 // create a new field for result with view initialized to zero (views are initialized to
1077 // zero by default)
1078 FieldLHS resultField(field.get_mesh(), field.getLayout(), nghost);
1079
1080 // List of quadrature weights
1082 this->quadrature_m.getWeightsForRefElement();
1083
1084 // List of quadrature nodes
1086 this->quadrature_m.getIntegrationNodesForRefElement();
1087
1088 // TODO move outside of evaluateAx (I think it is possible for other problems as well)
1089 // Gradients of the basis functions for the DOF at the quadrature nodes
1090 Vector<Vector<point_t, numElementDOFs>, QuadratureType::numElementNodes> grad_b_q;
1091 for (size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
1092 for (size_t i = 0; i < numElementDOFs; ++i) {
1093 grad_b_q[k][i] = this->evaluateRefElementShapeFunctionGradient(i, q[k]);
1094 }
1095 }
1096
1097 // Make local element matrix -- does not change through the element mesh
1098 // Element matrix
1100
1101 // 1. Compute the Galerkin element matrix A_K
1102 for (size_t i = 0; i < numElementDOFs; ++i) {
1103 A_K_diag[i] = 0.0;
1104 for (size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
1105 A_K_diag[i] += w[k] * evalFunction(i, i, grad_b_q[k]);
1106 }
1107 }
1108
1109 // Get field data and atomic result data,
1110 // since it will be added to during the kokkos loop
1111 ViewType view = field.getView();
1112 AtomicViewType resultView = resultField.getView();
1113
1114 // Get boundary conditions from field
1115 BConds<FieldLHS, Dim>& bcField = field.getFieldBC();
1116 FieldBC bcType = bcField[0]->getBCType();
1117
1118 // Get domain information
1119 auto ldom = (field.getLayout()).getLocalNDIndex();
1120
1121 using exec_space = typename Kokkos::View<const size_t*>::execution_space;
1122 using policy_type = Kokkos::RangePolicy<exec_space>;
1123
1124 // start a timer
1125 static IpplTimings::TimerRef outer_loop = IpplTimings::getTimer("evaluateAx: outer loop");
1126 IpplTimings::startTimer(outer_loop);
1127
1128 // Loop over elements to compute contributions
1129 Kokkos::parallel_for(
1130 "Loop over elements", policy_type(0, elementIndices.extent(0)),
1131 KOKKOS_CLASS_LAMBDA(const size_t index) {
1132 const size_t elementIndex = elementIndices(index);
1133 const Vector<size_t, numElementDOFs> global_dofs =
1134 this->LagrangeSpace::getGlobalDOFIndices(elementIndex);
1135 Vector<indices_t, numElementDOFs> global_dof_ndindices;
1136
1137 for (size_t i = 0; i < numElementDOFs; ++i) {
1138 global_dof_ndindices[i] = this->getMeshVertexNDIndex(global_dofs[i]);
1139 }
1140
1141 // local DOF indices
1142 size_t i;
1143
1144 // global DOF n-dimensional indices (Vector of N indices representing indices in
1145 // each dimension)
1146 indices_t I_nd;
1147
1148 // 2. Compute the contribution to resultAx = A*x with A_K
1149 for (i = 0; i < numElementDOFs; ++i) {
1150 I_nd = global_dof_ndindices[i];
1151
1152 // Handle boundary DOFs
1153 // If Zero Dirichlet BCs, skip this DOF
1154 // If Constant Dirichlet BCs, identity
1155 if ((bcType == CONSTANT_FACE) && (this->isDOFOnBoundary(I_nd))) {
1156 for (unsigned d = 0; d < Dim; ++d) {
1157 I_nd[d] = I_nd[d] - ldom[d].first() + nghost;
1158 }
1159 apply(resultView, I_nd) = apply(view, I_nd);
1160 continue;
1161 } else if ((bcType == ZERO_FACE) && (this->isDOFOnBoundary(I_nd))) {
1162 continue;
1163 }
1164
1165 // get the appropriate index for the Kokkos view of the field
1166 for (unsigned d = 0; d < Dim; ++d) {
1167 I_nd[d] = I_nd[d] - ldom[d].first() + nghost;
1168 }
1169 apply(resultView, I_nd) += A_K_diag[i] * apply(view, I_nd);
1170 }
1171 });
1172 IpplTimings::stopTimer(outer_loop);
1173
1174 if (bcType == PERIODIC_FACE) {
1175 resultField.accumulateHalo();
1176 bcField.apply(resultField);
1177 bcField.assignGhostToPhysical(resultField);
1178 } else {
1179 resultField.accumulateHalo_noghost();
1180 }
1181
1182 IpplTimings::stopTimer(evalAx);
1183
1184 return resultField;
1185 }
1186
1187 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
1188 typename QuadratureType, typename FieldLHS, typename FieldRHS>
1189 template <typename F>
1190 FieldLHS LagrangeSpace<T, Dim, Order, ElementType, QuadratureType, FieldLHS,
1191 FieldRHS>::evaluateAx_lift(FieldLHS& field, F& evalFunction) const {
1192 Inform m("");
1193
1194 // get number of ghost cells in field
1195 const int nghost = field.getNghost();
1196
1197 // create a new field for result with view initialized to zero (views are initialized to
1198 // zero by default)
1199 FieldLHS resultField(field.get_mesh(), field.getLayout(), nghost);
1200
1201 // List of quadrature weights
1203 this->quadrature_m.getWeightsForRefElement();
1204
1205 // List of quadrature nodes
1207 this->quadrature_m.getIntegrationNodesForRefElement();
1208
1209 // TODO move outside of evaluateAx (I think it is possible for other problems as well)
1210 // Gradients of the basis functions for the DOF at the quadrature nodes
1211 Vector<Vector<point_t, numElementDOFs>, QuadratureType::numElementNodes> grad_b_q;
1212 for (size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
1213 for (size_t i = 0; i < numElementDOFs; ++i) {
1214 grad_b_q[k][i] = this->evaluateRefElementShapeFunctionGradient(i, q[k]);
1215 }
1216 }
1217
1218 // Make local element matrix -- does not change through the element mesh
1219 // Element matrix
1221
1222 // 1. Compute the Galerkin element matrix A_K
1223 for (size_t i = 0; i < numElementDOFs; ++i) {
1224 for (size_t j = 0; j < numElementDOFs; ++j) {
1225 A_K[i][j] = 0.0;
1226 for (size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
1227 A_K[i][j] += w[k] * evalFunction(i, j, grad_b_q[k]);
1228 }
1229 }
1230 }
1231
1232 // Get field data and atomic result data,
1233 // since it will be added to during the kokkos loop
1234 ViewType view = field.getView();
1235 AtomicViewType resultView = resultField.getView();
1236
1237 // Get domain information
1238 auto ldom = (field.getLayout()).getLocalNDIndex();
1239
1240 using exec_space = typename Kokkos::View<const size_t*>::execution_space;
1241 using policy_type = Kokkos::RangePolicy<exec_space>;
1242
1243 // Loop over elements to compute contributions
1244 Kokkos::parallel_for(
1245 "Loop over elements", policy_type(0, elementIndices.extent(0)),
1246 KOKKOS_CLASS_LAMBDA(const size_t index) {
1247 const size_t elementIndex = elementIndices(index);
1248 const Vector<size_t, numElementDOFs> global_dofs =
1249 this->LagrangeSpace::getGlobalDOFIndices(elementIndex);
1250 Vector<indices_t, numElementDOFs> global_dof_ndindices;
1251
1252 for (size_t i = 0; i < numElementDOFs; ++i) {
1253 global_dof_ndindices[i] = this->getMeshVertexNDIndex(global_dofs[i]);
1254 }
1255
1256 // local DOF indices (both i and j go from 0 to numDOFs-1 in the element)
1257 size_t i, j;
1258
1259 // global DOF n-dimensional indices (Vector of N indices representing indices in
1260 // each dimension)
1261 indices_t I_nd, J_nd;
1262
1263 // 2. Compute the contribution to resultAx = A*x with A_K
1264 for (i = 0; i < numElementDOFs; ++i) {
1265 I_nd = global_dof_ndindices[i];
1266
1267 // Skip if on a row of the matrix
1268 if (this->isDOFOnBoundary(I_nd)) {
1269 continue;
1270 }
1271
1272 // get the appropriate index for the Kokkos view of the field
1273 for (unsigned d = 0; d < Dim; ++d) {
1274 I_nd[d] = I_nd[d] - ldom[d].first() + nghost;
1275 }
1276
1277 for (j = 0; j < numElementDOFs; ++j) {
1278 J_nd = global_dof_ndindices[j];
1279
1280 // Contribute to lifting only if on a boundary DOF
1281 if (this->isDOFOnBoundary(J_nd)) {
1282 // get the appropriate index for the Kokkos view of the field
1283 for (unsigned d = 0; d < Dim; ++d) {
1284 J_nd[d] = J_nd[d] - ldom[d].first() + nghost;
1285 }
1286 apply(resultView, I_nd) += A_K[i][j] * apply(view, J_nd);
1287 continue;
1288 }
1289
1290 }
1291 }
1292 });
1293 resultField.accumulateHalo();
1294
1295 return resultField;
1296 }
1297
1298 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
1299 typename QuadratureType, typename FieldLHS, typename FieldRHS>
1300 void LagrangeSpace<T, Dim, Order, ElementType, QuadratureType, FieldLHS,
1301 FieldRHS>::evaluateLoadVector(FieldRHS& field) const {
1302 Inform m("");
1303
1304 // start a timer
1305 static IpplTimings::TimerRef evalLoadV = IpplTimings::getTimer("evaluateLoadVector");
1306 IpplTimings::startTimer(evalLoadV);
1307
1308 // List of quadrature weights
1310 this->quadrature_m.getWeightsForRefElement();
1311
1312 // List of quadrature nodes
1314 this->quadrature_m.getIntegrationNodesForRefElement();
1315
1316 const indices_t zeroNdIndex = Vector<size_t, Dim>(0);
1317
1318 // Evaluate the basis functions for the DOF at the quadrature nodes
1319 Vector<Vector<T, numElementDOFs>, QuadratureType::numElementNodes> basis_q;
1320 for (size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
1321 for (size_t i = 0; i < numElementDOFs; ++i) {
1322 basis_q[k][i] = this->evaluateRefElementShapeFunction(i, q[k]);
1323 }
1324 }
1325
1326 // Absolute value of det Phi_K
1327 const T absDetDPhi = Kokkos::abs(this->ref_element_m.getDeterminantOfTransformationJacobian(
1328 this->getElementMeshVertexPoints(zeroNdIndex)));
1329
1330 // Get domain information and ghost cells
1331 auto ldom = (field.getLayout()).getLocalNDIndex();
1332 const int nghost = field.getNghost();
1333
1334 // Get boundary conditions from field
1335 BConds<FieldRHS, Dim>& bcField = field.getFieldBC();
1336 FieldBC bcType = bcField[0]->getBCType();
1337
1338 FieldRHS temp_field(field.get_mesh(), field.getLayout(), nghost);
1339 temp_field.setFieldBC(bcField);
1340
1341 // Get field data and make it atomic,
1342 // since it will be added to during the kokkos loop
1343 // We work with a temporary field since we need to use field
1344 // to evaluate the load vector; then we assign temp to RHS field
1345 AtomicViewType atomic_view = temp_field.getView();
1346
1347 using exec_space = typename Kokkos::View<const size_t*>::execution_space;
1348 using policy_type = Kokkos::RangePolicy<exec_space>;
1349
1350 // start a timer
1351 static IpplTimings::TimerRef outer_loop =
1352 IpplTimings::getTimer("evaluateLoadVec: outer loop");
1353 IpplTimings::startTimer(outer_loop);
1354
1355 // Loop over elements to compute contributions
1356 Kokkos::parallel_for(
1357 "Loop over elements", policy_type(0, elementIndices.extent(0)),
1358 KOKKOS_CLASS_LAMBDA(size_t index) {
1359 const size_t elementIndex = elementIndices(index);
1360 const Vector<size_t, numElementDOFs> global_dofs =
1361 this->LagrangeSpace::getGlobalDOFIndices(elementIndex);
1362
1363 size_t i, I;
1364
1365 // 1. Compute b_K
1366 for (i = 0; i < numElementDOFs; ++i) {
1367 I = global_dofs[i];
1368
1369 // TODO fix for higher order
1370 auto dof_ndindex_I = this->getMeshVertexNDIndex(I);
1371
1372 // Skip boundary DOFs (Zero and Constant Dirichlet BCs)
1373 if (((bcType == ZERO_FACE) || (bcType == CONSTANT_FACE))
1374 && (this->isDOFOnBoundary(dof_ndindex_I))) {
1375 continue;
1376 }
1377
1378 // calculate the contribution of this element
1379 T contrib = 0;
1380 for (size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
1381 T val = 0;
1382 for (size_t j = 0; j < numElementDOFs; ++j) {
1383 // get field index corresponding to this DOF
1384 size_t J = global_dofs[j];
1385 auto dof_ndindex_J = this->getMeshVertexNDIndex(J);
1386 for (unsigned d = 0; d < Dim; ++d) {
1387 dof_ndindex_J[d] = dof_ndindex_J[d] - ldom[d].first() + nghost;
1388 }
1389
1390 // get field value at DOF and interpolate to q_k
1391 val += basis_q[k][j] * apply(field, dof_ndindex_J);
1392 }
1393
1394 contrib += w[k] * basis_q[k][i] * absDetDPhi * val;
1395 }
1396
1397 // get the appropriate index for the Kokkos view of the field
1398 for (unsigned d = 0; d < Dim; ++d) {
1399 dof_ndindex_I[d] = dof_ndindex_I[d] - ldom[d].first() + nghost;
1400 }
1401
1402 // add the contribution of the element to the field
1403 apply(atomic_view, dof_ndindex_I) += contrib;
1404
1405 }
1406 });
1407 IpplTimings::stopTimer(outer_loop);
1408
1409 temp_field.accumulateHalo();
1410
1411 if ((bcType == PERIODIC_FACE) || (bcType == CONSTANT_FACE)) {
1412 bcField.apply(temp_field);
1413 bcField.assignGhostToPhysical(temp_field);
1414 }
1415
1416 field = temp_field;
1417
1418 IpplTimings::stopTimer(evalLoadV);
1419 }
1420
1424
1425 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
1426 typename QuadratureType, typename FieldLHS, typename FieldRHS>
1427 template <typename F>
1429 const FieldLHS& u_h, const F& u_sol) const {
1430 if (this->quadrature_m.getOrder() < (2 * Order + 1)) {
1431 // throw exception
1432 throw IpplException(
1433 "LagrangeSpace::computeErrorL2()",
1434 "Order of quadrature rule for error computation should be > 2*p + 1");
1435 }
1436
1437 // List of quadrature weights
1439 this->quadrature_m.getWeightsForRefElement();
1440
1441 // List of quadrature nodes
1443 this->quadrature_m.getIntegrationNodesForRefElement();
1444
1445 // Evaluate the basis functions for the DOF at the quadrature nodes
1446 Vector<Vector<T, numElementDOFs>, QuadratureType::numElementNodes> basis_q;
1447 for (size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
1448 for (size_t i = 0; i < numElementDOFs; ++i) {
1449 basis_q[k][i] = this->evaluateRefElementShapeFunction(i, q[k]);
1450 }
1451 }
1452
1453 const indices_t zeroNdIndex = Vector<size_t, Dim>(0);
1454
1455 // Absolute value of det Phi_K
1456 const T absDetDPhi = Kokkos::abs(this->ref_element_m.getDeterminantOfTransformationJacobian(
1457 this->getElementMeshVertexPoints(zeroNdIndex)));
1458
1459 // Variable to sum the error to
1460 T error = 0;
1461
1462 // Get domain information and ghost cells
1463 auto ldom = (u_h.getLayout()).getLocalNDIndex();
1464 const int nghost = u_h.getNghost();
1465
1466 using exec_space = typename Kokkos::View<const size_t*>::execution_space;
1467 using policy_type = Kokkos::RangePolicy<exec_space>;
1468
1469 // Loop over elements to compute contributions
1470 Kokkos::parallel_reduce(
1471 "Compute error over elements", policy_type(0, elementIndices.extent(0)),
1472 KOKKOS_CLASS_LAMBDA(size_t index, double& local) {
1473 const size_t elementIndex = elementIndices(index);
1474 const Vector<size_t, numElementDOFs> global_dofs =
1475 this->LagrangeSpace::getGlobalDOFIndices(elementIndex);
1476
1477 // contribution of this element to the error
1478 T contrib = 0;
1479 for (size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
1480 T val_u_sol = u_sol(this->ref_element_m.localToGlobal(
1481 this->getElementMeshVertexPoints(this->getElementNDIndex(elementIndex)),
1482 q[k]));
1483
1484 T val_u_h = 0;
1485 for (size_t i = 0; i < numElementDOFs; ++i) {
1486 // get field index corresponding to this DOF
1487 size_t I = global_dofs[i];
1488 auto dof_ndindex_I = this->getMeshVertexNDIndex(I);
1489 for (unsigned d = 0; d < Dim; ++d) {
1490 dof_ndindex_I[d] = dof_ndindex_I[d] - ldom[d].first() + nghost;
1491 }
1492
1493 // get field value at DOF and interpolate to q_k
1494 val_u_h += basis_q[k][i] * apply(u_h, dof_ndindex_I);
1495 }
1496
1497 contrib += w[k] * Kokkos::pow(val_u_sol - val_u_h, 2) * absDetDPhi;
1498 }
1499 local += contrib;
1500 },
1501 Kokkos::Sum<double>(error));
1502
1503 // MPI reduce
1504 T global_error = 0.0;
1505 Comm->allreduce(error, global_error, 1, std::plus<T>());
1506
1507 return Kokkos::sqrt(global_error);
1508 }
1509
1510 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
1511 typename QuadratureType, typename FieldLHS, typename FieldRHS>
1513 const FieldLHS& u_h) const {
1514 if (this->quadrature_m.getOrder() < (2 * Order + 1)) {
1515 // throw exception
1516 throw IpplException(
1517 "LagrangeSpace::computeAvg()",
1518 "Order of quadrature rule for error computation should be > 2*p + 1");
1519 }
1520
1521 // List of quadrature weights
1523 this->quadrature_m.getWeightsForRefElement();
1524
1525 // List of quadrature nodes
1527 this->quadrature_m.getIntegrationNodesForRefElement();
1528
1529 // Evaluate the basis functions for the DOF at the quadrature nodes
1530 Vector<Vector<T, numElementDOFs>, QuadratureType::numElementNodes> basis_q;
1531 for (size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
1532 for (size_t i = 0; i < numElementDOFs; ++i) {
1533 basis_q[k][i] = this->evaluateRefElementShapeFunction(i, q[k]);
1534 }
1535 }
1536
1537 const indices_t zeroNdIndex = Vector<size_t, Dim>(0);
1538
1539 // Absolute value of det Phi_K
1540 const T absDetDPhi = Kokkos::abs(this->ref_element_m.getDeterminantOfTransformationJacobian(
1541 this->getElementMeshVertexPoints(zeroNdIndex)));
1542
1543 // Variable to sum the error to
1544 T avg = 0;
1545
1546 // Get domain information and ghost cells
1547 auto ldom = (u_h.getLayout()).getLocalNDIndex();
1548 const int nghost = u_h.getNghost();
1549
1550 using exec_space = typename Kokkos::View<const size_t*>::execution_space;
1551 using policy_type = Kokkos::RangePolicy<exec_space>;
1552
1553 // Loop over elements to compute contributions
1554 Kokkos::parallel_reduce(
1555 "Compute average over elements", policy_type(0, elementIndices.extent(0)),
1556 KOKKOS_CLASS_LAMBDA(size_t index, double& local) {
1557 const size_t elementIndex = elementIndices(index);
1558 const Vector<size_t, numElementDOFs> global_dofs =
1559 this->LagrangeSpace::getGlobalDOFIndices(elementIndex);
1560
1561 // contribution of this element to the error
1562 T contrib = 0;
1563 for (size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
1564 T val_u_h = 0;
1565 for (size_t i = 0; i < numElementDOFs; ++i) {
1566 // get field index corresponding to this DOF
1567 size_t I = global_dofs[i];
1568 auto dof_ndindex_I = this->getMeshVertexNDIndex(I);
1569 for (unsigned d = 0; d < Dim; ++d) {
1570 dof_ndindex_I[d] = dof_ndindex_I[d] - ldom[d].first() + nghost;
1571 }
1572
1573 // get field value at DOF and interpolate to q_k
1574 val_u_h += basis_q[k][i] * apply(u_h, dof_ndindex_I);
1575 }
1576
1577 contrib += w[k] * val_u_h * absDetDPhi;
1578 }
1579 local += contrib;
1580 },
1581 Kokkos::Sum<double>(avg));
1582
1583 // MPI reduce
1584 T global_avg = 0.0;
1585 Comm->allreduce(avg, global_avg, 1, std::plus<T>());
1586
1587 return global_avg;
1588 }
1589
1593
1594 // Function to return the device struct of this Lagrange Space
1595 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
1596 typename QuadratureType, typename FieldLHS, typename FieldRHS>
1598 DeviceStruct
1600 getDeviceMirror() const {
1601 DeviceStruct space_mirror;
1602 space_mirror.nr_m = this->nr_m;
1603 space_mirror.ref_element_m = this->ref_element_m;
1604 return space_mirror;
1605 }
1606
1607 // I don't know how to avoid code duplication here...
1608 // Make sure that any changes in getLocalDOFIndex, getGlobalDOFIndices,
1609 // evaluateRefElementShapeFunction, and getMeshVertexNDIndex from the
1610 // parent class FiniteElementSpace get propagated here.
1611
1612 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
1613 typename QuadratureType, typename FieldLHS, typename FieldRHS>
1614 KOKKOS_FUNCTION size_t
1616 DeviceStruct::getLocalDOFIndex(const indices_t& elementNDIndex,
1617 const size_t& globalDOFIndex) const {
1618
1619 static_assert(Dim == 1 || Dim == 2 || Dim == 3, "Dim must be 1, 2 or 3");
1620 // TODO fix not order independent, only works for order 1
1621 static_assert(Order == 1, "Only order 1 is supported at the moment");
1622
1623 // Get all the global DOFs for the element
1624 const Vector<size_t, numElementDOFs> global_dofs =
1625 this->getGlobalDOFIndices(elementNDIndex);
1626
1627 // Find the global DOF in the vector and return the local DOF index
1628 // Note: It is important that this only works because the global_dofs
1629 // are already arranged in the correct order from getGlobalDOFIndices
1630 for (size_t i = 0; i < global_dofs.dim; ++i) {
1631 if (global_dofs[i] == globalDOFIndex) {
1632 return i;
1633 }
1634 }
1635 // commented this due to this being on device
1636 // however, it would be good to throw an error in this case
1637 //throw IpplException("LagrangeSpace::getLocalDOFIndex()",
1638 // "FEM Lagrange Space: Global DOF not found in specified element");
1639 return 0;
1640 }
1641
1642 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
1643 typename QuadratureType, typename FieldLHS, typename FieldRHS>
1644 KOKKOS_FUNCTION Vector<size_t, LagrangeSpace<T, Dim, Order, ElementType, QuadratureType,
1645 FieldLHS, FieldRHS>::DeviceStruct::numElementDOFs>
1647 DeviceStruct::getGlobalDOFIndices(const indices_t& elementNDIndex) const {
1648
1649 Vector<size_t, numElementDOFs> globalDOFs(0);
1650
1651 // Compute the vector to multiply the ndindex with
1653 for (size_t d = 1; d < dim; ++d) {
1654 for (size_t d2 = d; d2 < Dim; ++d2) {
1655 vec[d2] *= this->nr_m[d - 1];
1656 }
1657 }
1658 vec *= Order; // Multiply each dimension by the order
1659 size_t smallestGlobalDOF = elementNDIndex.dot(vec);
1660
1661 // Add the vertex DOFs
1662 globalDOFs[0] = smallestGlobalDOF;
1663 globalDOFs[1] = smallestGlobalDOF + Order;
1664
1665 if constexpr (Dim >= 2) {
1666 globalDOFs[2] = globalDOFs[1] + this->nr_m[0] * Order;
1667 globalDOFs[3] = globalDOFs[0] + this->nr_m[0] * Order;
1668 }
1669 if constexpr (Dim >= 3) {
1670 globalDOFs[4] = globalDOFs[0] + this->nr_m[1] * this->nr_m[0] * Order;
1671 globalDOFs[5] = globalDOFs[1] + this->nr_m[1] * this->nr_m[0] * Order;
1672 globalDOFs[6] = globalDOFs[2] + this->nr_m[1] * this->nr_m[0] * Order;
1673 globalDOFs[7] = globalDOFs[3] + this->nr_m[1] * this->nr_m[0] * Order;
1674 }
1675
1676 if constexpr (Order > 1) {
1677 // If the order is greater than 1, there are edge and face DOFs, otherwise the work is
1678 // done
1679
1680 // Add the edge DOFs
1681 if constexpr (Dim >= 2) {
1682 for (size_t i = 0; i < Order - 1; ++i) {
1683 globalDOFs[8 + i] = globalDOFs[0] + i + 1;
1684 globalDOFs[8 + Order - 1 + i] = globalDOFs[1] + (i + 1) * this->nr_m[1];
1685 globalDOFs[8 + 2 * (Order - 1) + i] = globalDOFs[2] - (i + 1);
1686 globalDOFs[8 + 3 * (Order - 1) + i] = globalDOFs[3] - (i + 1) * this->nr_m[1];
1687 }
1688 }
1689 if constexpr (Dim >= 3) {
1690 // TODO
1691 }
1692
1693 // Add the face DOFs
1694 if constexpr (Dim >= 2) {
1695 for (size_t i = 0; i < Order - 1; ++i) {
1696 for (size_t j = 0; j < Order - 1; ++j) {
1697 // TODO CHECK
1698 globalDOFs[8 + 4 * (Order - 1) + i * (Order - 1) + j] =
1699 globalDOFs[0] + (i + 1) + (j + 1) * this->nr_m[1];
1700 globalDOFs[8 + 4 * (Order - 1) + (Order - 1) * (Order - 1) + i * (Order - 1)
1701 + j] = globalDOFs[1] + (i + 1) + (j + 1) * this->nr_m[1];
1702 globalDOFs[8 + 4 * (Order - 1) + 2 * (Order - 1) * (Order - 1)
1703 + i * (Order - 1) + j] =
1704 globalDOFs[2] - (i + 1) + (j + 1) * this->nr_m[1];
1705 globalDOFs[8 + 4 * (Order - 1) + 3 * (Order - 1) * (Order - 1)
1706 + i * (Order - 1) + j] =
1707 globalDOFs[3] - (i + 1) + (j + 1) * this->nr_m[1];
1708 }
1709 }
1710 }
1711 }
1712
1713 return globalDOFs;
1714 }
1715
1716 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
1717 typename QuadratureType, typename FieldLHS, typename FieldRHS>
1718 KOKKOS_FUNCTION T
1721 const LagrangeSpace<T, Dim, Order, ElementType, QuadratureType, FieldLHS,
1722 FieldRHS>::point_t& localPoint) const {
1723
1724 static_assert(Order == 1, "Only order 1 is supported at the moment");
1725 // Assert that the local vertex index is valid.
1726 assert(localDOF < DeviceStruct::numElementDOFs
1727 && "The local vertex index is invalid"); // TODO assumes 1st order Lagrange
1728
1729 assert(this->ref_element_m.isPointInRefElement(localPoint)
1730 && "Point is not in reference element");
1731
1732 // Get the local vertex indices for the local vertex index.
1733 // TODO fix not order independent, only works for order 1
1734 const point_t ref_element_point = this->ref_element_m.getLocalVertices()[localDOF];
1735
1736 // The variable that accumulates the product of the shape functions.
1737 T product = 1;
1738
1739 for (size_t d = 0; d < Dim; d++) {
1740 if (localPoint[d] < ref_element_point[d]) {
1741 product *= localPoint[d];
1742 } else {
1743 product *= 1.0 - localPoint[d];
1744 }
1745 }
1746
1747 return product;
1748 }
1749
1750 template <typename T, unsigned Dim, unsigned Order, typename ElementType,
1751 typename QuadratureType, typename FieldLHS, typename FieldRHS>
1752 KOKKOS_FUNCTION typename LagrangeSpace<T, Dim, Order, ElementType, QuadratureType,
1753 FieldLHS, FieldRHS>::indices_t
1755 DeviceStruct::getMeshVertexNDIndex(const size_t& vertex_index) const {
1756 // Copy the vertex index to the index variable we can alter during the computation.
1757 size_t index = vertex_index;
1758
1759 // Create a vector to store the vertex indices in each dimension for the corresponding
1760 // vertex.
1761 indices_t vertex_indices;
1762
1763 // This is the number of vertices in each dimension.
1764 Vector<size_t, Dim> vertices_per_dim = nr_m;
1765
1766 // The number_of_lower_dim_vertices is the product of the number of vertices per
1767 // dimension, it will get divided by the current dimensions number to get the index in
1768 // that dimension
1769 size_t remaining_number_of_vertices = 1;
1770 for (const size_t num_vertices : vertices_per_dim) {
1771 remaining_number_of_vertices *= num_vertices;
1772 }
1773
1774 for (int d = Dim - 1; d >= 0; --d) {
1775 remaining_number_of_vertices /= vertices_per_dim[d];
1776 vertex_indices[d] = index / remaining_number_of_vertices;
1777 index -= vertex_indices[d] * remaining_number_of_vertices;
1778 }
1779
1780 return vertex_indices;
1781 };
1782
1783} // namespace ippl
constexpr unsigned Dim
constexpr unsigned getLagrangeNumElementDOFs(unsigned Dim, unsigned Order)
constexpr KOKKOS_INLINE_FUNCTION auto first()
Definition AbsorbingBC.h:10
Definition Archive.h:20
KOKKOS_INLINE_FUNCTION constexpr decltype(auto) apply(const View &view, const Coords &coords)
void parallel_for(const std::string &name, const ExecPolicy &policy, const FunctorType &functor)
FieldBC
Definition BcTypes.h:33
@ ZERO_FACE
Definition BcTypes.h:36
@ PERIODIC_FACE
Definition BcTypes.h:34
@ CONSTANT_FACE
Definition BcTypes.h:35
std::unique_ptr< mpi::Communicator > Comm
Definition Ippl.h:22
FiniteElementSpace(UniformCartesian< T, Dim > &mesh, ElementType &ref_element, const QuadratureType &quadrature)
KOKKOS_FUNCTION indices_t getElementNDIndex(const size_t &elementIndex) const
A class representing a Lagrange space for finite element methods on a structured, rectilinear grid.
FiniteElementSpace< Tlhs, Dim, numElementDOFs, ElementType, QuadratureType, FieldLHS, FieldLHS >::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.
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 ///////////////////////////////////////////////.
static constexpr unsigned numElementDOFs
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...
Kokkos::View< size_t * > elementIndices
FieldLHS evaluateAx_diag(FieldLHS &field, F &evalFunction) const
FiniteElementSpace< Tlhs, Dim, numElementDOFs, ElementType, QuadratureType, FieldLHS, FieldLHS >::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
void assignGhostToPhysical(Field &field)
Definition BConds.hpp:38
void apply(Field &field)
Definition BConds.hpp:29
const NDIndex_t & getLocalNDIndex() const
KOKKOS_INLINE_FUNCTION unsigned size() const noexcept
Definition NDIndex.hpp:43
static constexpr unsigned dim
Definition Vector.h:26
KOKKOS_INLINE_FUNCTION T dot(const Vector< T, Dim > &rhs) const
Definition Vector.hpp:182
Timing::TimerRef TimerRef
static TimerRef getTimer(const char *nm)
static void stopTimer(TimerRef t)
static void startTimer(TimerRef t)
::ippl::Vector< index_type, Dim > index_array_type