]> https://gitweb.dealii.org/ - dealii.git/commitdiff
add_boundary_indices
authorJiaqi Zhang <jiaqi2@clemson.edu>
Mon, 2 May 2022 19:15:53 +0000 (15:15 -0400)
committerJiaqi Zhang <jiaqi2@clemson.edu>
Mon, 2 May 2022 19:15:53 +0000 (15:15 -0400)
doc/news/changes/minor/20220502JiaqiZhang [new file with mode: 0644]
include/deal.II/multigrid/mg_constrained_dofs.h
tests/multigrid/constrained_dofs_01.cc
tests/multigrid/constrained_dofs_01.with_p4est=true.mpirun=1.output
tests/multigrid/constrained_dofs_01.with_p4est=true.mpirun=3.output

diff --git a/doc/news/changes/minor/20220502JiaqiZhang b/doc/news/changes/minor/20220502JiaqiZhang
new file mode 100644 (file)
index 0000000..5d812f9
--- /dev/null
@@ -0,0 +1,4 @@
+New: MGConstrainedDoFs now has a function add_boundary_indices()
+to add boundary indices on levels.
+<br>
+(Jiaqi Zhang, 2021/05/02)
index 9d9efbd31767aad3959081028f4b1f4d85aee30d..5b2605d371a975cf748d194d577cfad67d524772 100644 (file)
@@ -95,6 +95,18 @@ public:
     const std::set<types::boundary_id> &boundary_ids,
     const ComponentMask &               component_mask = ComponentMask());
 
+  /**
+   * Fill the internal data structures with information
+   * about Dirichlet boundary dofs on level @p level.
+   * The indices are restricted to the set of locally relevant
+   * level dofs.
+   */
+  template <int dim, int spacedim>
+  void
+  add_boundary_indices(const DoFHandler<dim, spacedim> &dof,
+                       const unsigned int               level,
+                       const IndexSet &                 boundary_indices);
+
   /**
    * Add user defined constraints to be used on level @p level.
    *
@@ -355,6 +367,25 @@ MGConstrainedDoFs::make_zero_boundary_constraints(
 
 
 
+template <int dim, int spacedim>
+inline void
+MGConstrainedDoFs::add_boundary_indices(const DoFHandler<dim, spacedim> &dof,
+                                        const unsigned int               level,
+                                        const IndexSet &level_boundary_indices)
+{
+  const unsigned int n_levels = dof.get_triangulation().n_global_levels();
+  if (boundary_indices.size() == 0)
+    {
+      boundary_indices.resize(n_levels);
+      for (unsigned int i = 0; i < n_levels; ++i)
+        boundary_indices[i] = IndexSet(dof.n_dofs(i));
+    }
+  AssertDimension(boundary_indices.size(), n_levels);
+  boundary_indices[level].add_indices(level_boundary_indices);
+}
+
+
+
 template <int dim, int spacedim>
 inline void
 MGConstrainedDoFs::make_no_normal_flux_constraints(
index 843aa2fc51b98659e7190f2dfada870c1f09bca0..d8cf2578da9f6fcee38df07a761f486dc51fb0ab 100644 (file)
@@ -166,6 +166,80 @@ check_fe(FiniteElement<dim> &fe)
               "FAIL ")
         << std::endl;
     }
+
+  {
+    deallog << "Test add_boundary_indices: " << std::endl;
+    // used to test add_boundary_indices()
+    MGConstrainedDoFs mg_constrained_dofs_1;
+    mg_constrained_dofs_1.initialize(dofh);
+
+    MGLevelObject<AffineConstraints<double>> mg_boundary_constraints;
+    mg_boundary_constraints.resize(0, n_levels - 1);
+
+    std::vector<IndexSet> boundary_indices;
+    boundary_indices.resize(n_levels);
+    MGTools::make_boundary_list(dofh, {0}, boundary_indices);
+
+    for (unsigned int level = 0; level < n_levels; ++level)
+      {
+        deallog << "Level " << level << ':' << std::endl;
+
+        IndexSet relevant;
+        DoFTools::extract_locally_relevant_level_dofs(dofh, level, relevant);
+        mg_boundary_constraints[level].reinit(relevant);
+        mg_boundary_constraints[level].add_lines(boundary_indices[level]);
+
+        mg_constrained_dofs_1.add_boundary_indices(dofh,
+                                                   level,
+                                                   boundary_indices[level]);
+        const auto &bi = mg_constrained_dofs_1.get_boundary_indices(level);
+
+        deallog << "get_boundary_indices test:" << std::endl;
+        bi.print(deallog);
+
+        deallog << "relevant:" << std::endl;
+        relevant.print(deallog);
+
+        deallog << ((bi ==
+                     (relevant &
+                      mg_constrained_dofs_ref.get_boundary_indices(level))) ?
+                      "ok " :
+                      "FAIL test")
+                << std::endl;
+      }
+
+    // extract boundary indices from constraint matrices
+    // this is probably how we would use the function add_boundary_indices()
+    mg_constrained_dofs_1.clear();
+    mg_constrained_dofs_1.initialize(dofh);
+    for (unsigned int level = 0; level < n_levels; ++level)
+      {
+        deallog << "Level " << level << ':' << std::endl;
+
+        IndexSet level_boundary_indices(dofh.n_dofs(level));
+        for (auto line : mg_boundary_constraints[level].get_lines())
+          {
+            if (line.entries.size() == 0)
+              {
+                level_boundary_indices.add_index(line.index);
+              }
+          }
+        mg_constrained_dofs_1.add_boundary_indices(dofh,
+                                                   level,
+                                                   level_boundary_indices);
+        const auto &bi = mg_constrained_dofs_1.get_boundary_indices(level);
+
+        IndexSet relevant;
+        DoFTools::extract_locally_relevant_level_dofs(dofh, level, relevant);
+
+        deallog << ((bi ==
+                     (relevant &
+                      mg_constrained_dofs_ref.get_boundary_indices(level))) ?
+                      "ok " :
+                      "FAIL test")
+                << std::endl;
+      }
+  }
 }
 
 
index c22d3d278e11a9bccb16f281e8651a505a652614..d2902b294d479548c7e9414d078578240c89c2be 100644 (file)
@@ -32,3 +32,36 @@ DEAL:0::{[0,1], [4,5], 8}
 DEAL:0::relevant:
 DEAL:0::{[0,8]}
 DEAL:0::ok ok 
+DEAL:0::Test add_boundary_indices: 
+DEAL:0::Level 0:
+DEAL:0::get_boundary_indices test:
+DEAL:0::{[0,3]}
+DEAL:0::relevant:
+DEAL:0::{[0,3]}
+DEAL:0::ok 
+DEAL:0::Level 1:
+DEAL:0::get_boundary_indices test:
+DEAL:0::{[0,2], [4,8]}
+DEAL:0::relevant:
+DEAL:0::{[0,8]}
+DEAL:0::ok 
+DEAL:0::Level 2:
+DEAL:0::get_boundary_indices test:
+DEAL:0::{[0,2], 4, 6, 9, [11,12], [14,15], [18,20], [22,24]}
+DEAL:0::relevant:
+DEAL:0::{[0,24]}
+DEAL:0::ok 
+DEAL:0::Level 3:
+DEAL:0::get_boundary_indices test:
+DEAL:0::{[0,1], [4,5], 8}
+DEAL:0::relevant:
+DEAL:0::{[0,8]}
+DEAL:0::ok 
+DEAL:0::Level 0:
+DEAL:0::ok 
+DEAL:0::Level 1:
+DEAL:0::ok 
+DEAL:0::Level 2:
+DEAL:0::ok 
+DEAL:0::Level 3:
+DEAL:0::ok 
index e45b5c97b37f951d330e434971b5559eff2daac1..6bdcab37cc26fd38e879eef927bd341a8a445398 100644 (file)
@@ -32,6 +32,39 @@ DEAL:0::{[0,1]}
 DEAL:0::relevant:
 DEAL:0::{[0,3], [6,7]}
 DEAL:0::ok ok 
+DEAL:0::Test add_boundary_indices: 
+DEAL:0::Level 0:
+DEAL:0::get_boundary_indices test:
+DEAL:0::{[0,3]}
+DEAL:0::relevant:
+DEAL:0::{[0,3]}
+DEAL:0::ok 
+DEAL:0::Level 1:
+DEAL:0::get_boundary_indices test:
+DEAL:0::{[0,2], [4,8]}
+DEAL:0::relevant:
+DEAL:0::{[0,8]}
+DEAL:0::ok 
+DEAL:0::Level 2:
+DEAL:0::get_boundary_indices test:
+DEAL:0::{[0,2], 4, 6, 9, [11,12], [14,15]}
+DEAL:0::relevant:
+DEAL:0::{[0,17], 21}
+DEAL:0::ok 
+DEAL:0::Level 3:
+DEAL:0::get_boundary_indices test:
+DEAL:0::{[0,1]}
+DEAL:0::relevant:
+DEAL:0::{[0,3], [6,7]}
+DEAL:0::ok 
+DEAL:0::Level 0:
+DEAL:0::ok 
+DEAL:0::Level 1:
+DEAL:0::ok 
+DEAL:0::Level 2:
+DEAL:0::ok 
+DEAL:0::Level 3:
+DEAL:0::ok 
 
 DEAL:1::FE_Q<2>(1)
 DEAL:1::Level 0:
@@ -66,6 +99,39 @@ DEAL:1::{[0,1], [4,5], 8}
 DEAL:1::relevant:
 DEAL:1::{[0,8]}
 DEAL:1::ok ok 
+DEAL:1::Test add_boundary_indices: 
+DEAL:1::Level 0:
+DEAL:1::get_boundary_indices test:
+DEAL:1::{}
+DEAL:1::relevant:
+DEAL:1::{}
+DEAL:1::ok 
+DEAL:1::Level 1:
+DEAL:1::get_boundary_indices test:
+DEAL:1::{1, [4,5]}
+DEAL:1::relevant:
+DEAL:1::{1, [3,5]}
+DEAL:1::ok 
+DEAL:1::Level 2:
+DEAL:1::get_boundary_indices test:
+DEAL:1::{1, 4, 9, [11,12], 14, 22}
+DEAL:1::relevant:
+DEAL:1::{1, [3,5], [7,14], [16,17], [21,22]}
+DEAL:1::ok 
+DEAL:1::Level 3:
+DEAL:1::get_boundary_indices test:
+DEAL:1::{[0,1], [4,5], 8}
+DEAL:1::relevant:
+DEAL:1::{[0,8]}
+DEAL:1::ok 
+DEAL:1::Level 0:
+DEAL:1::ok 
+DEAL:1::Level 1:
+DEAL:1::ok 
+DEAL:1::Level 2:
+DEAL:1::ok 
+DEAL:1::Level 3:
+DEAL:1::ok 
 
 
 DEAL:2::FE_Q<2>(1)
@@ -101,4 +167,37 @@ DEAL:2::{}
 DEAL:2::relevant:
 DEAL:2::{}
 DEAL:2::ok ok 
+DEAL:2::Test add_boundary_indices: 
+DEAL:2::Level 0:
+DEAL:2::get_boundary_indices test:
+DEAL:2::{[0,3]}
+DEAL:2::relevant:
+DEAL:2::{[0,3]}
+DEAL:2::ok 
+DEAL:2::Level 1:
+DEAL:2::get_boundary_indices test:
+DEAL:2::{[0,2], [4,8]}
+DEAL:2::relevant:
+DEAL:2::{[0,8]}
+DEAL:2::ok 
+DEAL:2::Level 2:
+DEAL:2::get_boundary_indices test:
+DEAL:2::{2, 6, 12, [14,15], [18,20], [22,24]}
+DEAL:2::relevant:
+DEAL:2::{[2,3], [5,8], 10, [12,24]}
+DEAL:2::ok 
+DEAL:2::Level 3:
+DEAL:2::get_boundary_indices test:
+DEAL:2::{}
+DEAL:2::relevant:
+DEAL:2::{}
+DEAL:2::ok 
+DEAL:2::Level 0:
+DEAL:2::ok 
+DEAL:2::Level 1:
+DEAL:2::ok 
+DEAL:2::Level 2:
+DEAL:2::ok 
+DEAL:2::Level 3:
+DEAL:2::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.