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 PetscErrorCode ierr; 20 21 /*load matrix*/ 22 ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 23 comm = PETSC_COMM_WORLD; 24 CHKERRMPI(MPI_Comm_size(comm, &size)); 25 CHKERRQ(PetscOptionsGetString(NULL,NULL,"-fin",filein,sizeof(filein),&flg)); 26 if (flg) { 27 PetscViewer view; 28 CHKERRQ(PetscViewerBinaryOpen(comm,filein,FILE_MODE_READ,&view)); 29 CHKERRQ(MatCreate(comm,&A)); 30 CHKERRQ(MatLoad(A,view)); 31 CHKERRQ(PetscViewerDestroy(&view)); 32 } 33 34 /*partition matrix*/ 35 CHKERRQ(MatPartitioningCreate(comm,&part)); 36 CHKERRQ(MatPartitioningSetAdjacency(part, A)); 37 CHKERRQ(MatPartitioningSetFromOptions(part)); 38 CHKERRQ(MatPartitioningApply(part, &partis)); 39 CHKERRQ(MatGetOwnershipRanges(A, &ranges)); 40 CHKERRQ(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 CHKERRQ(PetscPrintf(comm, "ranges: ")); 49 for (p = 0; p <= size; ++p) { 50 if (p > 0) CHKERRQ(PetscPrintf(comm, ", ")); 51 CHKERRQ(PetscPrintf(comm, "%D", ranges[p])); 52 } 53 CHKERRQ(PetscPrintf(comm, "\n")); 54 CHKERRQ(PetscPrintf(comm, "max:%.0lf min:%.0lf balance:%.11lf\n", (double) max,(double) min,(double) balance)); 55 CHKERRQ(PetscObjectViewFromOptions((PetscObject)partis,NULL,"-partition_view")); 56 CHKERRQ(MatPartitioningDestroy(&part)); 57 CHKERRQ(ISDestroy(&partis)); 58 CHKERRQ(MatDestroy(&A)); 59 ierr = PetscFinalize(); 60 return ierr; 61 62 } 63