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);
--- /dev/null
+
+#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;
+}