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