]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add test.
authorWolfgang Bangerth <bangerth@colostate.edu>
Sun, 2 Jul 2023 02:13:55 +0000 (20:13 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Sun, 2 Jul 2023 02:17:10 +0000 (20:17 -0600)
tests/petsc/solution_transfer_01.cc [new file with mode: 0644]
tests/petsc/solution_transfer_01.output [new file with mode: 0644]
tests/petsc/solution_transfer_02.cc [new file with mode: 0644]
tests/petsc/solution_transfer_02.output [new file with mode: 0644]

diff --git a/tests/petsc/solution_transfer_01.cc b/tests/petsc/solution_transfer_01.cc
new file mode 100644 (file)
index 0000000..5c0b92d
--- /dev/null
@@ -0,0 +1,100 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2023 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.
+//
+// ---------------------------------------------------------------------
+
+
+
+// SolutionTransfer writes into vectors, but used to forget to
+// compress the resulting vectors appropriately. This didn't matter
+// most of the time because the (sequential) SolutionTransfer was
+// apparently used only with our own vector types, for which
+// compress() is a no-op. But it failed with PETSc vectors.
+//
+// This variation of the function deals with the
+// SolutionTransfer::interpolate() variant that takes multiple vectors
+// to transfer.
+
+
+#include <deal.II/dofs/dof_accessor.h>
+#include <deal.II/dofs/dof_handler.h>
+
+#include <deal.II/fe/fe_q.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+
+#include <deal.II/lac/petsc_vector.h>
+
+#include <deal.II/numerics/solution_transfer.h>
+
+#include <iostream>
+
+#include "../tests.h"
+
+
+template <int dim>
+void
+test()
+{
+  Triangulation<dim> triangulation;
+
+  GridGenerator::hyper_cube(triangulation, 0, 1);
+  triangulation.refine_global(2);
+
+  FE_Q<dim>       fe(1);
+  DoFHandler<dim> dof_handler(triangulation);
+  dof_handler.distribute_dofs(fe);
+
+
+  PETScWrappers::MPI::Vector solution;
+  solution.reinit(dof_handler.locally_owned_dofs(), MPI_COMM_SELF);
+
+  // Do a fake transfer
+  SolutionTransfer<dim, PETScWrappers::MPI::Vector> solution_trans(dof_handler);
+
+  PETScWrappers::MPI::Vector previous_solution;
+  previous_solution = solution;
+  triangulation.prepare_coarsening_and_refinement();
+  solution_trans.prepare_for_coarsening_and_refinement(previous_solution);
+
+  triangulation.execute_coarsening_and_refinement();
+
+  solution_trans.interpolate(previous_solution, solution);
+
+  deallog << "OK" << std::endl;
+}
+
+
+
+int
+main(int argc, char **argv)
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(
+    argc, argv, testing_max_num_threads());
+
+  mpi_initlog();
+  deallog << std::setprecision(2);
+
+  {
+    deallog.push("2d");
+    test<2>();
+    deallog.pop();
+  }
+
+  {
+    deallog.push("3d");
+    test<3>();
+    deallog.pop();
+  }
+}
diff --git a/tests/petsc/solution_transfer_01.output b/tests/petsc/solution_transfer_01.output
new file mode 100644 (file)
index 0000000..25966e9
--- /dev/null
@@ -0,0 +1,3 @@
+
+DEAL:2d::OK
+DEAL:3d::OK
diff --git a/tests/petsc/solution_transfer_02.cc b/tests/petsc/solution_transfer_02.cc
new file mode 100644 (file)
index 0000000..5c0b92d
--- /dev/null
@@ -0,0 +1,100 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2023 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.
+//
+// ---------------------------------------------------------------------
+
+
+
+// SolutionTransfer writes into vectors, but used to forget to
+// compress the resulting vectors appropriately. This didn't matter
+// most of the time because the (sequential) SolutionTransfer was
+// apparently used only with our own vector types, for which
+// compress() is a no-op. But it failed with PETSc vectors.
+//
+// This variation of the function deals with the
+// SolutionTransfer::interpolate() variant that takes multiple vectors
+// to transfer.
+
+
+#include <deal.II/dofs/dof_accessor.h>
+#include <deal.II/dofs/dof_handler.h>
+
+#include <deal.II/fe/fe_q.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+
+#include <deal.II/lac/petsc_vector.h>
+
+#include <deal.II/numerics/solution_transfer.h>
+
+#include <iostream>
+
+#include "../tests.h"
+
+
+template <int dim>
+void
+test()
+{
+  Triangulation<dim> triangulation;
+
+  GridGenerator::hyper_cube(triangulation, 0, 1);
+  triangulation.refine_global(2);
+
+  FE_Q<dim>       fe(1);
+  DoFHandler<dim> dof_handler(triangulation);
+  dof_handler.distribute_dofs(fe);
+
+
+  PETScWrappers::MPI::Vector solution;
+  solution.reinit(dof_handler.locally_owned_dofs(), MPI_COMM_SELF);
+
+  // Do a fake transfer
+  SolutionTransfer<dim, PETScWrappers::MPI::Vector> solution_trans(dof_handler);
+
+  PETScWrappers::MPI::Vector previous_solution;
+  previous_solution = solution;
+  triangulation.prepare_coarsening_and_refinement();
+  solution_trans.prepare_for_coarsening_and_refinement(previous_solution);
+
+  triangulation.execute_coarsening_and_refinement();
+
+  solution_trans.interpolate(previous_solution, solution);
+
+  deallog << "OK" << std::endl;
+}
+
+
+
+int
+main(int argc, char **argv)
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(
+    argc, argv, testing_max_num_threads());
+
+  mpi_initlog();
+  deallog << std::setprecision(2);
+
+  {
+    deallog.push("2d");
+    test<2>();
+    deallog.pop();
+  }
+
+  {
+    deallog.push("3d");
+    test<3>();
+    deallog.pop();
+  }
+}
diff --git a/tests/petsc/solution_transfer_02.output b/tests/petsc/solution_transfer_02.output
new file mode 100644 (file)
index 0000000..25966e9
--- /dev/null
@@ -0,0 +1,3 @@
+
+DEAL:2d::OK
+DEAL:3d::OK

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.