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