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