]> https://gitweb.dealii.org/ - dealii.git/commitdiff
AffineConstraints: remove the block sparsity constraint function.
authorDavid Wells <drwells@email.unc.edu>
Sun, 6 Nov 2022 13:13:25 +0000 (08:13 -0500)
committerDavid Wells <drwells@email.unc.edu>
Sat, 12 Nov 2022 20:22:02 +0000 (15:22 -0500)
This function predates BlockSparsityPatternBase<T>::add_entries(), which does
the same thing (essentially just updating the individual blocks after computing
offsets). For good measure I added a test which verifies that we get identical
output with constraints and a DoF mask table before and after the switch.

include/deal.II/lac/affine_constraints.h
include/deal.II/lac/affine_constraints.templates.h
source/lac/affine_constraints.inst.in
tests/lac/is_block_matrix.cc
tests/lac/is_block_matrix.output
tests/sparsity/block_sparsity_pattern_04.cc [new file with mode: 0644]
tests/sparsity/block_sparsity_pattern_04.output [new file with mode: 0644]

index 0fa3fafc51724c2eeee64fb7b808b6e1ba7e4946..90a723e8ee0d629f73345e25339564dc9e7fa02a 100644 (file)
@@ -1972,30 +1972,6 @@ private:
                              const bool use_inhomogeneities_for_rhs,
                              const std::integral_constant<bool, true>) const;
 
-  /**
-   * This function actually implements the local_to_global function for
-   * standard (non-block) sparsity types.
-   */
-  template <typename SparsityPatternType>
-  void
-  add_entries_local_to_global(const std::vector<size_type> &local_dof_indices,
-                              SparsityPatternType &         sparsity_pattern,
-                              const bool            keep_constrained_entries,
-                              const Table<2, bool> &dof_mask,
-                              const std::integral_constant<bool, false>) const;
-
-  /**
-   * This function actually implements the local_to_global function for block
-   * sparsity types.
-   */
-  template <typename SparsityPatternType>
-  void
-  add_entries_local_to_global(const std::vector<size_type> &local_dof_indices,
-                              SparsityPatternType &         sparsity_pattern,
-                              const bool            keep_constrained_entries,
-                              const Table<2, bool> &dof_mask,
-                              const std::integral_constant<bool, true>) const;
-
   /**
    * Internal helper function for distribute_local_to_global function.
    *
@@ -2454,49 +2430,6 @@ namespace internal
     template <typename MatrixType>
     const bool IsBlockMatrix<MatrixType>::value;
 
-
-    /**
-     * A class that can be used to determine whether a given type is a block
-     * sparsity pattern type or not. In this, it matches the IsBlockMatrix
-     * class.
-     *
-     * @see
-     * @ref GlossBlockLA "Block (linear algebra)"
-     */
-    template <typename MatrixType>
-    struct IsBlockSparsityPattern
-    {
-    private:
-      /**
-       * Overload returning true if the class is derived from
-       * BlockSparsityPatternBase, which is what block sparsity patterns do.
-       */
-      template <typename T>
-      static std::true_type
-      check(const BlockSparsityPatternBase<T> *);
-
-      /**
-       * Catch all for all other potential types that are then apparently not
-       * block sparsity patterns.
-       */
-      static std::false_type
-      check(...);
-
-    public:
-      /**
-       * A statically computable value that indicates whether the template
-       * argument to this class is a block sparsity pattern (in fact whether the
-       * type is derived from BlockSparsityPatternBase<T>).
-       */
-      static const bool value =
-        std::is_same<decltype(check(std::declval<MatrixType *>())),
-                     std::true_type>::value;
-    };
-
-    // instantiation of the static member
-    template <typename MatrixType>
-    const bool IsBlockSparsityPattern<MatrixType>::value;
-
   } // namespace AffineConstraints
 } // namespace internal
 #endif
@@ -2570,29 +2503,6 @@ AffineConstraints<number>::distribute_local_to_global(
 
 
 
-template <typename number>
-template <typename SparsityPatternType>
-inline void
-AffineConstraints<number>::add_entries_local_to_global(
-  const std::vector<size_type> &local_dof_indices,
-  SparsityPatternType &         sparsity_pattern,
-  const bool                    keep_constrained_entries,
-  const Table<2, bool> &        dof_mask) const
-{
-  // enter the internal function with the respective block information set,
-  // the actual implementation follows in the cm.templates.h file.
-  add_entries_local_to_global(
-    local_dof_indices,
-    sparsity_pattern,
-    keep_constrained_entries,
-    dof_mask,
-    std::integral_constant<bool,
-                           internal::AffineConstraints::IsBlockSparsityPattern<
-                             SparsityPatternType>::value>());
-}
-
-
-
 template <typename number>
 inline AffineConstraints<number>::ConstraintLine::ConstraintLine(
   const size_type &                                                  index,
index 1bcd7f2c41f12778fdd1a8c47e2778ee795c0550..d6ce30a995a4fd2fcf560a5b685a9ff21c84a48b 100644 (file)
@@ -4281,8 +4281,7 @@ AffineConstraints<number>::add_entries_local_to_global(
   const std::vector<size_type> &local_dof_indices,
   SparsityPatternType &         sparsity_pattern,
   const bool                    keep_constrained_entries,
-  const Table<2, bool> &        dof_mask,
-  const std::integral_constant<bool, false>) const
+  const Table<2, bool> &        dof_mask) const
 {
   Assert(sparsity_pattern.n_rows() == sparsity_pattern.n_cols(),
          ExcNotQuadratic());
@@ -4448,140 +4447,6 @@ AffineConstraints<number>::add_entries_local_to_global(
                               dof_mask);
 }
 
-
-
-template <typename number>
-template <typename SparsityPatternType>
-void
-AffineConstraints<number>::add_entries_local_to_global(
-  const std::vector<size_type> &local_dof_indices,
-  SparsityPatternType &         sparsity_pattern,
-  const bool                    keep_constrained_entries,
-  const Table<2, bool> &        dof_mask,
-  const std::integral_constant<bool, true>) const
-{
-  // just as the other add_entries_local_to_global function, but now
-  // specialized for block matrices.
-  Assert(sparsity_pattern.n_rows() == sparsity_pattern.n_cols(),
-         ExcNotQuadratic());
-  Assert(sparsity_pattern.n_block_rows() == sparsity_pattern.n_block_cols(),
-         ExcNotQuadratic());
-
-  const size_type n_local_dofs = local_dof_indices.size();
-  const size_type num_blocks   = sparsity_pattern.n_block_rows();
-
-  typename internal::AffineConstraints::ScratchDataAccessor<number>
-    scratch_data(this->scratch_data);
-
-  const bool dof_mask_is_active = (dof_mask.n_rows() == n_local_dofs);
-  if (dof_mask_is_active == true)
-    {
-      AssertDimension(dof_mask.n_cols(), n_local_dofs);
-    }
-  else
-    {
-      std::vector<size_type> &actual_dof_indices = scratch_data->columns;
-      actual_dof_indices.resize(n_local_dofs);
-      make_sorted_row_list(local_dof_indices, actual_dof_indices);
-      const size_type n_actual_dofs = actual_dof_indices.size();
-      (void)n_actual_dofs;
-
-      // additional construct that also takes care of block indices.
-      std::vector<size_type> &block_starts = scratch_data->block_starts;
-      block_starts.resize(num_blocks + 1);
-      internal::AffineConstraints::make_block_starts(sparsity_pattern,
-                                                     actual_dof_indices,
-                                                     block_starts);
-
-      for (size_type block = 0; block < num_blocks; ++block)
-        {
-          const size_type next_block = block_starts[block + 1];
-          for (size_type i = block_starts[block]; i < next_block; ++i)
-            {
-              Assert(i < n_actual_dofs, ExcInternalError());
-              const size_type row = actual_dof_indices[i];
-              Assert(row < sparsity_pattern.block(block, 0).n_rows(),
-                     ExcInternalError());
-              std::vector<size_type>::iterator index_it =
-                actual_dof_indices.begin();
-              for (size_type block_col = 0; block_col < num_blocks; ++block_col)
-                {
-                  const size_type next_block_col = block_starts[block_col + 1];
-                  sparsity_pattern.block(block, block_col)
-                    .add_entries(row,
-                                 index_it,
-                                 actual_dof_indices.begin() + next_block_col,
-                                 true);
-                  index_it = actual_dof_indices.begin() + next_block_col;
-                }
-            }
-        }
-
-      for (size_type i = 0; i < n_local_dofs; ++i)
-        if (is_constrained(local_dof_indices[i]))
-          {
-            if (keep_constrained_entries == true)
-              for (size_type j = 0; j < n_local_dofs; ++j)
-                {
-                  sparsity_pattern.add(local_dof_indices[i],
-                                       local_dof_indices[j]);
-                  sparsity_pattern.add(local_dof_indices[j],
-                                       local_dof_indices[i]);
-                }
-            else
-              sparsity_pattern.add(local_dof_indices[i], local_dof_indices[i]);
-          }
-
-      return;
-    }
-
-  // difficult case with dof_mask, similar to the distribute_local_to_global
-  // function for block matrices
-  internal::AffineConstraints::GlobalRowsFromLocal<number> &global_rows =
-    scratch_data->global_rows;
-  global_rows.reinit(n_local_dofs);
-  make_sorted_row_list(local_dof_indices, global_rows);
-  const size_type n_actual_dofs = global_rows.size();
-
-  // additional construct that also takes care of block indices.
-  std::vector<size_type> &block_starts = scratch_data->block_starts;
-  block_starts.resize(num_blocks + 1);
-  internal::AffineConstraints::make_block_starts(sparsity_pattern,
-                                                 global_rows,
-                                                 block_starts);
-
-  std::vector<size_type> &cols = scratch_data->columns;
-  cols.resize(n_actual_dofs);
-
-  // the basic difference to the non-block variant from now onwards is that we
-  // go through the blocks of the matrix separately.
-  for (size_type block = 0; block < num_blocks; ++block)
-    {
-      const size_type next_block = block_starts[block + 1];
-      for (size_type i = block_starts[block]; i < next_block; ++i)
-        {
-          const size_type row = global_rows.global_row(i);
-          for (size_type block_col = 0; block_col < num_blocks; ++block_col)
-            {
-              const size_type begin_block = block_starts[block_col],
-                              end_block   = block_starts[block_col + 1];
-              std::vector<size_type>::iterator col_ptr = cols.begin();
-              internal::AffineConstraints::resolve_matrix_row(
-                global_rows, i, begin_block, end_block, dof_mask, col_ptr);
-
-              sparsity_pattern.block(block, block_col)
-                .add_entries(row, cols.begin(), col_ptr, true);
-            }
-        }
-    }
-
-  internal::AffineConstraints::set_sparsity_diagonals(global_rows,
-                                                      local_dof_indices,
-                                                      dof_mask,
-                                                      keep_constrained_entries,
-                                                      sparsity_pattern);
-}
-
 DEAL_II_NAMESPACE_CLOSE
 
 #endif
index 2b5a77657f11dc83356938d901ee58afc11fe837..b45ba98756b560b7b16cff9dc118a202278ce24e 100644 (file)
@@ -263,8 +263,7 @@ for (S : REAL_AND_COMPLEX_SCALARS; SP : AFFINE_CONSTRAINTS_SP)
       const std::vector<AffineConstraints<S>::size_type> &,
       SP &,
       const bool,
-      const Table<2, bool> &,
-      std::integral_constant<bool, false>) const;
+      const Table<2, bool> &) const;
 
     template void AffineConstraints<S>::add_entries_local_to_global<SP>(
       const std::vector<AffineConstraints<S>::size_type> &,
@@ -288,8 +287,7 @@ for (S : REAL_AND_COMPLEX_SCALARS; SP : AFFINE_CONSTRAINTS_SP_BLOCK)
       const std::vector<AffineConstraints<S>::size_type> &,
       SP &,
       const bool,
-      const Table<2, bool> &,
-      std::integral_constant<bool, true>) const;
+      const Table<2, bool> &) const;
 
     template void AffineConstraints<S>::add_entries_local_to_global<SP>(
       const std::vector<AffineConstraints<S>::size_type> &,
index e2156365f40e7ad4d8f901b7734b644b01d99ccc..c2137c670e17d2ce9dba64b66ca400b4098bd679 100644 (file)
@@ -58,20 +58,6 @@ test()
           << internal::AffineConstraints::IsBlockMatrix<
                BlockSparseMatrixEZ<float>>::value
           << std::endl;
-
-  deallog << internal::AffineConstraints::IsBlockSparsityPattern<
-               SparsityPattern>::value
-          << ' '
-          << internal::AffineConstraints::IsBlockSparsityPattern<
-               DynamicSparsityPattern>::value
-          << std::endl;
-
-  deallog << internal::AffineConstraints::IsBlockSparsityPattern<
-               BlockSparsityPattern>::value
-          << ' '
-          << internal::AffineConstraints::IsBlockSparsityPattern<
-               BlockDynamicSparsityPattern>::value
-          << std::endl;
 }
 
 
index 3f7223107b3350efb3199ab9d45d2115b07440d4..c7c5f754bf7cf46bc1a424c06177b97d83ef44b5 100644 (file)
@@ -1,5 +1,3 @@
 
 DEAL::0 0 0 0
 DEAL::1 1 1 1
-DEAL::0 0
-DEAL::1 1
diff --git a/tests/sparsity/block_sparsity_pattern_04.cc b/tests/sparsity/block_sparsity_pattern_04.cc
new file mode 100644 (file)
index 0000000..d0a1522
--- /dev/null
@@ -0,0 +1,85 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2022 - 2022 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+// Verify that AffineConstraints::add_entries_local_to_global() works correctly.
+// In particular, there used to be (written before
+// BlockSparsityPattern::add_entries() existed) a block-specific version of that
+// function. This test verifies that the non-block version, combined with the
+// block version of add_entries(), produces the same output.
+
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/dofs/dof_tools.h>
+
+#include <deal.II/fe/fe_dgq.h>
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_system.h>
+#include <deal.II/fe/mapping_cartesian.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+
+#include <deal.II/lac/affine_constraints.h>
+#include <deal.II/lac/block_sparsity_pattern.h>
+
+#include <deal.II/numerics/vector_tools_boundary.h>
+
+#include <vector>
+
+#include "../tests.h"
+
+int
+main()
+{
+  initlog();
+
+  Triangulation<2> tria;
+  GridGenerator::hyper_cube(tria);
+  tria.refine_global(1);
+  tria.begin_active()->set_refine_flag();
+  tria.execute_coarsening_and_refinement();
+
+  FESystem<2>        fe(FE_Q<2>(2), 2, FE_DGQ<2>(1), 1, FE_DGQ<2>(0), 2);
+  const unsigned int n_components = fe.n_components();
+  const unsigned int n_blocks     = fe.n_blocks();
+  DoFHandler<2>      dof_handler(tria);
+  dof_handler.distribute_dofs(fe);
+
+  AffineConstraints<double> constraints;
+  VectorTools::interpolate_boundary_values(
+    MappingCartesian<2>(),
+    dof_handler,
+    0,
+    Functions::ConstantFunction<2>(42.0, n_components),
+    constraints,
+    ComponentMask({true, true, false, false, false}));
+  DoFTools::make_hanging_node_constraints(dof_handler, constraints);
+
+  Table<2, DoFTools::Coupling> coupling(n_components, n_components);
+  for (unsigned int c = 0; c < n_components; ++c)
+    for (unsigned int d = 0; d < n_components; ++d)
+      // This is completely arbitrary for this test so just flip back and forth
+      coupling[c][d] = (c + d) % 2 == 0 ? DoFTools::always : DoFTools::none;
+
+  const std::vector<types::global_dof_index> dofs_per_block =
+    DoFTools::count_dofs_per_fe_block(dof_handler);
+
+  BlockDynamicSparsityPattern dsp(dofs_per_block, dofs_per_block);
+  DoFTools::make_sparsity_pattern(
+    dof_handler, coupling, dsp, constraints, false);
+
+  BlockSparsityPattern sp(n_blocks, n_blocks);
+  sp.copy_from(dsp);
+  sp.print(deallog.get_file_stream());
+}
diff --git a/tests/sparsity/block_sparsity_pattern_04.output b/tests/sparsity/block_sparsity_pattern_04.output
new file mode 100644 (file)
index 0000000..3ec84dc
--- /dev/null
@@ -0,0 +1,129 @@
+
+[0,0,4,8,24,34,43,66,70,90,92,94,95,96,97,99,104,120,122,123,124,125,127]
+[1,1,5,9,25,35,44,67,71,91,93,98,105,121,126]
+[2,2,45]
+[3,3,46]
+[4,0,4,8,14,16,18,19,20,21,23,24,32,34,38,40,41,42,43,45,47,52,54,55,56,57,59,66,70,74,90,92,94,95,96,97,99,104,108,110,111,112,113,115,120,122,123,124,125,127]
+[5,1,5,9,15,17,22,25,33,35,39,44,48,53,58,67,71,75,91,93,98,105,109,114,121,126]
+[6,6,49]
+[7,7,50]
+[8,0,4,8,14,16,18,19,20,21,23,24,34,51,66,70,90,92,94,95,96,97,99,104,120,122,123,124,125,127]
+[9,1,5,9,15,17,22,25,35,52,67,71,91,93,98,105,121,126]
+[10,10,53]
+[11,11,54]
+[12,12,55]
+[13,13,56]
+[14,4,8,14,16,18,19,20,21,23,32,52,54,55,56,57,59]
+[15,5,9,15,17,22,33,53,58]
+[16,4,8,14,16,18,19,20,21,23,59]
+[17,5,9,15,17,22,60]
+[18,4,8,14,16,18,19,20,21,23,61]
+[19,4,8,14,16,18,19,20,21,23,62]
+[20,4,8,14,16,18,19,20,21,23,63]
+[21,4,8,14,16,18,19,20,21,23,64]
+[22,5,9,15,17,22,65]
+[23,4,8,14,16,18,19,20,21,23,66]
+[24,0,4,8,24,34,66,67,74,90,104,108,110,111,112,113,115,120,122,123,124,125,127]
+[25,1,5,9,25,35,67,68,75,91,105,109,114,121,126]
+[26,26,69]
+[27,27,70]
+[28,28,71]
+[29,29,72]
+[30,30,73]
+[31,31,74]
+[32,4,14,32,34,38,40,41,42,43,45,52,54,55,56,57,59,75]
+[33,5,15,33,35,39,44,53,58,76]
+[34,0,4,8,24,32,34,38,40,41,42,43,45,66,74,77,90,104,108,110,111,112,113,115,120,122,123,124,125,127]
+[35,1,5,9,25,33,35,39,44,67,75,78,91,105,109,114,121,126]
+[36,36,79]
+[37,37,80]
+[38,4,32,34,38,40,41,42,43,45,81]
+[39,5,33,35,39,44,82]
+[40,4,32,34,38,40,41,42,43,45,83]
+[41,4,32,34,38,40,41,42,43,45,84]
+[42,4,32,34,38,40,41,42,43,45,85]
+[43,0,4,32,34,38,40,41,42,43,45]
+[44,1,5,33,35,39,44]
+[45,2,4,32,34,38,40,41,42,43,45]
+[46,3,46]
+[47,4,47]
+[48,5,48]
+[49,6,49]
+[50,7,50]
+[51,8,51]
+[52,4,9,14,32,52,54,55,56,57,59]
+[53,5,10,15,33,53,58]
+[54,4,11,14,32,52,54,55,56,57,59]
+[55,4,12,14,32,52,54,55,56,57,59]
+[56,4,13,14,32,52,54,55,56,57,59]
+[57,4,14,32,52,54,55,56,57,59]
+[58,5,15,33,53,58]
+[59,4,14,16,32,52,54,55,56,57,59]
+[60,17,60]
+[61,18,61]
+[62,19,62]
+[63,20,63]
+[64,21,64]
+[65,22,65]
+[66,0,4,8,23,24,34,66,70,74,76,78,79,80,81,83,90,92,94,95,96,97,99,104,108,110,111,112,113,115,120,122,123,124,125,127]
+[67,1,5,9,24,25,35,67,71,75,77,82,91,93,98,105,109,114,121,126]
+[68,25,68]
+[69,26,69]
+[70,0,4,8,27,66,70,74,76,78,79,80,81,83,90,92,94,95,96,97,99]
+[71,1,5,9,28,67,71,75,77,82,91,93,98]
+[72,29,72]
+[73,30,73]
+[74,4,24,31,34,66,70,74,76,78,79,80,81,83,104,108,110,111,112,113,115]
+[75,5,25,32,35,67,71,75,77,82,105,109,114]
+[76,33,66,70,74,76,78,79,80,81,83]
+[77,34,67,71,75,77,82]
+[78,35,66,70,74,76,78,79,80,81,83]
+[79,36,66,70,74,76,78,79,80,81,83]
+[80,37,66,70,74,76,78,79,80,81,83]
+[81,38,66,70,74,76,78,79,80,81,83]
+[82,39,67,71,75,77,82]
+[83,40,66,70,74,76,78,79,80,81,83]
+[84,41,84]
+[85,42,85]
+[86,86]
+[87,87]
+[88,88]
+[89,89]
+[90,0,4,8,24,34,66,70,90,92,94,95,96,97,99,104,120,122,123,124,125,127]
+[91,1,5,9,25,35,67,71,91,93,98,105,121,126]
+[92,0,4,8,66,70,90,92,94,95,96,97,99]
+[93,1,5,9,67,71,91,93,98]
+[94,0,4,8,66,70,90,92,94,95,96,97,99]
+[95,0,4,8,66,70,90,92,94,95,96,97,99]
+[96,0,4,8,66,70,90,92,94,95,96,97,99]
+[97,0,4,8,66,70,90,92,94,95,96,97,99]
+[98,1,5,9,67,71,91,93,98]
+[99,0,4,8,66,70,90,92,94,95,96,97,99]
+[100,100]
+[101,101]
+[102,102]
+[103,103]
+[104,0,4,8,24,34,66,74,90,104,108,110,111,112,113,115,120,122,123,124,125,127]
+[105,1,5,9,25,35,67,75,91,105,109,114,121,126]
+[106,106]
+[107,107]
+[108,4,24,34,66,74,104,108,110,111,112,113,115]
+[109,5,25,35,67,75,105,109,114]
+[110,4,24,34,66,74,104,108,110,111,112,113,115]
+[111,4,24,34,66,74,104,108,110,111,112,113,115]
+[112,4,24,34,66,74,104,108,110,111,112,113,115]
+[113,4,24,34,66,74,104,108,110,111,112,113,115]
+[114,5,25,35,67,75,105,109,114,121]
+[115,4,24,34,66,74,104,108,110,111,112,113,115,122]
+[116,116,123]
+[117,117,124]
+[118,118,125]
+[119,119,126]
+[120,0,4,8,24,34,66,90,104,120,122,123,124,125,127]
+[121,1,5,9,25,35,67,91,105,114,121,126]
+[122,0,4,8,24,34,66,90,104,115,120,122,123,124,125,127]
+[123,0,4,8,24,34,66,90,104,116,120,122,123,124,125,127]
+[124,0,4,8,24,34,66,90,104,117,120,122,123,124,125,127]
+[125,0,4,8,24,34,66,90,104,118,120,122,123,124,125,127]
+[126,1,5,9,25,35,67,91,105,119,121,126]
+[127,0,4,8,24,34,66,90,104,120,122,123,124,125,127]

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.