xref: /petsc/src/dm/impls/plex/tests/ex28.c (revision 63a3b9bc7a1f24f247904ccba9383635fe6abade)
1c4762a1bSJed Brown static char help[] = "Compare parallel partitioning strategies using matrix graphs\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown int main(int argc, char **args)
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   MatPartitioning part;
8c4762a1bSJed Brown   IS              partis;
9c4762a1bSJed Brown   Mat             A        = NULL;
10c4762a1bSJed Brown   PetscInt        max      = -1;
11c4762a1bSJed Brown   PetscInt        min      = -1;
12c4762a1bSJed Brown   PetscReal       balance  = 0.0;
13c4762a1bSJed Brown   const PetscInt *ranges  = NULL;
14c4762a1bSJed Brown   char            filein[PETSC_MAX_PATH_LEN];
15c4762a1bSJed Brown   MPI_Comm        comm;
16c4762a1bSJed Brown   PetscMPIInt     size;
17c4762a1bSJed Brown   PetscInt        p;
18c4762a1bSJed Brown   PetscBool       flg;
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   /*load matrix*/
219566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
22c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL,NULL,"-fin",filein,sizeof(filein),&flg));
25c4762a1bSJed Brown   if (flg) {
26c4762a1bSJed Brown     PetscViewer view;
279566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(comm,filein,FILE_MODE_READ,&view));
289566063dSJacob Faibussowitsch     PetscCall(MatCreate(comm,&A));
299566063dSJacob Faibussowitsch     PetscCall(MatLoad(A,view));
309566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&view));
31c4762a1bSJed Brown   }
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   /*partition matrix*/
349566063dSJacob Faibussowitsch   PetscCall(MatPartitioningCreate(comm,&part));
359566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetAdjacency(part, A));
369566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetFromOptions(part));
379566063dSJacob Faibussowitsch   PetscCall(MatPartitioningApply(part, &partis));
389566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRanges(A, &ranges));
399566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &min, NULL));
40c4762a1bSJed Brown   for (p = 0; p < size; ++p) {
41c4762a1bSJed Brown     const PetscInt partsize = ranges[p+1]-ranges[p];
42c4762a1bSJed Brown 
43c4762a1bSJed Brown     max = PetscMax(max, partsize);
44c4762a1bSJed Brown     min = PetscMin(min, partsize);
45c4762a1bSJed Brown   }
46c4762a1bSJed Brown   balance = ((PetscReal) max)/min;
479566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "ranges: "));
48c4762a1bSJed Brown   for (p = 0; p <= size; ++p) {
499566063dSJacob Faibussowitsch     if (p > 0) PetscCall(PetscPrintf(comm, ", "));
50*63a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "%" PetscInt_FMT, ranges[p]));
51c4762a1bSJed Brown   }
529566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "\n"));
539566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "max:%.0lf min:%.0lf balance:%.11lf\n", (double) max,(double) min,(double) balance));
549566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)partis,NULL,"-partition_view"));
559566063dSJacob Faibussowitsch   PetscCall(MatPartitioningDestroy(&part));
569566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&partis));
579566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
589566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
59b122ec5aSJacob Faibussowitsch   return 0;
60c4762a1bSJed Brown 
61c4762a1bSJed Brown }
62