6 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
7 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
12 QuadratureType, FieldLHS, FieldRHS>(mesh, ref_element, quadrature) {
14 static_assert(
Dim >= 1 &&
Dim <= 3,
15 "Finite Element space only supports 1D, 2D and 3D meshes");
22 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
23 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
27 QuadratureType, FieldLHS, FieldRHS>(mesh, ref_element, quadrature) {
29 static_assert(
Dim >= 1 &&
Dim <= 3,
30 "Finite Element space only supports 1D, 2D and 3D meshes");
36 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
37 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
42 QuadratureType, FieldLHS, FieldRHS>::setMesh(mesh);
45 initializeElementIndices(layout);
49 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
50 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
52 FieldRHS>::initializeElementIndices(
const Layout_t& layout) {
54 int npoints = ldom.
size();
55 auto first = ldom.first();
56 auto last = ldom.last();
59 for (
size_t d = 0; d <
Dim; ++d) {
60 bounds[d] = this->nr_m[d] - 1;
63 int upperBoundaryPoints = -1;
67 Kokkos::View<size_t*> points(
"npoints", npoints);
69 "ComputePoints", npoints,
70 KOKKOS_CLASS_LAMBDA(
const int i,
int& local) {
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);
78 if (val[d] == bounds[d]) {
82 points(i) = (!isBoundary) * (this->getElementIndex(val));
85 Kokkos::Sum<int>(upperBoundaryPoints));
90 int elementsPerRank = npoints - upperBoundaryPoints;
91 elementIndices = Kokkos::View<size_t*>(
"i", elementsPerRank);
92 Kokkos::View<size_t> index(
"index");
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);
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;
116 return num_global_dofs;
119 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
120 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
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");
126 static_assert(Order == 1,
"Only order 1 is supported at the moment");
135 for (
size_t i = 0; i < global_dofs.
dim; ++i) {
136 if (global_dofs[i] == globalDOFIndex) {
147 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
148 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
149 KOKKOS_FUNCTION
size_t
151 FieldRHS>::getGlobalDOFIndex(
const size_t& elementIndex,
152 const size_t& localDOFIndex)
const {
155 return global_dofs[localDOFIndex];
158 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
159 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
161 FieldLHS, FieldRHS>::numElementDOFs>
163 FieldRHS>::getLocalDOFIndices()
const {
166 for (
size_t dof = 0; dof < numElementDOFs; ++dof) {
167 localDOFs[dof] = dof;
173 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
174 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
176 FieldLHS, FieldRHS>::numElementDOFs>
178 FieldRHS>::getGlobalDOFIndices(
const size_t& elementIndex)
const {
182 indices_t elementPos = this->getElementNDIndex(elementIndex);
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];
192 size_t smallestGlobalDOF = elementPos.
dot(vec);
195 globalDOFs[0] = smallestGlobalDOF;
196 globalDOFs[1] = smallestGlobalDOF + Order;
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;
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;
209 if constexpr (Order > 1) {
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];
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) {
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];
253 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
254 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
258 const size_t& localDOF,
260 FieldRHS>::point_t& localPoint)
const {
261 static_assert(Order == 1,
"Only order 1 is supported at the moment");
263 assert(localDOF < numElementDOFs
264 &&
"The local vertex index is invalid");
266 assert(this->ref_element_m.isPointInRefElement(localPoint)
267 &&
"Point is not in reference element");
271 const point_t ref_element_point = this->ref_element_m.getLocalVertices()[localDOF];
276 for (
size_t d = 0; d <
Dim; d++) {
277 if (localPoint[d] < ref_element_point[d]) {
278 product *= localPoint[d];
280 product *= 1.0 - localPoint[d];
287 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
288 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
293 const size_t& localDOF,
295 FieldRHS>::point_t& localPoint)
const {
297 static_assert(Order == 1 &&
"Only order 1 is supported at the moment");
300 assert(localDOF < numElementDOFs &&
"The local vertex index is invalid");
302 assert(this->ref_element_m.isPointInRefElement(localPoint)
303 &&
"Point is not in reference element");
306 const vertex_points_t local_vertex_points = this->ref_element_m.getLocalVertices();
308 const point_t& local_vertex_point = local_vertex_points[localDOF];
316 for (
size_t d = 0; d <
Dim; d++) {
320 for (
size_t d2 = 0; d2 <
Dim; d2++) {
322 if (localPoint[d] < local_vertex_point[d]) {
328 if (localPoint[d2] < local_vertex_point[d2]) {
329 product *= localPoint[d2];
331 product *= 1.0 - localPoint[d2];
336 gradient[d] = product;
346 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
347 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
348 template <
typename F>
350 FieldRHS>::evaluateAx(FieldLHS& field, F& evalFunction)
const {
358 const int nghost = field.getNghost();
362 FieldLHS resultField(field.get_mesh(), field.getLayout(), nghost);
366 this->quadrature_m.getWeightsForRefElement();
370 this->quadrature_m.getIntegrationNodesForRefElement();
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]);
386 for (
size_t i = 0; i < numElementDOFs; ++i) {
387 for (
size_t j = 0; j < numElementDOFs; ++j) {
389 for (
size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
390 A_K[i][j] += w[k] * evalFunction(i, j, grad_b_q[k]);
402 FieldBC bcType = bcField[0]->getBCType();
405 auto ldom = (field.getLayout()).getLocalNDIndex();
407 using exec_space =
typename Kokkos::View<const size_t*>::execution_space;
408 using policy_type = Kokkos::RangePolicy<exec_space>;
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);
423 for (
size_t i = 0; i < numElementDOFs; ++i) {
424 global_dof_ndindices[i] = this->getMeshVertexNDIndex(global_dofs[i]);
435 for (i = 0; i < numElementDOFs; ++i) {
436 I_nd = global_dof_ndindices[i];
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;
447 }
else if ((bcType ==
ZERO_FACE) && (this->isDOFOnBoundary(I_nd))) {
452 for (
unsigned d = 0; d <
Dim; ++d) {
453 I_nd[d] = I_nd[d] - ldom[d].first() + nghost;
456 for (j = 0; j < numElementDOFs; ++j) {
457 J_nd = global_dof_ndindices[j];
461 && this->isDOFOnBoundary(J_nd)) {
466 for (
unsigned d = 0; d <
Dim; ++d) {
467 J_nd[d] = J_nd[d] - ldom[d].first() + nghost;
470 apply(resultView, I_nd) += A_K[i][j] *
apply(view, J_nd);
477 resultField.accumulateHalo();
478 bcField.
apply(resultField);
481 resultField.accumulateHalo_noghost();
489 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
490 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
491 template <
typename F>
493 FieldRHS>::evaluateAx_lower(FieldLHS& field, F& evalFunction)
const {
501 const int nghost = field.getNghost();
505 FieldLHS resultField(field.get_mesh(), field.getLayout(), nghost);
509 this->quadrature_m.getWeightsForRefElement();
513 this->quadrature_m.getIntegrationNodesForRefElement();
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]);
529 for (
size_t i = 0; i < numElementDOFs; ++i) {
530 for (
size_t j = 0; j < numElementDOFs; ++j) {
532 for (
size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
533 A_K[i][j] += w[k] * evalFunction(i, j, grad_b_q[k]);
545 FieldBC bcType = bcField[0]->getBCType();
548 auto ldom = (field.getLayout()).getLocalNDIndex();
550 using exec_space =
typename Kokkos::View<const size_t*>::execution_space;
551 using policy_type = Kokkos::RangePolicy<exec_space>;
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);
566 for (
size_t i = 0; i < numElementDOFs; ++i) {
567 global_dof_ndindices[i] = this->getMeshVertexNDIndex(global_dofs[i]);
578 for (i = 0; i < numElementDOFs; ++i) {
579 I_nd = global_dof_ndindices[i];
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;
590 }
else if ((bcType ==
ZERO_FACE) && (this->isDOFOnBoundary(I_nd))) {
595 for (
unsigned d = 0; d <
Dim; ++d) {
596 I_nd[d] = I_nd[d] - ldom[d].first() + nghost;
599 for (j = 0; j < numElementDOFs; ++j) {
600 J_nd = global_dof_ndindices[j];
602 if (global_dofs[i] >= global_dofs[j]) {
608 && this->isDOFOnBoundary(J_nd)) {
613 for (
unsigned d = 0; d <
Dim; ++d) {
614 J_nd[d] = J_nd[d] - ldom[d].first() + nghost;
617 apply(resultView, I_nd) += A_K[i][j] *
apply(view, J_nd);
624 resultField.accumulateHalo();
625 bcField.
apply(resultField);
628 resultField.accumulateHalo_noghost();
636 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
637 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
638 template <
typename F>
640 FieldRHS>::evaluateAx_upper(FieldLHS& field, F& evalFunction)
const {
648 const int nghost = field.getNghost();
652 FieldLHS resultField(field.get_mesh(), field.getLayout(), nghost);
656 this->quadrature_m.getWeightsForRefElement();
660 this->quadrature_m.getIntegrationNodesForRefElement();
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]);
676 for (
size_t i = 0; i < numElementDOFs; ++i) {
677 for (
size_t j = 0; j < numElementDOFs; ++j) {
679 for (
size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
680 A_K[i][j] += w[k] * evalFunction(i, j, grad_b_q[k]);
692 FieldBC bcType = bcField[0]->getBCType();
695 auto ldom = (field.getLayout()).getLocalNDIndex();
697 using exec_space =
typename Kokkos::View<const size_t*>::execution_space;
698 using policy_type = Kokkos::RangePolicy<exec_space>;
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);
713 for (
size_t i = 0; i < numElementDOFs; ++i) {
714 global_dof_ndindices[i] = this->getMeshVertexNDIndex(global_dofs[i]);
725 for (i = 0; i < numElementDOFs; ++i) {
726 I_nd = global_dof_ndindices[i];
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;
737 }
else if ((bcType ==
ZERO_FACE) && (this->isDOFOnBoundary(I_nd))) {
742 for (
unsigned d = 0; d <
Dim; ++d) {
743 I_nd[d] = I_nd[d] - ldom[d].first() + nghost;
746 for (j = 0; j < numElementDOFs; ++j) {
747 J_nd = global_dof_ndindices[j];
749 if (global_dofs[i] <= global_dofs[j]) {
755 && this->isDOFOnBoundary(J_nd)) {
760 for (
unsigned d = 0; d <
Dim; ++d) {
761 J_nd[d] = J_nd[d] - ldom[d].first() + nghost;
764 apply(resultView, I_nd) += A_K[i][j] *
apply(view, J_nd);
771 resultField.accumulateHalo();
772 bcField.
apply(resultField);
775 resultField.accumulateHalo_noghost();
783 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
784 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
785 template <
typename F>
787 FieldRHS>::evaluateAx_upperlower(FieldLHS& field, F& evalFunction)
const {
795 const int nghost = field.getNghost();
799 FieldLHS resultField(field.get_mesh(), field.getLayout(), nghost);
803 this->quadrature_m.getWeightsForRefElement();
807 this->quadrature_m.getIntegrationNodesForRefElement();
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]);
823 for (
size_t i = 0; i < numElementDOFs; ++i) {
824 for (
size_t j = 0; j < numElementDOFs; ++j) {
826 for (
size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
827 A_K[i][j] += w[k] * evalFunction(i, j, grad_b_q[k]);
839 FieldBC bcType = bcField[0]->getBCType();
842 auto ldom = (field.getLayout()).getLocalNDIndex();
844 using exec_space =
typename Kokkos::View<const size_t*>::execution_space;
845 using policy_type = Kokkos::RangePolicy<exec_space>;
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);
860 for (
size_t i = 0; i < numElementDOFs; ++i) {
861 global_dof_ndindices[i] = this->getMeshVertexNDIndex(global_dofs[i]);
872 for (i = 0; i < numElementDOFs; ++i) {
873 I_nd = global_dof_ndindices[i];
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;
884 }
else if ((bcType ==
ZERO_FACE) && (this->isDOFOnBoundary(I_nd))) {
889 for (
unsigned d = 0; d <
Dim; ++d) {
890 I_nd[d] = I_nd[d] - ldom[d].first() + nghost;
893 for (j = 0; j < i; ++j) {
894 J_nd = global_dof_ndindices[j];
898 && this->isDOFOnBoundary(J_nd)) {
903 for (
unsigned d = 0; d <
Dim; ++d) {
904 J_nd[d] = J_nd[d] - ldom[d].first() + nghost;
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);
915 resultField.accumulateHalo();
916 bcField.
apply(resultField);
919 resultField.accumulateHalo_noghost();
927 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
928 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
929 template <
typename F>
931 FieldRHS>::evaluateAx_inversediag(FieldLHS& field, F& evalFunction)
const {
939 const int nghost = field.getNghost();
943 FieldLHS resultField(field.get_mesh(), field.getLayout(), nghost);
947 this->quadrature_m.getWeightsForRefElement();
951 this->quadrature_m.getIntegrationNodesForRefElement();
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]);
967 for (
size_t i = 0; i < numElementDOFs; ++i) {
969 for (
size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
970 A_K_diag[i] += w[k] * evalFunction(i, i, grad_b_q[k]);
981 FieldBC bcType = bcField[0]->getBCType();
984 auto ldom = (field.getLayout()).getLocalNDIndex();
986 using exec_space =
typename Kokkos::View<const size_t*>::execution_space;
987 using policy_type = Kokkos::RangePolicy<exec_space>;
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);
1002 for (
size_t i = 0; i < numElementDOFs; ++i) {
1003 global_dof_ndindices[i] = this->getMeshVertexNDIndex(global_dofs[i]);
1014 for (i = 0; i < numElementDOFs; ++i) {
1015 I_nd = global_dof_ndindices[i];
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;
1024 apply(resultView, I_nd) = 1.0;
1026 }
else if ((bcType ==
ZERO_FACE) && (this->isDOFOnBoundary(I_nd))) {
1031 for (
unsigned d = 0; d <
Dim; ++d) {
1032 I_nd[d] = I_nd[d] - ldom[d].first() + nghost;
1035 apply(resultView, I_nd) += A_K_diag[i];
1040 resultField.accumulateHalo();
1041 bcField.
apply(resultField);
1044 resultField.accumulateHalo_noghost();
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);
1062 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
1063 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
1064 template <
typename F>
1066 FieldRHS>::evaluateAx_diag(FieldLHS& field, F& evalFunction)
const {
1074 const int nghost = field.getNghost();
1078 FieldLHS resultField(field.get_mesh(), field.getLayout(), nghost);
1082 this->quadrature_m.getWeightsForRefElement();
1086 this->quadrature_m.getIntegrationNodesForRefElement();
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]);
1102 for (
size_t i = 0; i < numElementDOFs; ++i) {
1104 for (
size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
1105 A_K_diag[i] += w[k] * evalFunction(i, i, grad_b_q[k]);
1116 FieldBC bcType = bcField[0]->getBCType();
1119 auto ldom = (field.getLayout()).getLocalNDIndex();
1121 using exec_space =
typename Kokkos::View<const size_t*>::execution_space;
1122 using policy_type = Kokkos::RangePolicy<exec_space>;
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);
1137 for (
size_t i = 0; i < numElementDOFs; ++i) {
1138 global_dof_ndindices[i] = this->getMeshVertexNDIndex(global_dofs[i]);
1149 for (i = 0; i < numElementDOFs; ++i) {
1150 I_nd = global_dof_ndindices[i];
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;
1161 }
else if ((bcType ==
ZERO_FACE) && (this->isDOFOnBoundary(I_nd))) {
1166 for (
unsigned d = 0; d <
Dim; ++d) {
1167 I_nd[d] = I_nd[d] - ldom[d].first() + nghost;
1169 apply(resultView, I_nd) += A_K_diag[i] *
apply(view, I_nd);
1175 resultField.accumulateHalo();
1176 bcField.
apply(resultField);
1179 resultField.accumulateHalo_noghost();
1187 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
1188 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
1189 template <
typename F>
1191 FieldRHS>::evaluateAx_lift(FieldLHS& field, F& evalFunction)
const {
1195 const int nghost = field.getNghost();
1199 FieldLHS resultField(field.get_mesh(), field.getLayout(), nghost);
1203 this->quadrature_m.getWeightsForRefElement();
1207 this->quadrature_m.getIntegrationNodesForRefElement();
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]);
1223 for (
size_t i = 0; i < numElementDOFs; ++i) {
1224 for (
size_t j = 0; j < numElementDOFs; ++j) {
1226 for (
size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
1227 A_K[i][j] += w[k] * evalFunction(i, j, grad_b_q[k]);
1238 auto ldom = (field.getLayout()).getLocalNDIndex();
1240 using exec_space =
typename Kokkos::View<const size_t*>::execution_space;
1241 using policy_type = Kokkos::RangePolicy<exec_space>;
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);
1252 for (
size_t i = 0; i < numElementDOFs; ++i) {
1253 global_dof_ndindices[i] = this->getMeshVertexNDIndex(global_dofs[i]);
1264 for (i = 0; i < numElementDOFs; ++i) {
1265 I_nd = global_dof_ndindices[i];
1268 if (this->isDOFOnBoundary(I_nd)) {
1273 for (
unsigned d = 0; d <
Dim; ++d) {
1274 I_nd[d] = I_nd[d] - ldom[d].first() + nghost;
1277 for (j = 0; j < numElementDOFs; ++j) {
1278 J_nd = global_dof_ndindices[j];
1281 if (this->isDOFOnBoundary(J_nd)) {
1283 for (
unsigned d = 0; d <
Dim; ++d) {
1284 J_nd[d] = J_nd[d] - ldom[d].first() + nghost;
1286 apply(resultView, I_nd) += A_K[i][j] *
apply(view, J_nd);
1293 resultField.accumulateHalo();
1298 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
1299 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
1301 FieldRHS>::evaluateLoadVector(FieldRHS& field)
const {
1310 this->quadrature_m.getWeightsForRefElement();
1314 this->quadrature_m.getIntegrationNodesForRefElement();
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]);
1327 const T absDetDPhi = Kokkos::abs(this->ref_element_m.getDeterminantOfTransformationJacobian(
1328 this->getElementMeshVertexPoints(zeroNdIndex)));
1331 auto ldom = (field.getLayout()).getLocalNDIndex();
1332 const int nghost = field.getNghost();
1336 FieldBC bcType = bcField[0]->getBCType();
1338 FieldRHS temp_field(field.get_mesh(), field.getLayout(), nghost);
1339 temp_field.setFieldBC(bcField);
1347 using exec_space =
typename Kokkos::View<const size_t*>::execution_space;
1348 using policy_type = Kokkos::RangePolicy<exec_space>;
1357 "Loop over elements", policy_type(0, elementIndices.extent(0)),
1358 KOKKOS_CLASS_LAMBDA(
size_t index) {
1359 const size_t elementIndex = elementIndices(index);
1366 for (i = 0; i < numElementDOFs; ++i) {
1370 auto dof_ndindex_I = this->getMeshVertexNDIndex(I);
1374 && (this->isDOFOnBoundary(dof_ndindex_I))) {
1380 for (
size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
1382 for (
size_t j = 0; j < numElementDOFs; ++j) {
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;
1391 val += basis_q[k][j] *
apply(field, dof_ndindex_J);
1394 contrib += w[k] * basis_q[k][i] * absDetDPhi * val;
1398 for (
unsigned d = 0; d <
Dim; ++d) {
1399 dof_ndindex_I[d] = dof_ndindex_I[d] - ldom[d].first() + nghost;
1403 apply(atomic_view, dof_ndindex_I) += contrib;
1409 temp_field.accumulateHalo();
1412 bcField.
apply(temp_field);
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)) {
1433 "LagrangeSpace::computeErrorL2()",
1434 "Order of quadrature rule for error computation should be > 2*p + 1");
1439 this->quadrature_m.getWeightsForRefElement();
1443 this->quadrature_m.getIntegrationNodesForRefElement();
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]);
1456 const T absDetDPhi = Kokkos::abs(this->ref_element_m.getDeterminantOfTransformationJacobian(
1457 this->getElementMeshVertexPoints(zeroNdIndex)));
1463 auto ldom = (u_h.getLayout()).getLocalNDIndex();
1464 const int nghost = u_h.getNghost();
1466 using exec_space =
typename Kokkos::View<const size_t*>::execution_space;
1467 using policy_type = Kokkos::RangePolicy<exec_space>;
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);
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)),
1485 for (
size_t i = 0; i < numElementDOFs; ++i) {
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;
1494 val_u_h += basis_q[k][i] *
apply(u_h, dof_ndindex_I);
1497 contrib += w[k] * Kokkos::pow(val_u_sol - val_u_h, 2) * absDetDPhi;
1501 Kokkos::Sum<double>(error));
1504 T global_error = 0.0;
1505 Comm->allreduce(error, global_error, 1, std::plus<T>());
1507 return Kokkos::sqrt(global_error);
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)) {
1517 "LagrangeSpace::computeAvg()",
1518 "Order of quadrature rule for error computation should be > 2*p + 1");
1523 this->quadrature_m.getWeightsForRefElement();
1527 this->quadrature_m.getIntegrationNodesForRefElement();
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]);
1540 const T absDetDPhi = Kokkos::abs(this->ref_element_m.getDeterminantOfTransformationJacobian(
1541 this->getElementMeshVertexPoints(zeroNdIndex)));
1547 auto ldom = (u_h.getLayout()).getLocalNDIndex();
1548 const int nghost = u_h.getNghost();
1550 using exec_space =
typename Kokkos::View<const size_t*>::execution_space;
1551 using policy_type = Kokkos::RangePolicy<exec_space>;
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);
1563 for (
size_t k = 0; k < QuadratureType::numElementNodes; ++k) {
1565 for (
size_t i = 0; i < numElementDOFs; ++i) {
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;
1574 val_u_h += basis_q[k][i] *
apply(u_h, dof_ndindex_I);
1577 contrib += w[k] * val_u_h * absDetDPhi;
1581 Kokkos::Sum<double>(avg));
1585 Comm->allreduce(avg, global_avg, 1, std::plus<T>());
1595 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
1596 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
1602 space_mirror.
nr_m = this->nr_m;
1604 return space_mirror;
1612 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
1613 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
1614 KOKKOS_FUNCTION
size_t
1617 const size_t& globalDOFIndex)
const {
1619 static_assert(
Dim == 1 ||
Dim == 2 ||
Dim == 3,
"Dim must be 1, 2 or 3");
1621 static_assert(Order == 1,
"Only order 1 is supported at the moment");
1625 this->getGlobalDOFIndices(elementNDIndex);
1630 for (
size_t i = 0; i < global_dofs.
dim; ++i) {
1631 if (global_dofs[i] == globalDOFIndex) {
1642 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
1643 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
1645 FieldLHS, FieldRHS>::DeviceStruct::numElementDOFs>
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];
1659 size_t smallestGlobalDOF = elementNDIndex.
dot(vec);
1662 globalDOFs[0] = smallestGlobalDOF;
1663 globalDOFs[1] = smallestGlobalDOF + Order;
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;
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;
1676 if constexpr (Order > 1) {
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];
1689 if constexpr (
Dim >= 3) {
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) {
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];
1716 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
1717 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
1722 FieldRHS>::point_t& localPoint)
const {
1724 static_assert(Order == 1,
"Only order 1 is supported at the moment");
1726 assert(localDOF < DeviceStruct::numElementDOFs
1727 &&
"The local vertex index is invalid");
1729 assert(this->ref_element_m.isPointInRefElement(localPoint)
1730 &&
"Point is not in reference element");
1734 const point_t ref_element_point = this->ref_element_m.getLocalVertices()[localDOF];
1739 for (
size_t d = 0; d <
Dim; d++) {
1740 if (localPoint[d] < ref_element_point[d]) {
1741 product *= localPoint[d];
1743 product *= 1.0 - localPoint[d];
1750 template <
typename T,
unsigned Dim,
unsigned Order,
typename ElementType,
1751 typename QuadratureType,
typename FieldLHS,
typename FieldRHS>
1753 FieldLHS, FieldRHS>::indices_t
1757 size_t index = vertex_index;
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;
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;
1780 return vertex_indices;
constexpr unsigned getLagrangeNumElementDOFs(unsigned Dim, unsigned Order)
constexpr KOKKOS_INLINE_FUNCTION auto first()
Implementations for FFT constructor/destructor and transforms.
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)
void parallel_reduce(const std::string &name, const ExecPolicy &policy, const FunctorType &functor, ReducerArgument &&... reducer)
std::unique_ptr< mpi::Communicator > Comm
The FiniteElementSpace class handles the mesh index mapping to vertices and elements and is the base ...
A class representing a Lagrange space for finite element methods on a structured, rectilinear grid.
T computeAvg(const FieldLHS &u_h) const
Given a field, compute the average.
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.
detail::ViewType< T, Dim >::view_type ViewType
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.
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 T evaluateRefElementShapeFunction(const size_t &localDOF, const point_t &localPoint) const
Basis functions and gradients /////////////////////////////////////.
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...
void initializeElementIndices(const Layout_t &layout)
Initialize a Kokkos view containing the element indices.
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
ElementType ref_element_m
Vector< size_t, Dim > nr_m
KOKKOS_FUNCTION indices_t getMeshVertexNDIndex(const size_t &vertex_index) const
KOKKOS_FUNCTION T evaluateRefElementShapeFunction(const size_t &localDOF, const point_t &localPoint) const
void assignGhostToPhysical(Field &field)
const NDIndex_t & getLocalNDIndex() const
KOKKOS_INLINE_FUNCTION unsigned size() const noexcept
static constexpr unsigned dim
KOKKOS_INLINE_FUNCTION T dot(const Vector< T, Dim > &rhs) const
Timing::TimerRef TimerRef
static TimerRef getTimer(const char *nm)
static void stopTimer(TimerRef t)
static void startTimer(TimerRef t)