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