--- /dev/null
+// ------------------------------------------------------------------------
+//
+// SPDX-License-Identifier: LGPL-2.1-or-later
+// Copyright (C) 2018 - 2023 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// Part of the source code is dual licensed under Apache-2.0 WITH
+// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
+// governing the source code and code contributions can be found in
+// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
+//
+// ------------------------------------------------------------------------
+
+
+
+// test AffineConstraints<double>::make_consistent_in_parallel for a
+// case where each process owns one DoF and on each process p>0 we
+// know about one constraint that constrains X_p=0.5*X_{p-1}. After making
+// things consistent in parallel, this needs to all reduce to X_p=factor*X1.
+
+#include <deal.II/base/mpi.h>
+
+#include <deal.II/lac/affine_constraints.h>
+
+#include "../tests.h"
+
+int
+main(int argc, char **argv)
+{
+ Utilities::MPI::MPI_InitFinalize mpi_init(argc, argv);
+ MPILogInitAll all;
+
+ const unsigned int my_proc = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD);
+ const unsigned int n_procs = Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD);
+
+ const unsigned int n_indices = n_procs;
+ IndexSet owned_elements(n_indices);
+ owned_elements.add_index(my_proc);
+ IndexSet relevant_elements(n_indices);
+ relevant_elements.add_index(my_proc);
+ if (my_proc > 0)
+ relevant_elements.add_index(my_proc - 1);
+ AffineConstraints<double> constraints(owned_elements, relevant_elements);
+
+ if (my_proc > 0)
+ constraints.add_constraint(my_proc, {{my_proc - 1, 0.5}});
+
+ deallog << "Constraints as created:" << std::endl;
+ constraints.print(deallog.get_file_stream());
+
+ constraints.make_consistent_in_parallel(owned_elements,
+ relevant_elements,
+ MPI_COMM_WORLD);
+
+ deallog << "Constraints as made consistent:" << std::endl;
+ constraints.print(deallog.get_file_stream());
+
+ return 0;
+}
--- /dev/null
+
+DEAL:0::Constraints as created:
+DEAL:0::Constraints as made consistent:
+
+DEAL:1::Constraints as created:
+ 1 0: 0.500000
+DEAL:1::Constraints as made consistent:
+ 1 0: 0.500000
+
+
+DEAL:2::Constraints as created:
+ 2 1: 0.500000
+DEAL:2::Constraints as made consistent:
+ 1 0: 0.500000
+ 2 0: 0.250000
+
+
+DEAL:3::Constraints as created:
+ 3 2: 0.500000
+DEAL:3::Constraints as made consistent:
+ 1 0: 0.500000
+ 2 0: 0.250000
+ 3 0: 0.125000
+
+
+DEAL:4::Constraints as created:
+ 4 3: 0.500000
+DEAL:4::Constraints as made consistent:
+ 1 0: 0.500000
+ 2 0: 0.250000
+ 3 0: 0.125000
+ 4 0: 0.0625000
+
--- /dev/null
+
+DEAL:0::Constraints as created:
+DEAL:0::Constraints as made consistent:
+
+DEAL:1::Constraints as created:
+ 1 0: 0.500000
+DEAL:1::Constraints as made consistent:
+ 1 0: 0.500000
+
+
+DEAL:2::Constraints as created:
+ 2 1: 0.500000
+DEAL:2::Constraints as made consistent:
+ 1 0: 0.500000
+ 2 0: 0.250000
+
+
+DEAL:3::Constraints as created:
+ 3 2: 0.500000
+DEAL:3::Constraints as made consistent:
+ 1 0: 0.500000
+ 2 0: 0.250000
+ 3 0: 0.125000
+
+
+DEAL:4::Constraints as created:
+ 4 3: 0.500000
+DEAL:4::Constraints as made consistent:
+ 1 0: 0.500000
+ 2 0: 0.250000
+ 3 0: 0.125000
+ 4 0: 0.0625000
+
+
+DEAL:5::Constraints as created:
+ 5 4: 0.500000
+DEAL:5::Constraints as made consistent:
+ 1 0: 0.500000
+ 2 0: 0.250000
+ 3 0: 0.125000
+ 4 0: 0.0625000
+ 5 0: 0.0312500
+
+
+DEAL:6::Constraints as created:
+ 6 5: 0.500000
+DEAL:6::Constraints as made consistent:
+ 1 0: 0.500000
+ 2 0: 0.250000
+ 3 0: 0.125000
+ 4 0: 0.0625000
+ 5 0: 0.0312500
+ 6 0: 0.0156250
+