]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Another rectangular version of add_entries_local_to_global() 13415/head
authorLuca Heltai <luca.heltai@sissa.it>
Fri, 18 Feb 2022 08:47:30 +0000 (08:47 +0000)
committerLuca Heltai <luca.heltai@sissa.it>
Fri, 18 Feb 2022 08:50:16 +0000 (08:50 +0000)
doc/news/changes/minor/20220218LucaHeltai [new file with mode: 0644]
include/deal.II/lac/affine_constraints.h
include/deal.II/lac/affine_constraints.templates.h
source/lac/affine_constraints.inst.in
tests/lac/constraints_add_entries_local_to_global_rect.cc [new file with mode: 0644]
tests/lac/constraints_add_entries_local_to_global_rect.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20220218LucaHeltai b/doc/news/changes/minor/20220218LucaHeltai
new file mode 100644 (file)
index 0000000..68975e0
--- /dev/null
@@ -0,0 +1,4 @@
+New: AffineConstraints::add_entries_local_to_global() now supports also the case
+where columns don't use the same constraints as rows. 
+<br> 
+(Luca Heltai, 2022/02/18)
index eaa2d95da87218b81756099fd36a0acc49b34f2b..2a1ecec635ddc219de65ba7efd811b40926369c5 100644 (file)
@@ -1525,6 +1525,20 @@ public:
     const bool                    keep_constrained_entries = true,
     const Table<2, bool> &        dof_mask = Table<2, bool>()) const;
 
+  /**
+   * Similar to the other function, but for non-quadratic sparsity patterns, and
+   * for different constraints in the column space.
+   */
+  template <typename SparsityPatternType>
+  void
+  add_entries_local_to_global(
+    const std::vector<size_type> &   row_indices,
+    const AffineConstraints<number> &col_constraints,
+    const std::vector<size_type> &   col_indices,
+    SparsityPatternType &            sparsity_pattern,
+    const bool                       keep_constrained_entries = true,
+    const Table<2, bool> &           dof_mask = Table<2, bool>()) const;
+
   /**
    * This function imports values from a global vector (@p global_vector) by
    * applying the constraints to a vector of local values, expressed in
index 76c6d8e68af199fcca3a91f5c59f517d1fc26628..489624c50c340d7e5c7078a62723f8b66ea1e724 100644 (file)
@@ -4356,11 +4356,12 @@ template <typename number>
 template <typename SparsityPatternType>
 void
 AffineConstraints<number>::add_entries_local_to_global(
-  const std::vector<size_type> &row_indices,
-  const std::vector<size_type> &col_indices,
-  SparsityPatternType &         sparsity_pattern,
-  const bool                    keep_constrained_entries,
-  const Table<2, bool> &        dof_mask) const
+  const std::vector<size_type> &   row_indices,
+  const AffineConstraints<number> &col_constraints,
+  const std::vector<size_type> &   col_indices,
+  SparsityPatternType &            sparsity_pattern,
+  const bool                       keep_constrained_entries,
+  const Table<2, bool> &           dof_mask) const
 {
   const size_type n_local_rows = row_indices.size();
   const size_type n_local_cols = col_indices.size();
@@ -4374,7 +4375,7 @@ AffineConstraints<number>::add_entries_local_to_global(
           for (const size_type col_index : col_indices)
             sparsity_pattern.add(row_index, col_index);
       for (const size_type col_index : col_indices)
-        if (is_constrained(col_index))
+        if (col_constraints.is_constrained(col_index))
           for (const size_type row_index : row_indices)
             sparsity_pattern.add(row_index, col_index);
     }
@@ -4390,7 +4391,7 @@ AffineConstraints<number>::add_entries_local_to_global(
       std::vector<size_type> actual_row_indices(n_local_rows);
       std::vector<size_type> actual_col_indices(n_local_cols);
       make_sorted_row_list(row_indices, actual_row_indices);
-      make_sorted_row_list(col_indices, actual_col_indices);
+      col_constraints.make_sorted_row_list(col_indices, actual_col_indices);
       const size_type n_actual_rows = actual_row_indices.size();
 
       // now add the indices we collected above to the sparsity pattern. Very
@@ -4409,6 +4410,27 @@ AffineConstraints<number>::add_entries_local_to_global(
 
 
 
+template <typename number>
+template <typename SparsityPatternType>
+void
+AffineConstraints<number>::add_entries_local_to_global(
+  const std::vector<size_type> &row_indices,
+  const std::vector<size_type> &col_indices,
+  SparsityPatternType &         sparsity_pattern,
+  const bool                    keep_constrained_entries,
+  const Table<2, bool> &        dof_mask) const
+{
+  // Call the function with the same name that takes a column constraint as well
+  add_entries_local_to_global(row_indices,
+                              *this,
+                              col_indices,
+                              sparsity_pattern,
+                              keep_constrained_entries,
+                              dof_mask);
+}
+
+
+
 template <typename number>
 template <typename SparsityPatternType>
 void
index 4fcb197a8a549d9e452258f7f3be1fe525d96007..100236c2705beff936bfa24748799011d4424e67 100644 (file)
@@ -272,6 +272,14 @@ for (S : REAL_AND_COMPLEX_SCALARS; SP : AFFINE_CONSTRAINTS_SP)
       SP &,
       const bool,
       const Table<2, bool> &) const;
+
+    template void AffineConstraints<S>::add_entries_local_to_global<SP>(
+      const std::vector<AffineConstraints<S>::size_type> &,
+      const AffineConstraints<S> &,
+      const std::vector<AffineConstraints<S>::size_type> &,
+      SP &,
+      const bool,
+      const Table<2, bool> &) const;
   }
 
 for (S : REAL_AND_COMPLEX_SCALARS; SP : AFFINE_CONSTRAINTS_SP_BLOCK)
diff --git a/tests/lac/constraints_add_entries_local_to_global_rect.cc b/tests/lac/constraints_add_entries_local_to_global_rect.cc
new file mode 100644 (file)
index 0000000..78bb9d2
--- /dev/null
@@ -0,0 +1,65 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2017 - 2018 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.
+//
+// ---------------------------------------------------------------------
+
+
+
+// this function tests the correctness of the implementation of
+// AffineConstraints<double>::add_entries_local_to_global for row and column
+// indices with different constraints on the rows and columns.
+
+#include <deal.II/base/function.h>
+
+#include <deal.II/lac/affine_constraints.h>
+#include <deal.II/lac/dynamic_sparsity_pattern.h>
+
+#include <iostream>
+
+#include "../tests.h"
+
+void
+test()
+{
+  DynamicSparsityPattern    dsp1(4, 6), dsp2(4, 6);
+  AffineConstraints<double> cm1, cm2;
+  cm1.add_line(2);
+  cm1.add_entry(2, 1, 0.5);
+  cm1.add_entry(2, 3, 0.5);
+  cm1.close();
+  cm2.add_line(4);
+  cm2.add_entry(4, 2, 0.3);
+  cm2.add_entry(4, 5, 0.7);
+  cm2.close();
+
+  std::vector<types::global_dof_index> indices1(2), indices2(2);
+  indices1[0] = 1;
+  indices1[1] = 2;
+  indices2[0] = 4;
+  indices2[1] = 5;
+  cm1.add_entries_local_to_global(indices1, indices2, dsp1);
+  deallog << "Same constraints: " << std::endl;
+  dsp1.print(deallog.get_file_stream());
+  cm1.add_entries_local_to_global(indices1, cm2, indices2, dsp2);
+  deallog << "Different constraints: " << std::endl;
+  dsp2.print(deallog.get_file_stream());
+}
+
+
+int
+main()
+{
+  initlog();
+
+  test();
+}
diff --git a/tests/lac/constraints_add_entries_local_to_global_rect.output b/tests/lac/constraints_add_entries_local_to_global_rect.output
new file mode 100644 (file)
index 0000000..c9cf61d
--- /dev/null
@@ -0,0 +1,11 @@
+
+DEAL::Same constraints: 
+[0]
+[1,4,5]
+[2,4,5]
+[3,4,5]
+DEAL::Different constraints: 
+[0]
+[1,2,4,5]
+[2,4,5]
+[3,2,5]

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.