--- /dev/null
+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)
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
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();
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);
}
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
+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
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)
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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();
+}
--- /dev/null
+
+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]