]> https://gitweb.dealii.org/ - dealii.git/commitdiff
internal::MatrixFreeFunctions::ConstraintInfo: switch data type 14279/head
authorPeter Munch <peterrmuench@gmail.com>
Sun, 18 Sep 2022 12:10:57 +0000 (14:10 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Sun, 18 Sep 2022 12:10:57 +0000 (14:10 +0200)
include/deal.II/matrix_free/constraint_info.h
include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h
tests/matrix_free/constraint_info_01.cc

index 44a13c9f65055a81ab4cb4b7dd82de88f40513de..91c735c36f7cbae6391930bc8e1b15ead8fb5f58 100644 (file)
@@ -77,13 +77,13 @@ namespace internal
 
       template <typename T, typename VectorType>
       void
-      read_write_operation(const T &              operation,
-                           VectorType &           global_vector,
-                           AlignedVector<Number> &local_vector,
-                           const unsigned int     first_cell,
-                           const unsigned int     n_cells,
-                           const unsigned int     n_dofs_per_cell,
-                           const bool             apply_constraints) const;
+      read_write_operation(const T &          operation,
+                           VectorType &       global_vector,
+                           Number *           local_vector,
+                           const unsigned int first_cell,
+                           const unsigned int n_cells,
+                           const unsigned int n_dofs_per_cell,
+                           const bool         apply_constraints) const;
 
       void
       apply_hanging_node_constraints(
@@ -457,13 +457,13 @@ namespace internal
     template <typename T, typename VectorType>
     inline void
     ConstraintInfo<dim, Number>::read_write_operation(
-      const T &              operation,
-      VectorType &           global_vector,
-      AlignedVector<Number> &local_vector,
-      const unsigned int     first_cell,
-      const unsigned int     n_cells,
-      const unsigned int     n_dofs_per_cell,
-      const bool             apply_constraints) const
+      const T &          operation,
+      VectorType &       global_vector,
+      Number *           local_vector,
+      const unsigned int first_cell,
+      const unsigned int n_cells,
+      const unsigned int n_dofs_per_cell,
+      const bool         apply_constraints) const
     {
       if (apply_constraints == false)
         {
index ce715bcf46629a98c2552c4df3e7a5f93a739da5..08de0842e91ce660f4914efcc14ed0cdb3a152f1 100644 (file)
@@ -2557,7 +2557,7 @@ MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
           internal::VectorReader<Number, VectorizedArrayType> reader;
           constraint_info.read_write_operation(reader,
                                                *vec_coarse_ptr,
-                                               evaluation_data_coarse,
+                                               evaluation_data_coarse.data(),
                                                cell_counter,
                                                n_lanes_filled,
                                                scheme.n_dofs_per_cell_coarse,
@@ -2761,7 +2761,7 @@ MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
             cell_counter, n_lanes_filled, true, evaluation_data_coarse);
           constraint_info.read_write_operation(writer,
                                                *vec_coarse_ptr,
-                                               evaluation_data_coarse,
+                                               evaluation_data_coarse.data(),
                                                cell_counter,
                                                n_lanes_filled,
                                                scheme.n_dofs_per_cell_coarse,
@@ -2894,7 +2894,7 @@ MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
           internal::VectorSetter<Number, VectorizedArrayType> writer;
           constraint_info.read_write_operation(writer,
                                                *vec_coarse_ptr,
-                                               evaluation_data_coarse,
+                                               evaluation_data_coarse.data(),
                                                cell_counter,
                                                n_lanes_filled,
                                                scheme.n_dofs_per_cell_coarse,
index 2cc5a69ec43c2e63f3d8f769ef70bee1c307ea5b..5047b6cc7de5baf11135b21e1dd73f5e6164af29 100644 (file)
@@ -133,8 +133,13 @@ main()
                                    const unsigned int n_cells,
                                    const unsigned int n_dofs_per_cell) {
     internal::VectorReader<Number, VectorizedArrayType> reader;
-    constraint_info.read_write_operation(
-      reader, src, local_vector, first_cell, n_cells, n_dofs_per_cell, true);
+    constraint_info.read_write_operation(reader,
+                                         src,
+                                         local_vector.data(),
+                                         first_cell,
+                                         n_cells,
+                                         n_dofs_per_cell,
+                                         true);
     constraint_info.apply_hanging_node_constraints(first_cell,
                                                    n_cells,
                                                    false,
@@ -151,24 +156,39 @@ main()
                                                      n_cells,
                                                      true,
                                                      local_vector);
-      constraint_info.read_write_operation(
-        writer, dst, local_vector, first_cell, n_cells, n_dofs_per_cell, true);
+      constraint_info.read_write_operation(writer,
+                                           dst,
+                                           local_vector.data(),
+                                           first_cell,
+                                           n_cells,
+                                           n_dofs_per_cell,
+                                           true);
     };
 
   const auto read_dof_values_plain = [&](const unsigned int first_cell,
                                          const unsigned int n_cells,
                                          const unsigned int n_dofs_per_cell) {
     internal::VectorReader<Number, VectorizedArrayType> reader;
-    constraint_info.read_write_operation(
-      reader, src, local_vector, first_cell, n_cells, n_dofs_per_cell, false);
+    constraint_info.read_write_operation(reader,
+                                         src,
+                                         local_vector.data(),
+                                         first_cell,
+                                         n_cells,
+                                         n_dofs_per_cell,
+                                         false);
   };
 
   const auto set_dof_values_plain = [&](const unsigned int first_cell,
                                         const unsigned int n_cells,
                                         const unsigned int n_dofs_per_cell) {
     internal::VectorSetter<Number, VectorizedArrayType> writer;
-    constraint_info.read_write_operation(
-      writer, dst, local_vector, first_cell, n_cells, n_dofs_per_cell, false);
+    constraint_info.read_write_operation(writer,
+                                         dst,
+                                         local_vector.data(),
+                                         first_cell,
+                                         n_cells,
+                                         n_dofs_per_cell,
+                                         false);
   };
 
   unsigned int       first_cell      = 0;

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.