]> https://gitweb.dealii.org/ - dealii.git/commitdiff
suppress PHG warning saying: PHG_EDGE_SIZE_THRESHOLD is low
authorvishalkenchan <vishal.boddu@fau.de>
Wed, 10 Jan 2018 10:42:45 +0000 (11:42 +0100)
committervishalkenchan <vishal.boddu@fau.de>
Wed, 10 Jan 2018 10:42:45 +0000 (11:42 +0100)
source/lac/sparsity_tools.cc
tests/zoltan/tria_zoltan_01.cc [new file with mode: 0644]
tests/zoltan/tria_zoltan_01.mpirun=2.output [new file with mode: 0644]
tests/zoltan/tria_zoltan_01.output [new file with mode: 0644]

index 5bf24eb8a27fee7996e5c4c462cc8af41305bc6a..ffadb21f64ffa96c9b69bbd3e62f5f5ed01864d2 100644 (file)
@@ -245,6 +245,21 @@ namespace SparsityTools
       zz->Set_Param( "LB_METHOD", "GRAPH" );  //graph based partition method (LB-load balancing)
       zz->Set_Param( "NUM_LOCAL_PARTS", std::to_string(n_partitions) ); //set number of partitions
 
+      // The PHG partitioner is a hypergraph partitioner that Zoltan could use
+      // for graph partitioning.
+      // If number of vertices in hyperedge divided by total vertices in
+      // hypergraph exceeds PHG_EDGE_SIZE_THRESHOLD,
+      // then the hyperedge will be omitted as such (dense) edges will likely
+      // incur high communication costs regardless of the partition.
+      // PHG_EDGE_SIZE_THRESHOLD value is raised to 0.5 from the default
+      // value of 0.25 so that the PHG partitioner doesn't throw warning saying
+      // "PHG_EDGE_SIZE_THRESHOLD is low ..." after removing all dense edges.
+      // For instance, in two dimensions if the triangulation being partitioned
+      // is two quadrilaterals sharing an edge and if PHG_EDGE_SIZE_THRESHOLD
+      // value is set to 0.25, PHG will remove all the edges throwing the
+      // above warning.
+      zz->Set_Param( "PHG_EDGE_SIZE_THRESHOLD", "0.5" );
+
       //Need a non-const object equal to sparsity_pattern
       SparsityPattern graph;
       graph.copy_from(sparsity_pattern);
diff --git a/tests/zoltan/tria_zoltan_01.cc b/tests/zoltan/tria_zoltan_01.cc
new file mode 100644 (file)
index 0000000..ece6ffd
--- /dev/null
@@ -0,0 +1,70 @@
+
+#include "../tests.h"
+
+#include <deal.II/base/mpi.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/distributed/shared_tria.h>
+
+#include <iostream>
+
+using namespace dealii;
+
+// Test to check whether Zoltan PHG throws warning because
+// PHG_EDGE_SIZE_THRESHOLD value is low.
+
+template <int dim>
+void test (const MPI_Comm &mpi_communicator)
+{
+  parallel::shared::Triangulation<dim>
+  triangulation (mpi_communicator,
+                 Triangulation<dim>::limit_level_difference_at_vertices);
+
+  GridGenerator::subdivided_hyper_cube (triangulation, 2, 0, 1);
+
+  // Some dummy output.
+  if (Utilities::MPI::this_mpi_process(mpi_communicator)==0)
+    std::cout << "Number of vertices: "
+              << triangulation.n_vertices()
+              << std::endl;
+}
+
+int main (int argc, char *argv[])
+{
+  try
+    {
+      dealii::Utilities::MPI::MPI_InitFinalize
+      mpi_initialization (argc,
+                          argv,
+                          dealii::numbers::invalid_unsigned_int);
+
+      test<1>(MPI_COMM_WORLD);
+      test<2>(MPI_COMM_WORLD);
+      test<3>(MPI_COMM_WORLD);
+
+    }
+  catch (std::exception &exc)
+    {
+      std::cerr << std::endl << std::endl
+                << "----------------------------------------------------"
+                << std::endl;
+      std::cerr << "Exception on processing: " << std::endl
+                << exc.what() << std::endl
+                << "Aborting!" << std::endl
+                << "----------------------------------------------------"
+                << std::endl;
+      throw;
+    }
+  catch (...)
+    {
+      std::cerr << std::endl << std::endl
+                << "----------------------------------------------------"
+                << std::endl;
+      std::cerr << "Unknown exception!" << std::endl
+                << "Aborting!" << std::endl
+                << "----------------------------------------------------"
+                << std::endl;
+      throw;
+    }
+
+  return 0;
+}
diff --git a/tests/zoltan/tria_zoltan_01.mpirun=2.output b/tests/zoltan/tria_zoltan_01.mpirun=2.output
new file mode 100644 (file)
index 0000000..744ee60
--- /dev/null
@@ -0,0 +1,3 @@
+Number of vertices: 3
+Number of vertices: 9
+Number of vertices: 27
diff --git a/tests/zoltan/tria_zoltan_01.output b/tests/zoltan/tria_zoltan_01.output
new file mode 100644 (file)
index 0000000..744ee60
--- /dev/null
@@ -0,0 +1,3 @@
+Number of vertices: 3
+Number of vertices: 9
+Number of vertices: 27

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.