]> https://gitweb.dealii.org/ - dealii.git/commitdiff
add test 10481/head
authortcclevenger <tcleven@clemson.edu>
Mon, 8 Jun 2020 18:28:17 +0000 (14:28 -0400)
committertcclevenger <tcleven@clemson.edu>
Thu, 18 Jun 2020 21:34:38 +0000 (17:34 -0400)
doc/news/changes/minor/20200605Roth
tests/multigrid/transfer_matrix_free_13.cc
tests/multigrid/transfer_matrix_free_13.output
tests/multigrid/transfer_prebuilt_05.cc [new file with mode: 0644]
tests/multigrid/transfer_prebuilt_05.output [new file with mode: 0644]

index c5086e6e8fabdcc1646aa8479d37970be4749f31..eb48afd0e35e4a0a492b5aaec7c35e1825ff1b24 100644 (file)
@@ -1,4 +1,5 @@
 New: Applying user constraints 
 before prolongation in MGTransferPrebuilt.
 <br>
-(Julian Roth, 2020/06/05)
+(Julian Roth and Conrad Clevenger, 2020/06/05)
+
index c3954e918d45249a41a1312845d8bed8c350f37d..f677a4e5f3378ac481af687e8909fee156d5d3b8 100644 (file)
@@ -69,20 +69,23 @@ check()
   user_constraints.close();
   mg_constrained_dofs.add_user_constraints(0, user_constraints);
 
-  // build matrix-free transfer
-  MGTransferMatrixFree<dim, double> transfer(mg_constrained_dofs);
-  transfer.build(mgdof);
+  MGTransferMatrixFree<dim, double> transfer_mf(mg_constrained_dofs);
+  transfer_mf.build(mgdof);
 
+  deallog << "SRC Vector" << std::endl;
   LinearAlgebra::distributed::Vector<double> src_level_0(mgdof.n_dofs(0));
   for (unsigned int i = 0; i < mgdof.n_dofs(0); ++i)
     deallog << src_level_0(i) << " ";
-  deallog << std::endl;
-
-  LinearAlgebra::distributed::Vector<double> dst_level_1(mgdof.n_dofs(1));
-  transfer.prolongate(1, dst_level_1, src_level_0);
-  for (unsigned int i = 0; i < mgdof.n_dofs(1); ++i)
-    deallog << dst_level_1(i) << " ";
-  deallog << std::endl;
+  deallog << std::endl << std::endl;
+
+  {
+    LinearAlgebra::distributed::Vector<double> dst_level_1(mgdof.n_dofs(1));
+    transfer_mf.prolongate(1, dst_level_1, src_level_0);
+    deallog << "DST Vector" << std::endl;
+    for (unsigned int i = 0; i < mgdof.n_dofs(1); ++i)
+      deallog << dst_level_1(i) << " ";
+    deallog << std::endl;
+  }
 }
 
 
index f51ed587324e17ea3a1de4f94ead93f54aa5f4b8..91b1651826d69dc09aac1788303156bd5078b455 100644 (file)
@@ -1,7 +1,13 @@
 
+DEAL::SRC Vector
 DEAL::0.00000 0.00000 0.00000 0.00000 
+DEAL::
+DEAL::DST Vector
 DEAL::5.00000 2.50000 5.00000 2.50000 0.00000 0.00000 5.00000 2.50000 0.00000 
 DEAL::
+DEAL::SRC Vector
 DEAL::0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 
+DEAL::
+DEAL::DST Vector
 DEAL::5.00000 2.50000 5.00000 2.50000 5.00000 2.50000 5.00000 2.50000 0.00000 0.00000 0.00000 0.00000 5.00000 2.50000 5.00000 2.50000 0.00000 0.00000 5.00000 2.50000 5.00000 2.50000 0.00000 0.00000 5.00000 2.50000 0.00000 
 DEAL::
diff --git a/tests/multigrid/transfer_prebuilt_05.cc b/tests/multigrid/transfer_prebuilt_05.cc
new file mode 100644 (file)
index 0000000..119091f
--- /dev/null
@@ -0,0 +1,98 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2016 - 2019 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.
+//
+// ---------------------------------------------------------------------
+
+
+// Check MGTransferPrebuilt with custom user constraints
+
+#include <deal.II/distributed/tria.h>
+
+#include <deal.II/fe/fe_q.h>
+
+#include <deal.II/grid/grid_generator.h>
+
+#include <deal.II/multigrid/mg_transfer.h>
+
+#include "../tests.h"
+
+template <int dim>
+void
+check()
+{
+  Triangulation<dim> tr(Triangulation<dim>::limit_level_difference_at_vertices);
+  GridGenerator::hyper_cube(tr, 0, 1, false);
+
+  typename Triangulation<dim>::active_cell_iterator cell = tr.begin_active();
+  cell->face(0)->set_all_boundary_ids(1);
+
+  tr.refine_global(1);
+
+  FE_Q<dim> fe = FE_Q<dim>(1);
+
+  DoFHandler<dim> mgdof(tr);
+  mgdof.distribute_dofs(fe);
+  mgdof.distribute_mg_dofs();
+
+  MGConstrainedDoFs mg_constrained_dofs;
+  mg_constrained_dofs.initialize(mgdof);
+
+  IndexSet relevant_dofs;
+  DoFTools::extract_locally_relevant_level_dofs(mgdof, 0, relevant_dofs);
+  AffineConstraints<double> user_constraints;
+  user_constraints.reinit(relevant_dofs);
+
+  typename DoFHandler<dim>::level_face_iterator face0 = mgdof.begin(0)->face(0);
+  std::vector<types::global_dof_index>          face_dofs(fe.dofs_per_face);
+  face0->get_mg_dof_indices(0, face_dofs);
+  for (unsigned int i = 0; i < face_dofs.size(); ++i)
+    {
+      if (user_constraints.can_store_line(face_dofs[i]))
+        {
+          user_constraints.add_line(face_dofs[i]);
+          user_constraints.set_inhomogeneity(face_dofs[i], 5.0);
+        }
+    }
+  user_constraints.close();
+  mg_constrained_dofs.add_user_constraints(0, user_constraints);
+
+  MGTransferPrebuilt<Vector<double>> transfer_pb(mg_constrained_dofs);
+  transfer_pb.build(mgdof);
+
+  deallog << "SRC Vector" << std::endl;
+  Vector<double> src_level_0(mgdof.n_dofs(0));
+  for (unsigned int i = 0; i < mgdof.n_dofs(0); ++i)
+    deallog << src_level_0(i) << " ";
+  deallog << std::endl << std::endl;
+
+  {
+    Vector<double> dst_level_1(mgdof.n_dofs(1));
+    transfer_pb.prolongate(1, dst_level_1, src_level_0);
+    deallog << "DST Vector" << std::endl;
+    for (unsigned int i = 0; i < mgdof.n_dofs(1); ++i)
+      deallog << dst_level_1(i) << " ";
+    deallog << std::endl;
+  }
+}
+
+
+int
+main()
+{
+  initlog();
+
+  check<2>();
+  deallog << std::endl;
+  check<3>();
+  deallog << std::endl;
+}
diff --git a/tests/multigrid/transfer_prebuilt_05.output b/tests/multigrid/transfer_prebuilt_05.output
new file mode 100644 (file)
index 0000000..91b1651
--- /dev/null
@@ -0,0 +1,13 @@
+
+DEAL::SRC Vector
+DEAL::0.00000 0.00000 0.00000 0.00000 
+DEAL::
+DEAL::DST Vector
+DEAL::5.00000 2.50000 5.00000 2.50000 0.00000 0.00000 5.00000 2.50000 0.00000 
+DEAL::
+DEAL::SRC Vector
+DEAL::0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 
+DEAL::
+DEAL::DST Vector
+DEAL::5.00000 2.50000 5.00000 2.50000 5.00000 2.50000 5.00000 2.50000 0.00000 0.00000 0.00000 0.00000 5.00000 2.50000 5.00000 2.50000 0.00000 0.00000 5.00000 2.50000 5.00000 2.50000 0.00000 0.00000 5.00000 2.50000 0.00000 
+DEAL::

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.