]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Portable MatrixFree: Some cleanup for fe_degree+1 != n_q_points_1d 17863/head
authorMartin Kronbichler <martin.kronbichler@rub.de>
Wed, 13 Nov 2024 14:07:06 +0000 (15:07 +0100)
committerMartin Kronbichler <martin.kronbichler@rub.de>
Wed, 13 Nov 2024 14:07:06 +0000 (15:07 +0100)
include/deal.II/matrix_free/portable_matrix_free.templates.h
include/deal.II/matrix_free/tools.h

index ed94802a6e691e299c939d683714dcc0a2fd48bf..422010cfe69609f719e350888b23e1e49017a118 100644 (file)
@@ -182,7 +182,7 @@ namespace Portable
             Kokkos::view_alloc("JxW_" + std::to_string(color),
                                Kokkos::WithoutInitializing),
             n_cells,
-            scalar_dofs_per_cell);
+            q_points_per_cell);
 
       if (update_flags & update_gradients)
         data->inv_jacobian[color] =
@@ -190,7 +190,7 @@ namespace Portable
             Kokkos::view_alloc("inv_jacobian_" + std::to_string(color),
                                Kokkos::WithoutInitializing),
             n_cells,
-            scalar_dofs_per_cell);
+            q_points_per_cell);
 
       // Initialize to zero, i.e., unconstrained cell
       data->constraint_mask[color] =
index a8903d8d5e4c53cbfdb6628434637d062a1f920d..c68490128d29981b2d6d58eaf788a1ce7ce332ed 100644 (file)
@@ -1374,7 +1374,11 @@ namespace MatrixFreeTools
   } // namespace internal
 
 
-  template <int dim, int fe_degree, typename Number, typename QuadOperation>
+  template <int dim,
+            int fe_degree,
+            int n_q_points_1d,
+            typename Number,
+            typename QuadOperation>
   class CellAction
   {
   public:
@@ -1393,7 +1397,7 @@ namespace MatrixFreeTools
                const Number *,
                Number *dst) const
     {
-      Portable::FEEvaluation<dim, fe_degree, fe_degree + 1, 1, Number> fe_eval(
+      Portable::FEEvaluation<dim, fe_degree, n_q_points_1d, 1, Number> fe_eval(
         gpu_data, shared_data);
       m_quad_operation.set_matrix_free_data(*gpu_data);
       m_quad_operation.set_cell(cell);
@@ -1490,8 +1494,8 @@ namespace MatrixFreeTools
     matrix_free.initialize_dof_vector(diagonal_global);
 
 
-    CellAction<dim, fe_degree, Number, QuadOperation> cell_action(
-      quad_operation, evaluation_flags, integration_flags);
+    CellAction<dim, fe_degree, n_q_points_1d, Number, QuadOperation>
+      cell_action(quad_operation, evaluation_flags, integration_flags);
     LinearAlgebra::distributed::Vector<Number, MemorySpace> dummy;
     matrix_free.cell_loop(cell_action, dummy, diagonal_global);
 

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.