xref: /petsc/src/dm/impls/plex/tests/ex28.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
1 static char help[] = "Compare parallel partitioning strategies using matrix graphs\n\n";
2 
3 #include <petscmat.h>
4 
main(int argc,char ** args)5 int main(int argc, char **args)
6 {
7   MatPartitioning part;
8   IS              partis;
9   Mat             A       = NULL;
10   PetscInt        max     = -1;
11   PetscInt        min     = -1;
12   PetscReal       balance = 0.0;
13   const PetscInt *ranges  = NULL;
14   char            filein[PETSC_MAX_PATH_LEN];
15   MPI_Comm        comm;
16   PetscMPIInt     size;
17   PetscInt        p;
18   PetscBool       flg;
19 
20   /*load matrix*/
21   PetscFunctionBeginUser;
22   PetscCall(PetscInitialize(&argc, &args, NULL, help));
23   comm = PETSC_COMM_WORLD;
24   PetscCallMPI(MPI_Comm_size(comm, &size));
25   PetscCall(PetscOptionsGetString(NULL, NULL, "-fin", filein, sizeof(filein), &flg));
26   if (flg) {
27     PetscViewer view;
28     PetscCall(PetscViewerBinaryOpen(comm, filein, FILE_MODE_READ, &view));
29     PetscCall(MatCreate(comm, &A));
30     PetscCall(MatLoad(A, view));
31     PetscCall(PetscViewerDestroy(&view));
32   }
33 
34   /*partition matrix*/
35   PetscCall(MatPartitioningCreate(comm, &part));
36   PetscCall(MatPartitioningSetAdjacency(part, A));
37   PetscCall(MatPartitioningSetFromOptions(part));
38   PetscCall(MatPartitioningApply(part, &partis));
39   PetscCall(MatGetOwnershipRanges(A, &ranges));
40   PetscCall(MatGetSize(A, &min, NULL));
41   for (p = 0; p < size; ++p) {
42     const PetscInt partsize = ranges[p + 1] - ranges[p];
43 
44     max = PetscMax(max, partsize);
45     min = PetscMin(min, partsize);
46   }
47   balance = ((PetscReal)max) / min;
48   PetscCall(PetscPrintf(comm, "ranges: "));
49   for (p = 0; p <= size; ++p) {
50     if (p > 0) PetscCall(PetscPrintf(comm, ", "));
51     PetscCall(PetscPrintf(comm, "%" PetscInt_FMT, ranges[p]));
52   }
53   PetscCall(PetscPrintf(comm, "\n"));
54   PetscCall(PetscPrintf(comm, "max:%.0lf min:%.0lf balance:%.11lf\n", (double)max, (double)min, (double)balance));
55   PetscCall(PetscObjectViewFromOptions((PetscObject)partis, NULL, "-partition_view"));
56   PetscCall(MatPartitioningDestroy(&part));
57   PetscCall(ISDestroy(&partis));
58   PetscCall(MatDestroy(&A));
59   PetscCall(PetscFinalize());
60   return 0;
61 }
62