]> https://gitweb.dealii.org/ - dealii.git/commitdiff
DoFInfo::read_dof_indices(): add additional index vector 12581/head
authorPeter Munch <peterrmuench@gmail.com>
Mon, 19 Jul 2021 10:29:59 +0000 (12:29 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Mon, 19 Jul 2021 12:17:19 +0000 (14:17 +0200)
include/deal.II/matrix_free/dof_info.h
include/deal.II/matrix_free/dof_info.templates.h
include/deal.II/matrix_free/matrix_free.templates.h
source/matrix_free/dof_info.cc

index cd5b0344e3ca5d0c763084d90dcb4edc48d5dfc3..90e346f0314f232f7dee9aa5b0f1913bef8adf0f 100644 (file)
@@ -184,18 +184,19 @@ namespace internal
                                     const bool with_constraints = true) const;
 
       /**
-       * This internal method takes the local indices on a cell and fills them
-       * into this class. It resolves the constraints and distributes the
-       * results. Ghost indices, i.e., indices that are located on another
-       * processor, get a temporary number by this function, and will later be
-       * assigned the correct index after all the ghost indices have been
-       * collected by the call to @p assign_ghosts.
+       * This internal method takes the local indices on a cell (two versions:
+       * hanging-node constraints resolved if possible and plain, i.e., not
+       * resolved) and fills them into this class. It resolves the constraints
+       * and distributes the results. Ghost indices, i.e., indices that are
+       * located on another processor, get a temporary number by this function,
+       * and will later be assigned the correct index after all the ghost
+       * indices have been collected by the call to @p assign_ghosts.
        */
       template <typename number>
       void
       read_dof_indices(
+        const std::vector<types::global_dof_index> &local_indices_resolved,
         const std::vector<types::global_dof_index> &local_indices,
-        const std::vector<unsigned int> &           lexicographic_inv,
         const dealii::AffineConstraints<number> &   constraints,
         const unsigned int                          cell_number,
         ConstraintValues<double> &                  constraint_values,
@@ -497,6 +498,12 @@ namespace internal
        */
       std::vector<unsigned int> dof_indices;
 
+      /**
+       * Masks indicating for each cell and component if the optimized
+       * hanging-node constraint is applicable and if yes which type.
+       */
+      std::vector<unsigned int> hanging_node_constraint_masks;
+
       /**
        * This variable describes the position of constraints in terms of the
        * local numbering of degrees of freedom on a cell. The first number
index 296e7d339b8600f677cbf6f4200bc8d5a1e76e6e..bc688a086bdfd652713c797db60844338beb5d90 100644 (file)
@@ -104,8 +104,8 @@ namespace internal
     template <typename number>
     void
     DoFInfo::read_dof_indices(
+      const std::vector<types::global_dof_index> &local_indices_resolved,
       const std::vector<types::global_dof_index> &local_indices,
-      const std::vector<unsigned int> &           lexicographic_inv,
       const dealii::AffineConstraints<number> &   constraints,
       const unsigned int                          cell_number,
       ConstraintValues<double> &                  constraint_values,
@@ -137,9 +137,8 @@ namespace internal
                i < component_dof_indices_offset[fe_index][comp + 1];
                i++)
             {
-              types::global_dof_index current_dof =
-                local_indices[lexicographic_inv[i]];
-              const auto *entries_ptr =
+              types::global_dof_index current_dof = local_indices_resolved[i];
+              const auto *            entries_ptr =
                 constraints.get_constraint_entries(current_dof);
 
               // dof is constrained
@@ -253,8 +252,7 @@ namespace internal
             {
               for (unsigned int i = 0; i < dofs_this_cell; ++i)
                 {
-                  types::global_dof_index current_dof =
-                    local_indices[lexicographic_inv[i]];
+                  types::global_dof_index current_dof = local_indices[i];
                   if (n_mpi_procs > 1 &&
                       (current_dof < first_owned || current_dof >= last_owned))
                     {
index 0325470e71c9a572e4032a117cd1f1ac04582491..66f560eb7c044cd92a364d0e823f96da4f9a84e5 100644 (file)
@@ -1050,7 +1050,8 @@ namespace internal
     AssertDimension(n_dof_handlers, locally_owned_dofs.size());
     AssertDimension(n_dof_handlers, constraint.size());
 
-    std::vector<types::global_dof_index>                local_dof_indices;
+    std::vector<types::global_dof_index> local_dof_indices_resolved;
+    std::vector<types::global_dof_index> local_dof_indices;
     std::vector<std::vector<std::vector<unsigned int>>> lexicographic(
       n_dof_handlers);
 
@@ -1196,10 +1197,18 @@ namespace internal
                     0;
                 if (dofh->get_fe_collection().size() > 1)
                   dof_info[no].cell_active_fe_index[counter] = fe_index;
-                local_dof_indices.resize(dof_info[no].dofs_per_cell[fe_index]);
-                cell_it->get_dof_indices(local_dof_indices);
+                local_dof_indices_resolved.resize(
+                  dof_info[no].dofs_per_cell[fe_index]);
+                cell_it->get_dof_indices(local_dof_indices_resolved);
+
+                local_dof_indices.resize(local_dof_indices_resolved.size());
+                for (unsigned int i = 0; i < local_dof_indices_resolved.size();
+                     ++i)
+                  local_dof_indices[i] =
+                    local_dof_indices_resolved[lexicographic[no][fe_index][i]];
+
                 dof_info[no].read_dof_indices(local_dof_indices,
-                                              lexicographic[no][fe_index],
+                                              local_dof_indices,
                                               *constraint[no],
                                               counter,
                                               constraint_values,
@@ -1223,10 +1232,18 @@ namespace internal
                   cell_level_index[counter].first,
                   cell_level_index[counter].second,
                   dofh);
-                local_dof_indices.resize(dof_info[no].dofs_per_cell[0]);
-                cell_it->get_mg_dof_indices(local_dof_indices);
+                local_dof_indices_resolved.resize(
+                  dof_info[no].dofs_per_cell[0]);
+                cell_it->get_mg_dof_indices(local_dof_indices_resolved);
+
+                local_dof_indices.resize(local_dof_indices_resolved.size());
+                for (unsigned int i = 0; i < local_dof_indices_resolved.size();
+                     ++i)
+                  local_dof_indices[i] =
+                    local_dof_indices_resolved[lexicographic[no][0][i]];
+
                 dof_info[no].read_dof_indices(local_dof_indices,
-                                              lexicographic[no][0],
+                                              local_dof_indices,
                                               *constraint[no],
                                               counter,
                                               constraint_values,
index 4b2d8b925b83bd0d55cf5f7051fda94b93ee2409..b6b25f2320719db304888c3ba9bb3d62bd1fc7b1 100644 (file)
@@ -1459,7 +1459,7 @@ namespace internal
     template void
     DoFInfo::read_dof_indices<double>(
       const std::vector<types::global_dof_index> &,
-      const std::vector<unsigned int> &,
+      const std::vector<types::global_dof_index> &,
       const dealii::AffineConstraints<double> &,
       const unsigned int,
       ConstraintValues<double> &,
@@ -1468,7 +1468,7 @@ namespace internal
     template void
     DoFInfo::read_dof_indices<float>(
       const std::vector<types::global_dof_index> &,
-      const std::vector<unsigned int> &,
+      const std::vector<types::global_dof_index> &,
       const dealii::AffineConstraints<float> &,
       const unsigned int,
       ConstraintValues<double> &,

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.