]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Compute DoFInfo::hanging_node_constraint_masks_comp once 13528/head
authorPeter Munch <peterrmuench@gmail.com>
Thu, 10 Mar 2022 08:28:30 +0000 (09:28 +0100)
committerPeter Munch <peterrmuench@gmail.com>
Thu, 10 Mar 2022 08:28:30 +0000 (09:28 +0100)
include/deal.II/matrix_free/dof_info.templates.h
include/deal.II/matrix_free/fe_evaluation.h
include/deal.II/matrix_free/matrix_free.templates.h
include/deal.II/matrix_free/tools.h
source/matrix_free/dof_info.cc

index 9e3571a1544cf400e1f5fe11327dea035258c72f..2b4d859e167a93c265d37ea9f43ba253421748d9 100644 (file)
@@ -283,18 +283,7 @@ namespace internal
       const TriaIterator<DoFCellAccessor<dim, dim, false>> &cell,
       std::vector<types::global_dof_index> &                dof_indices)
     {
-      // 1) check if finite elements support fast hanging-node algorithm
-      this->hanging_node_constraint_masks_comp =
-        hanging_nodes.compute_supported_components(
-          cell->get_dof_handler().get_fe_collection());
-
-      if ([](const auto &supported_components) {
-            return std::none_of(supported_components.begin(),
-                                supported_components.end(),
-                                [](const auto &a) {
-                                  return *std::max_element(a.begin(), a.end());
-                                });
-          }(hanging_node_constraint_masks_comp))
+      if (this->hanging_node_constraint_masks_comp.size() == 0)
         return false;
 
       // 2) determine the refinement configuration of the cell
index e5f58db72788b1e8b2901a08677bcc9d39938be4..1464223eb4e2b08ce7fcaedd36eb8a1ecbff5584 100644 (file)
@@ -3325,6 +3325,7 @@ FEEvaluationBase<dim, n_components_, Number, is_face, VectorizedArrayType>::
       for (unsigned int v = 0; v < n_lanes; ++v)
         if (cells[v] != numbers::invalid_unsigned_int &&
             this->dof_info->hanging_node_constraint_masks.size() > 0 &&
+            this->dof_info->hanging_node_constraint_masks_comp.size() > 0 &&
             this->dof_info->hanging_node_constraint_masks[cells[v]] !=
               internal::MatrixFreeFunctions::ConstraintKinds::unconstrained &&
             this->dof_info->hanging_node_constraint_masks_comp
@@ -3427,6 +3428,7 @@ FEEvaluationBase<dim, n_components_, Number, is_face, VectorizedArrayType>::
             has_constraints = true;
 
           if (this->dof_info->hanging_node_constraint_masks.size() > 0 &&
+              this->dof_info->hanging_node_constraint_masks_comp.size() > 0 &&
               this->dof_info->hanging_node_constraint_masks[cells[v]] !=
                 internal::MatrixFreeFunctions::ConstraintKinds::unconstrained &&
               this->dof_info->hanging_node_constraint_masks_comp
@@ -3515,6 +3517,7 @@ FEEvaluationBase<dim, n_components_, Number, is_face, VectorizedArrayType>::
              this->dof_info->row_starts[cell_dof_index + n_components_read]
                .second ||
            ((this->dof_info->hanging_node_constraint_masks.size() > 0 &&
+             this->dof_info->hanging_node_constraint_masks_comp.size() > 0 &&
              this->dof_info->hanging_node_constraint_masks[cell_index] !=
                internal::MatrixFreeFunctions::ConstraintKinds::unconstrained) &&
             this->dof_info->hanging_node_constraint_masks_comp
@@ -4155,6 +4158,7 @@ FEEvaluationBase<dim, n_components_, Number, is_face, VectorizedArrayType>::
 {
   if (this->dof_info == nullptr ||
       this->dof_info->hanging_node_constraint_masks.size() == 0 ||
+      this->dof_info->hanging_node_constraint_masks_comp.size() == 0 ||
       this->dof_info->hanging_node_constraint_masks_comp
           [this->active_fe_index][this->first_selected_component] == false)
     return; // nothing to do with faces
index 284da0c74ad6b507d2ae314e06a7b04f856d5951..a2f19515281b6ea4b8f6c509a5db92095c144ccb 100644 (file)
@@ -1241,7 +1241,23 @@ namespace internal
         hanging_nodes = std::make_unique<
           dealii::internal::MatrixFreeFunctions::HangingNodes<dim>>(tria);
         for (unsigned int no = 0; no < n_dof_handlers; ++no)
-          dof_info[no].hanging_node_constraint_masks.resize(n_active_cells);
+          {
+            dof_info[no].hanging_node_constraint_masks.resize(n_active_cells);
+
+            dof_info[no].hanging_node_constraint_masks_comp =
+              hanging_nodes->compute_supported_components(
+                dof_handler[no]->get_fe_collection());
+
+            if ([](const auto &supported_components) {
+                  return std::none_of(supported_components.begin(),
+                                      supported_components.end(),
+                                      [](const auto &a) {
+                                        return *std::max_element(a.begin(),
+                                                                 a.end());
+                                      });
+                }(dof_info[no].hanging_node_constraint_masks_comp))
+              dof_info[no].hanging_node_constraint_masks_comp = {};
+          }
       }
 
     for (unsigned int counter = 0; counter < n_active_cells; ++counter)
index 39d52486a36f173cc54a2bc6e9423374cb40e865..fc62fe3a3e9034d5bcaeeae90e4dbbcc5f9876a7 100644 (file)
@@ -448,6 +448,7 @@ namespace MatrixFreeTools
 
             // STEP 2c: apply hanging-node constraints
             if (dof_info.hanging_node_constraint_masks.size() > 0 &&
+                dof_info.hanging_node_constraint_masks_comp.size() > 0 &&
                 dof_info
                   .hanging_node_constraint_masks_comp[phi.get_active_fe_index()]
                                                      [first_selected_component])
index c9e2717eab504fe0b33d30ed85bd9560a9a1fbc5..38cfa797ef8d6f8acf73b6d7a7e7025f2d48973e 100644 (file)
@@ -104,6 +104,7 @@ namespace internal
           // one
           const bool has_constraints =
             (hanging_node_constraint_masks.size() != 0 &&
+             hanging_node_constraint_masks_comp.size() != 0 &&
              hanging_node_constraint_masks[cell * n_vectorization + v] !=
                ConstraintKinds::unconstrained &&
              hanging_node_constraint_masks_comp[fe_index][0 /*TODO*/]) ||
@@ -235,6 +236,7 @@ namespace internal
                   AssertIndexRange(fe_index, dofs_per_cell.size());
 
                   if (hanging_node_constraint_masks.size() > 0 &&
+                      hanging_node_constraint_masks_comp.size() > 0 &&
                       hanging_node_constraint_masks[boundary_cells[i]] !=
                         ConstraintKinds::unconstrained)
                     for (unsigned int comp = 0; comp < n_components; ++comp)
@@ -374,7 +376,8 @@ namespace internal
 
               bool has_hanging_nodes = false;
 
-              if (hanging_node_constraint_masks.size() > 0)
+              if (hanging_node_constraint_masks.size() > 0 &&
+                  hanging_node_constraint_masks_comp.size() > 0)
                 {
                   const auto mask =
                     hanging_node_constraint_masks[renumbering[position_cell +

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.