1*8be712e4SBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 2*8be712e4SBarry Smith 3*8be712e4SBarry Smith /* Logging support */ 4*8be712e4SBarry Smith PetscClassId MAT_PARTITIONING_CLASSID; 5*8be712e4SBarry Smith 6*8be712e4SBarry Smith /* 7*8be712e4SBarry Smith Simplest partitioning, keeps the current partitioning. 8*8be712e4SBarry Smith */ 9*8be712e4SBarry Smith static PetscErrorCode MatPartitioningApply_Current(MatPartitioning part, IS *partitioning) 10*8be712e4SBarry Smith { 11*8be712e4SBarry Smith PetscInt m; 12*8be712e4SBarry Smith PetscMPIInt rank, size; 13*8be712e4SBarry Smith 14*8be712e4SBarry Smith PetscFunctionBegin; 15*8be712e4SBarry Smith PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)part), &size)); 16*8be712e4SBarry Smith if (part->n != size) { 17*8be712e4SBarry Smith const char *prefix; 18*8be712e4SBarry Smith PetscCall(PetscObjectGetOptionsPrefix((PetscObject)part, &prefix)); 19*8be712e4SBarry Smith SETERRQ(PetscObjectComm((PetscObject)part), PETSC_ERR_SUP, "This is the DEFAULT NO-OP partitioner, it currently only supports one domain per processor\nuse -%smat_partitioning_type parmetis or chaco or ptscotch for more than one subdomain per processor", prefix ? prefix : ""); 20*8be712e4SBarry Smith } 21*8be712e4SBarry Smith PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)part), &rank)); 22*8be712e4SBarry Smith 23*8be712e4SBarry Smith PetscCall(MatGetLocalSize(part->adj, &m, NULL)); 24*8be712e4SBarry Smith PetscCall(ISCreateStride(PetscObjectComm((PetscObject)part), m, rank, 0, partitioning)); 25*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 26*8be712e4SBarry Smith } 27*8be712e4SBarry Smith 28*8be712e4SBarry Smith /* 29*8be712e4SBarry Smith partition an index to rebalance the computation 30*8be712e4SBarry Smith */ 31*8be712e4SBarry Smith static PetscErrorCode MatPartitioningApply_Average(MatPartitioning part, IS *partitioning) 32*8be712e4SBarry Smith { 33*8be712e4SBarry Smith PetscInt m, M, nparts, *indices, r, d, *parts, i, start, end, loc; 34*8be712e4SBarry Smith 35*8be712e4SBarry Smith PetscFunctionBegin; 36*8be712e4SBarry Smith PetscCall(MatGetSize(part->adj, &M, NULL)); 37*8be712e4SBarry Smith PetscCall(MatGetLocalSize(part->adj, &m, NULL)); 38*8be712e4SBarry Smith nparts = part->n; 39*8be712e4SBarry Smith PetscCall(PetscMalloc1(nparts, &parts)); 40*8be712e4SBarry Smith d = M / nparts; 41*8be712e4SBarry Smith for (i = 0; i < nparts; i++) parts[i] = d; 42*8be712e4SBarry Smith r = M % nparts; 43*8be712e4SBarry Smith for (i = 0; i < r; i++) parts[i] += 1; 44*8be712e4SBarry Smith for (i = 1; i < nparts; i++) parts[i] += parts[i - 1]; 45*8be712e4SBarry Smith PetscCall(PetscMalloc1(m, &indices)); 46*8be712e4SBarry Smith PetscCall(MatGetOwnershipRange(part->adj, &start, &end)); 47*8be712e4SBarry Smith for (i = start; i < end; i++) { 48*8be712e4SBarry Smith PetscCall(PetscFindInt(i, nparts, parts, &loc)); 49*8be712e4SBarry Smith if (loc < 0) loc = -(loc + 1); 50*8be712e4SBarry Smith else loc = loc + 1; 51*8be712e4SBarry Smith indices[i - start] = loc; 52*8be712e4SBarry Smith } 53*8be712e4SBarry Smith PetscCall(PetscFree(parts)); 54*8be712e4SBarry Smith PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)part), m, indices, PETSC_OWN_POINTER, partitioning)); 55*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 56*8be712e4SBarry Smith } 57*8be712e4SBarry Smith 58*8be712e4SBarry Smith static PetscErrorCode MatPartitioningApply_Square(MatPartitioning part, IS *partitioning) 59*8be712e4SBarry Smith { 60*8be712e4SBarry Smith PetscInt cell, n, N, p, rstart, rend, *color; 61*8be712e4SBarry Smith PetscMPIInt size; 62*8be712e4SBarry Smith 63*8be712e4SBarry Smith PetscFunctionBegin; 64*8be712e4SBarry Smith PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)part), &size)); 65*8be712e4SBarry Smith PetscCheck(part->n == size, PetscObjectComm((PetscObject)part), PETSC_ERR_SUP, "Currently only supports one domain per processor"); 66*8be712e4SBarry Smith p = (PetscInt)PetscSqrtReal((PetscReal)part->n); 67*8be712e4SBarry Smith PetscCheck(p * p == part->n, PetscObjectComm((PetscObject)part), PETSC_ERR_SUP, "Square partitioning requires \"perfect square\" number of domains"); 68*8be712e4SBarry Smith 69*8be712e4SBarry Smith PetscCall(MatGetSize(part->adj, &N, NULL)); 70*8be712e4SBarry Smith n = (PetscInt)PetscSqrtReal((PetscReal)N); 71*8be712e4SBarry Smith PetscCheck(n * n == N, PetscObjectComm((PetscObject)part), PETSC_ERR_SUP, "Square partitioning requires square domain"); 72*8be712e4SBarry Smith PetscCheck(n % p == 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "Square partitioning requires p to divide n"); 73*8be712e4SBarry Smith PetscCall(MatGetOwnershipRange(part->adj, &rstart, &rend)); 74*8be712e4SBarry Smith PetscCall(PetscMalloc1(rend - rstart, &color)); 75*8be712e4SBarry Smith /* for (int cell=rstart; cell<rend; cell++) color[cell-rstart] = ((cell%n) < (n/2)) + 2 * ((cell/n) < (n/2)); */ 76*8be712e4SBarry Smith for (cell = rstart; cell < rend; cell++) color[cell - rstart] = ((cell % n) / (n / p)) + p * ((cell / n) / (n / p)); 77*8be712e4SBarry Smith PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)part), rend - rstart, color, PETSC_OWN_POINTER, partitioning)); 78*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 79*8be712e4SBarry Smith } 80*8be712e4SBarry Smith 81*8be712e4SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Current(MatPartitioning part) 82*8be712e4SBarry Smith { 83*8be712e4SBarry Smith PetscFunctionBegin; 84*8be712e4SBarry Smith part->ops->apply = MatPartitioningApply_Current; 85*8be712e4SBarry Smith part->ops->view = NULL; 86*8be712e4SBarry Smith part->ops->destroy = NULL; 87*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 88*8be712e4SBarry Smith } 89*8be712e4SBarry Smith 90*8be712e4SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Average(MatPartitioning part) 91*8be712e4SBarry Smith { 92*8be712e4SBarry Smith PetscFunctionBegin; 93*8be712e4SBarry Smith part->ops->apply = MatPartitioningApply_Average; 94*8be712e4SBarry Smith part->ops->view = NULL; 95*8be712e4SBarry Smith part->ops->destroy = NULL; 96*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 97*8be712e4SBarry Smith } 98*8be712e4SBarry Smith 99*8be712e4SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Square(MatPartitioning part) 100*8be712e4SBarry Smith { 101*8be712e4SBarry Smith PetscFunctionBegin; 102*8be712e4SBarry Smith part->ops->apply = MatPartitioningApply_Square; 103*8be712e4SBarry Smith part->ops->view = NULL; 104*8be712e4SBarry Smith part->ops->destroy = NULL; 105*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 106*8be712e4SBarry Smith } 107*8be712e4SBarry Smith 108*8be712e4SBarry Smith /* gets as input the "sizes" array computed by ParMetis_*_NodeND and returns 109*8be712e4SBarry Smith seps[ 0 : 2*p) : the start and end node of each subdomain 110*8be712e4SBarry Smith seps[2*p : 2*p+2*(p-1)) : the start and end node of each separator 111*8be712e4SBarry Smith levels[ 0 : p-1) : level in the tree for each separator (-1 root, -2 and -3 first level and so on) 112*8be712e4SBarry Smith The arrays must be large enough 113*8be712e4SBarry Smith */ 114*8be712e4SBarry Smith PETSC_INTERN PetscErrorCode MatPartitioningSizesToSep_Private(PetscInt p, PetscInt sizes[], PetscInt seps[], PetscInt level[]) 115*8be712e4SBarry Smith { 116*8be712e4SBarry Smith PetscInt l2p, i, pTree, pStartTree; 117*8be712e4SBarry Smith 118*8be712e4SBarry Smith PetscFunctionBegin; 119*8be712e4SBarry Smith l2p = PetscLog2Real(p); 120*8be712e4SBarry Smith PetscCheck(!(l2p - (PetscInt)PetscLog2Real(p)), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "%" PetscInt_FMT " is not a power of 2", p); 121*8be712e4SBarry Smith if (!p) PetscFunctionReturn(PETSC_SUCCESS); 122*8be712e4SBarry Smith PetscCall(PetscArrayzero(seps, 2 * p - 2)); 123*8be712e4SBarry Smith PetscCall(PetscArrayzero(level, p - 1)); 124*8be712e4SBarry Smith seps[2 * p - 2] = sizes[2 * p - 2]; 125*8be712e4SBarry Smith pTree = p; 126*8be712e4SBarry Smith pStartTree = 0; 127*8be712e4SBarry Smith while (pTree != 1) { 128*8be712e4SBarry Smith for (i = pStartTree; i < pStartTree + pTree; i++) { 129*8be712e4SBarry Smith seps[i] += sizes[i]; 130*8be712e4SBarry Smith seps[pStartTree + pTree + (i - pStartTree) / 2] += seps[i]; 131*8be712e4SBarry Smith } 132*8be712e4SBarry Smith pStartTree += pTree; 133*8be712e4SBarry Smith pTree = pTree / 2; 134*8be712e4SBarry Smith } 135*8be712e4SBarry Smith seps[2 * p - 2] -= sizes[2 * p - 2]; 136*8be712e4SBarry Smith 137*8be712e4SBarry Smith pStartTree = 2 * p - 2; 138*8be712e4SBarry Smith pTree = 1; 139*8be712e4SBarry Smith while (pStartTree > 0) { 140*8be712e4SBarry Smith for (i = pStartTree; i < pStartTree + pTree; i++) { 141*8be712e4SBarry Smith PetscInt k = 2 * i - (pStartTree + 2 * pTree); 142*8be712e4SBarry Smith PetscInt n = seps[k + 1]; 143*8be712e4SBarry Smith 144*8be712e4SBarry Smith seps[k + 1] = seps[i] - sizes[k + 1]; 145*8be712e4SBarry Smith seps[k] = seps[k + 1] + sizes[k + 1] - n - sizes[k]; 146*8be712e4SBarry Smith level[i - p] = -pTree - i + pStartTree; 147*8be712e4SBarry Smith } 148*8be712e4SBarry Smith pTree *= 2; 149*8be712e4SBarry Smith pStartTree -= pTree; 150*8be712e4SBarry Smith } 151*8be712e4SBarry Smith /* I know there should be a formula */ 152*8be712e4SBarry Smith PetscCall(PetscSortIntWithArrayPair(p - 1, seps + p, sizes + p, level)); 153*8be712e4SBarry Smith for (i = 2 * p - 2; i >= 0; i--) { 154*8be712e4SBarry Smith seps[2 * i] = seps[i]; 155*8be712e4SBarry Smith seps[2 * i + 1] = seps[i] + PetscMax(sizes[i] - 1, 0); 156*8be712e4SBarry Smith } 157*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 158*8be712e4SBarry Smith } 159*8be712e4SBarry Smith 160*8be712e4SBarry Smith PetscFunctionList MatPartitioningList = NULL; 161*8be712e4SBarry Smith PetscBool MatPartitioningRegisterAllCalled = PETSC_FALSE; 162*8be712e4SBarry Smith 163*8be712e4SBarry Smith /*@C 164*8be712e4SBarry Smith MatPartitioningRegister - Adds a new sparse matrix partitioning to the matrix package. 165*8be712e4SBarry Smith 166*8be712e4SBarry Smith Not Collective 167*8be712e4SBarry Smith 168*8be712e4SBarry Smith Input Parameters: 169*8be712e4SBarry Smith + sname - name of partitioning (for example `MATPARTITIONINGCURRENT`) or `MATPARTITIONINGPARMETIS` 170*8be712e4SBarry Smith - function - function pointer that creates the partitioning type 171*8be712e4SBarry Smith 172*8be712e4SBarry Smith Level: developer 173*8be712e4SBarry Smith 174*8be712e4SBarry Smith Example Usage: 175*8be712e4SBarry Smith .vb 176*8be712e4SBarry Smith MatPartitioningRegister("my_part", MyPartCreate); 177*8be712e4SBarry Smith .ve 178*8be712e4SBarry Smith 179*8be712e4SBarry Smith Then, your partitioner can be chosen with the procedural interface via 180*8be712e4SBarry Smith $ MatPartitioningSetType(part, "my_part") 181*8be712e4SBarry Smith or at runtime via the option 182*8be712e4SBarry Smith $ -mat_partitioning_type my_part 183*8be712e4SBarry Smith 184*8be712e4SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatPartitioning`, `MatPartitioningType`, `MatPartitioningCreate()`, `MatPartitioningRegisterDestroy()`, `MatPartitioningRegisterAll()` 185*8be712e4SBarry Smith @*/ 186*8be712e4SBarry Smith PetscErrorCode MatPartitioningRegister(const char sname[], PetscErrorCode (*function)(MatPartitioning)) 187*8be712e4SBarry Smith { 188*8be712e4SBarry Smith PetscFunctionBegin; 189*8be712e4SBarry Smith PetscCall(MatInitializePackage()); 190*8be712e4SBarry Smith PetscCall(PetscFunctionListAdd(&MatPartitioningList, sname, function)); 191*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 192*8be712e4SBarry Smith } 193*8be712e4SBarry Smith 194*8be712e4SBarry Smith /*@C 195*8be712e4SBarry Smith MatPartitioningGetType - Gets the Partitioning method type and name (as a string) 196*8be712e4SBarry Smith from the partitioning context. 197*8be712e4SBarry Smith 198*8be712e4SBarry Smith Not Collective 199*8be712e4SBarry Smith 200*8be712e4SBarry Smith Input Parameter: 201*8be712e4SBarry Smith . partitioning - the partitioning context 202*8be712e4SBarry Smith 203*8be712e4SBarry Smith Output Parameter: 204*8be712e4SBarry Smith . type - partitioner type 205*8be712e4SBarry Smith 206*8be712e4SBarry Smith Level: intermediate 207*8be712e4SBarry Smith 208*8be712e4SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatPartitioning`, `MatPartitioningType`, `MatPartitioningCreate()`, `MatPartitioningRegisterDestroy()`, `MatPartitioningRegisterAll()` 209*8be712e4SBarry Smith @*/ 210*8be712e4SBarry Smith PetscErrorCode MatPartitioningGetType(MatPartitioning partitioning, MatPartitioningType *type) 211*8be712e4SBarry Smith { 212*8be712e4SBarry Smith PetscFunctionBegin; 213*8be712e4SBarry Smith PetscValidHeaderSpecific(partitioning, MAT_PARTITIONING_CLASSID, 1); 214*8be712e4SBarry Smith PetscAssertPointer(type, 2); 215*8be712e4SBarry Smith *type = ((PetscObject)partitioning)->type_name; 216*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 217*8be712e4SBarry Smith } 218*8be712e4SBarry Smith 219*8be712e4SBarry Smith /*@C 220*8be712e4SBarry Smith MatPartitioningSetNParts - Set how many partitions need to be created; 221*8be712e4SBarry Smith by default this is one per processor. Certain partitioning schemes may 222*8be712e4SBarry Smith in fact only support that option. 223*8be712e4SBarry Smith 224*8be712e4SBarry Smith Collective 225*8be712e4SBarry Smith 226*8be712e4SBarry Smith Input Parameters: 227*8be712e4SBarry Smith + part - the partitioning context 228*8be712e4SBarry Smith - n - the number of partitions 229*8be712e4SBarry Smith 230*8be712e4SBarry Smith Level: intermediate 231*8be712e4SBarry Smith 232*8be712e4SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatPartitioning`, `MatPartitioningCreate()`, `MatPartitioningApply()` 233*8be712e4SBarry Smith @*/ 234*8be712e4SBarry Smith PetscErrorCode MatPartitioningSetNParts(MatPartitioning part, PetscInt n) 235*8be712e4SBarry Smith { 236*8be712e4SBarry Smith PetscFunctionBegin; 237*8be712e4SBarry Smith part->n = n; 238*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 239*8be712e4SBarry Smith } 240*8be712e4SBarry Smith 241*8be712e4SBarry Smith /*@ 242*8be712e4SBarry Smith MatPartitioningApplyND - Gets a nested dissection partitioning for a matrix. 243*8be712e4SBarry Smith 244*8be712e4SBarry Smith Collective 245*8be712e4SBarry Smith 246*8be712e4SBarry Smith Input Parameter: 247*8be712e4SBarry Smith . matp - the matrix partitioning object 248*8be712e4SBarry Smith 249*8be712e4SBarry Smith Output Parameter: 250*8be712e4SBarry Smith . partitioning - the partitioning. For each local node, a positive value indicates the processor 251*8be712e4SBarry Smith number the node has been assigned to. Negative x values indicate the separator level -(x+1). 252*8be712e4SBarry Smith 253*8be712e4SBarry Smith Level: intermediate 254*8be712e4SBarry Smith 255*8be712e4SBarry Smith Note: 256*8be712e4SBarry Smith The user can define additional partitionings; see `MatPartitioningRegister()`. 257*8be712e4SBarry Smith 258*8be712e4SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatPartitioningRegister()`, `MatPartitioningCreate()`, 259*8be712e4SBarry Smith `MatPartitioningDestroy()`, `MatPartitioningSetAdjacency()`, `ISPartitioningToNumbering()`, 260*8be712e4SBarry Smith `ISPartitioningCount()` 261*8be712e4SBarry Smith @*/ 262*8be712e4SBarry Smith PetscErrorCode MatPartitioningApplyND(MatPartitioning matp, IS *partitioning) 263*8be712e4SBarry Smith { 264*8be712e4SBarry Smith PetscFunctionBegin; 265*8be712e4SBarry Smith PetscValidHeaderSpecific(matp, MAT_PARTITIONING_CLASSID, 1); 266*8be712e4SBarry Smith PetscAssertPointer(partitioning, 2); 267*8be712e4SBarry Smith PetscCheck(matp->adj->assembled, PetscObjectComm((PetscObject)matp), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 268*8be712e4SBarry Smith PetscCheck(!matp->adj->factortype, PetscObjectComm((PetscObject)matp), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 269*8be712e4SBarry Smith PetscCall(PetscLogEventBegin(MAT_PartitioningND, matp, 0, 0, 0)); 270*8be712e4SBarry Smith PetscUseTypeMethod(matp, applynd, partitioning); 271*8be712e4SBarry Smith PetscCall(PetscLogEventEnd(MAT_PartitioningND, matp, 0, 0, 0)); 272*8be712e4SBarry Smith 273*8be712e4SBarry Smith PetscCall(MatPartitioningViewFromOptions(matp, NULL, "-mat_partitioning_view")); 274*8be712e4SBarry Smith PetscCall(ISViewFromOptions(*partitioning, NULL, "-mat_partitioning_view")); 275*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 276*8be712e4SBarry Smith } 277*8be712e4SBarry Smith 278*8be712e4SBarry Smith /*@ 279*8be712e4SBarry Smith MatPartitioningApply - Gets a partitioning for the graph represented by a sparse matrix. 280*8be712e4SBarry Smith 281*8be712e4SBarry Smith Collective 282*8be712e4SBarry Smith 283*8be712e4SBarry Smith Input Parameter: 284*8be712e4SBarry Smith . matp - the matrix partitioning object 285*8be712e4SBarry Smith 286*8be712e4SBarry Smith Output Parameter: 287*8be712e4SBarry Smith . partitioning - the partitioning. For each local node this tells the processor 288*8be712e4SBarry Smith number that that node is assigned to. 289*8be712e4SBarry Smith 290*8be712e4SBarry Smith Options Database Keys: 291*8be712e4SBarry Smith + -mat_partitioning_type <type> - set the partitioning package or algorithm to use 292*8be712e4SBarry Smith - -mat_partitioning_view - display information about the partitioning object 293*8be712e4SBarry Smith 294*8be712e4SBarry Smith Level: beginner 295*8be712e4SBarry Smith 296*8be712e4SBarry Smith The user can define additional partitionings; see `MatPartitioningRegister()`. 297*8be712e4SBarry Smith 298*8be712e4SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatPartitioning`, `MatPartitioningType`, `MatPartitioningRegister()`, `MatPartitioningCreate()`, 299*8be712e4SBarry Smith `MatPartitioningDestroy()`, `MatPartitioningSetAdjacency()`, `ISPartitioningToNumbering()`, 300*8be712e4SBarry Smith `ISPartitioningCount()` 301*8be712e4SBarry Smith @*/ 302*8be712e4SBarry Smith PetscErrorCode MatPartitioningApply(MatPartitioning matp, IS *partitioning) 303*8be712e4SBarry Smith { 304*8be712e4SBarry Smith PetscBool viewbalance, improve; 305*8be712e4SBarry Smith 306*8be712e4SBarry Smith PetscFunctionBegin; 307*8be712e4SBarry Smith PetscValidHeaderSpecific(matp, MAT_PARTITIONING_CLASSID, 1); 308*8be712e4SBarry Smith PetscAssertPointer(partitioning, 2); 309*8be712e4SBarry Smith PetscCheck(matp->adj->assembled, PetscObjectComm((PetscObject)matp), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 310*8be712e4SBarry Smith PetscCheck(!matp->adj->factortype, PetscObjectComm((PetscObject)matp), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 311*8be712e4SBarry Smith PetscCall(PetscLogEventBegin(MAT_Partitioning, matp, 0, 0, 0)); 312*8be712e4SBarry Smith PetscUseTypeMethod(matp, apply, partitioning); 313*8be712e4SBarry Smith PetscCall(PetscLogEventEnd(MAT_Partitioning, matp, 0, 0, 0)); 314*8be712e4SBarry Smith 315*8be712e4SBarry Smith PetscCall(MatPartitioningViewFromOptions(matp, NULL, "-mat_partitioning_view")); 316*8be712e4SBarry Smith PetscCall(ISViewFromOptions(*partitioning, NULL, "-mat_partitioning_view")); 317*8be712e4SBarry Smith 318*8be712e4SBarry Smith PetscObjectOptionsBegin((PetscObject)matp); 319*8be712e4SBarry Smith viewbalance = PETSC_FALSE; 320*8be712e4SBarry Smith PetscCall(PetscOptionsBool("-mat_partitioning_view_imbalance", "Display imbalance information of a partition", NULL, PETSC_FALSE, &viewbalance, NULL)); 321*8be712e4SBarry Smith improve = PETSC_FALSE; 322*8be712e4SBarry Smith PetscCall(PetscOptionsBool("-mat_partitioning_improve", "Improve the quality of a partition", NULL, PETSC_FALSE, &improve, NULL)); 323*8be712e4SBarry Smith PetscOptionsEnd(); 324*8be712e4SBarry Smith 325*8be712e4SBarry Smith if (improve) PetscCall(MatPartitioningImprove(matp, partitioning)); 326*8be712e4SBarry Smith 327*8be712e4SBarry Smith if (viewbalance) PetscCall(MatPartitioningViewImbalance(matp, *partitioning)); 328*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 329*8be712e4SBarry Smith } 330*8be712e4SBarry Smith 331*8be712e4SBarry Smith /*@ 332*8be712e4SBarry Smith MatPartitioningImprove - Improves the quality of a given partition. 333*8be712e4SBarry Smith 334*8be712e4SBarry Smith Collective 335*8be712e4SBarry Smith 336*8be712e4SBarry Smith Input Parameters: 337*8be712e4SBarry Smith + matp - the matrix partitioning object 338*8be712e4SBarry Smith - partitioning - the original partitioning. For each local node this tells the processor 339*8be712e4SBarry Smith number that that node is assigned to. 340*8be712e4SBarry Smith 341*8be712e4SBarry Smith Options Database Key: 342*8be712e4SBarry Smith . -mat_partitioning_improve - improve the quality of the given partition 343*8be712e4SBarry Smith 344*8be712e4SBarry Smith Level: beginner 345*8be712e4SBarry Smith 346*8be712e4SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatPartitioning`, `MatPartitioningType`, `MatPartitioningApply()`, `MatPartitioningCreate()`, 347*8be712e4SBarry Smith `MatPartitioningDestroy()`, `MatPartitioningSetAdjacency()`, `ISPartitioningToNumbering()`, 348*8be712e4SBarry Smith `ISPartitioningCount()` 349*8be712e4SBarry Smith @*/ 350*8be712e4SBarry Smith PetscErrorCode MatPartitioningImprove(MatPartitioning matp, IS *partitioning) 351*8be712e4SBarry Smith { 352*8be712e4SBarry Smith PetscFunctionBegin; 353*8be712e4SBarry Smith PetscValidHeaderSpecific(matp, MAT_PARTITIONING_CLASSID, 1); 354*8be712e4SBarry Smith PetscAssertPointer(partitioning, 2); 355*8be712e4SBarry Smith PetscCheck(matp->adj->assembled, PetscObjectComm((PetscObject)matp), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 356*8be712e4SBarry Smith PetscCheck(!matp->adj->factortype, PetscObjectComm((PetscObject)matp), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 357*8be712e4SBarry Smith PetscCall(PetscLogEventBegin(MAT_Partitioning, matp, 0, 0, 0)); 358*8be712e4SBarry Smith PetscTryTypeMethod(matp, improve, partitioning); 359*8be712e4SBarry Smith PetscCall(PetscLogEventEnd(MAT_Partitioning, matp, 0, 0, 0)); 360*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 361*8be712e4SBarry Smith } 362*8be712e4SBarry Smith 363*8be712e4SBarry Smith /*@ 364*8be712e4SBarry Smith MatPartitioningViewImbalance - Display partitioning imbalance information. 365*8be712e4SBarry Smith 366*8be712e4SBarry Smith Collective 367*8be712e4SBarry Smith 368*8be712e4SBarry Smith Input Parameters: 369*8be712e4SBarry Smith + matp - the matrix partitioning object 370*8be712e4SBarry Smith - partitioning - the partitioning. For each local node this tells the processor 371*8be712e4SBarry Smith number that that node is assigned to. 372*8be712e4SBarry Smith 373*8be712e4SBarry Smith Options Database Key: 374*8be712e4SBarry Smith . -mat_partitioning_view_balance - view the balance information from the last partitioning 375*8be712e4SBarry Smith 376*8be712e4SBarry Smith Level: beginner 377*8be712e4SBarry Smith 378*8be712e4SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatPartitioning`, `MatPartitioningType`, `MatPartitioningApply()`, `MatPartitioningView()` 379*8be712e4SBarry Smith @*/ 380*8be712e4SBarry Smith PetscErrorCode MatPartitioningViewImbalance(MatPartitioning matp, IS partitioning) 381*8be712e4SBarry Smith { 382*8be712e4SBarry Smith PetscInt nparts, *subdomainsizes, *subdomainsizes_tmp, nlocal, i, maxsub, minsub, avgsub; 383*8be712e4SBarry Smith const PetscInt *indices; 384*8be712e4SBarry Smith PetscViewer viewer; 385*8be712e4SBarry Smith 386*8be712e4SBarry Smith PetscFunctionBegin; 387*8be712e4SBarry Smith PetscValidHeaderSpecific(matp, MAT_PARTITIONING_CLASSID, 1); 388*8be712e4SBarry Smith PetscValidHeaderSpecific(partitioning, IS_CLASSID, 2); 389*8be712e4SBarry Smith nparts = matp->n; 390*8be712e4SBarry Smith PetscCall(PetscCalloc2(nparts, &subdomainsizes, nparts, &subdomainsizes_tmp)); 391*8be712e4SBarry Smith PetscCall(ISGetLocalSize(partitioning, &nlocal)); 392*8be712e4SBarry Smith PetscCall(ISGetIndices(partitioning, &indices)); 393*8be712e4SBarry Smith for (i = 0; i < nlocal; i++) subdomainsizes_tmp[indices[i]] += matp->vertex_weights ? matp->vertex_weights[i] : 1; 394*8be712e4SBarry Smith PetscCall(MPIU_Allreduce(subdomainsizes_tmp, subdomainsizes, nparts, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)matp))); 395*8be712e4SBarry Smith PetscCall(ISRestoreIndices(partitioning, &indices)); 396*8be712e4SBarry Smith minsub = PETSC_MAX_INT, maxsub = PETSC_MIN_INT, avgsub = 0; 397*8be712e4SBarry Smith for (i = 0; i < nparts; i++) { 398*8be712e4SBarry Smith minsub = PetscMin(minsub, subdomainsizes[i]); 399*8be712e4SBarry Smith maxsub = PetscMax(maxsub, subdomainsizes[i]); 400*8be712e4SBarry Smith avgsub += subdomainsizes[i]; 401*8be712e4SBarry Smith } 402*8be712e4SBarry Smith avgsub /= nparts; 403*8be712e4SBarry Smith PetscCall(PetscFree2(subdomainsizes, subdomainsizes_tmp)); 404*8be712e4SBarry Smith PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)matp), &viewer)); 405*8be712e4SBarry Smith PetscCall(MatPartitioningView(matp, viewer)); 406*8be712e4SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "Partitioning Imbalance Info: Max %" PetscInt_FMT ", Min %" PetscInt_FMT ", Avg %" PetscInt_FMT ", R %g\n", maxsub, minsub, avgsub, (double)(maxsub / (PetscReal)minsub))); 407*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 408*8be712e4SBarry Smith } 409*8be712e4SBarry Smith 410*8be712e4SBarry Smith /*@ 411*8be712e4SBarry Smith MatPartitioningSetAdjacency - Sets the adjacency graph (matrix) of the thing to be 412*8be712e4SBarry Smith partitioned. 413*8be712e4SBarry Smith 414*8be712e4SBarry Smith Collective 415*8be712e4SBarry Smith 416*8be712e4SBarry Smith Input Parameters: 417*8be712e4SBarry Smith + part - the partitioning context 418*8be712e4SBarry Smith - adj - the adjacency matrix, this can be any `MatType` but the natural representation is `MATMPIADJ` 419*8be712e4SBarry Smith 420*8be712e4SBarry Smith Level: beginner 421*8be712e4SBarry Smith 422*8be712e4SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatPartitioning`, `MatPartitioningType`, `MatPartitioningCreate()` 423*8be712e4SBarry Smith @*/ 424*8be712e4SBarry Smith PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning part, Mat adj) 425*8be712e4SBarry Smith { 426*8be712e4SBarry Smith PetscFunctionBegin; 427*8be712e4SBarry Smith PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1); 428*8be712e4SBarry Smith PetscValidHeaderSpecific(adj, MAT_CLASSID, 2); 429*8be712e4SBarry Smith part->adj = adj; 430*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 431*8be712e4SBarry Smith } 432*8be712e4SBarry Smith 433*8be712e4SBarry Smith /*@ 434*8be712e4SBarry Smith MatPartitioningDestroy - Destroys the partitioning context. 435*8be712e4SBarry Smith 436*8be712e4SBarry Smith Collective 437*8be712e4SBarry Smith 438*8be712e4SBarry Smith Input Parameter: 439*8be712e4SBarry Smith . part - the partitioning context 440*8be712e4SBarry Smith 441*8be712e4SBarry Smith Level: beginner 442*8be712e4SBarry Smith 443*8be712e4SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatPartitioning`, `MatPartitioningType`, `MatPartitioningCreate()` 444*8be712e4SBarry Smith @*/ 445*8be712e4SBarry Smith PetscErrorCode MatPartitioningDestroy(MatPartitioning *part) 446*8be712e4SBarry Smith { 447*8be712e4SBarry Smith PetscFunctionBegin; 448*8be712e4SBarry Smith if (!*part) PetscFunctionReturn(PETSC_SUCCESS); 449*8be712e4SBarry Smith PetscValidHeaderSpecific((*part), MAT_PARTITIONING_CLASSID, 1); 450*8be712e4SBarry Smith if (--((PetscObject)(*part))->refct > 0) { 451*8be712e4SBarry Smith *part = NULL; 452*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 453*8be712e4SBarry Smith } 454*8be712e4SBarry Smith 455*8be712e4SBarry Smith if ((*part)->ops->destroy) PetscCall((*(*part)->ops->destroy)((*part))); 456*8be712e4SBarry Smith PetscCall(PetscFree((*part)->vertex_weights)); 457*8be712e4SBarry Smith PetscCall(PetscFree((*part)->part_weights)); 458*8be712e4SBarry Smith PetscCall(PetscHeaderDestroy(part)); 459*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 460*8be712e4SBarry Smith } 461*8be712e4SBarry Smith 462*8be712e4SBarry Smith /*@C 463*8be712e4SBarry Smith MatPartitioningSetVertexWeights - Sets the weights for vertices for a partitioning. 464*8be712e4SBarry Smith 465*8be712e4SBarry Smith Logically Collective 466*8be712e4SBarry Smith 467*8be712e4SBarry Smith Input Parameters: 468*8be712e4SBarry Smith + part - the partitioning context 469*8be712e4SBarry Smith - weights - the weights, on each process this array must have the same size as the number of local rows times the value passed with `MatPartitioningSetNumberVertexWeights()` or 470*8be712e4SBarry Smith 1 if that is not provided 471*8be712e4SBarry Smith 472*8be712e4SBarry Smith Level: beginner 473*8be712e4SBarry Smith 474*8be712e4SBarry Smith Notes: 475*8be712e4SBarry Smith The array weights is freed by PETSc so the user should not free the array. In C/C++ 476*8be712e4SBarry Smith the array must be obtained with a call to `PetscMalloc()`, not malloc(). 477*8be712e4SBarry Smith 478*8be712e4SBarry Smith The weights may not be used by some partitioners 479*8be712e4SBarry Smith 480*8be712e4SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatPartitioning`, `MatPartitioningCreate()`, `MatPartitioningSetType()`, `MatPartitioningSetPartitionWeights()`, `MatPartitioningSetNumberVertexWeights()` 481*8be712e4SBarry Smith @*/ 482*8be712e4SBarry Smith PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning part, const PetscInt weights[]) 483*8be712e4SBarry Smith { 484*8be712e4SBarry Smith PetscFunctionBegin; 485*8be712e4SBarry Smith PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1); 486*8be712e4SBarry Smith PetscCall(PetscFree(part->vertex_weights)); 487*8be712e4SBarry Smith part->vertex_weights = (PetscInt *)weights; 488*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 489*8be712e4SBarry Smith } 490*8be712e4SBarry Smith 491*8be712e4SBarry Smith /*@C 492*8be712e4SBarry Smith MatPartitioningSetPartitionWeights - Sets the weights for each partition. 493*8be712e4SBarry Smith 494*8be712e4SBarry Smith Logically Collective 495*8be712e4SBarry Smith 496*8be712e4SBarry Smith Input Parameters: 497*8be712e4SBarry Smith + part - the partitioning context 498*8be712e4SBarry Smith - weights - An array of size nparts that is used to specify the fraction of 499*8be712e4SBarry Smith vertex weight that should be distributed to each sub-domain for 500*8be712e4SBarry Smith the balance constraint. If all of the sub-domains are to be of 501*8be712e4SBarry Smith the same size, then each of the nparts elements should be set 502*8be712e4SBarry Smith to a value of 1/nparts. Note that the sum of all of the weights 503*8be712e4SBarry Smith should be one. 504*8be712e4SBarry Smith 505*8be712e4SBarry Smith Level: beginner 506*8be712e4SBarry Smith 507*8be712e4SBarry Smith Note: 508*8be712e4SBarry Smith The array weights is freed by PETSc so the user should not free the array. In C/C++ 509*8be712e4SBarry Smith the array must be obtained with a call to `PetscMalloc()`, not malloc(). 510*8be712e4SBarry Smith 511*8be712e4SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatPartitioning`, `MatPartitioningSetVertexWeights()`, `MatPartitioningCreate()`, `MatPartitioningSetType()` 512*8be712e4SBarry Smith @*/ 513*8be712e4SBarry Smith PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning part, const PetscReal weights[]) 514*8be712e4SBarry Smith { 515*8be712e4SBarry Smith PetscFunctionBegin; 516*8be712e4SBarry Smith PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1); 517*8be712e4SBarry Smith PetscCall(PetscFree(part->part_weights)); 518*8be712e4SBarry Smith part->part_weights = (PetscReal *)weights; 519*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 520*8be712e4SBarry Smith } 521*8be712e4SBarry Smith 522*8be712e4SBarry Smith /*@ 523*8be712e4SBarry Smith MatPartitioningSetUseEdgeWeights - Set a flag to indicate whether or not to use edge weights. 524*8be712e4SBarry Smith 525*8be712e4SBarry Smith Logically Collective 526*8be712e4SBarry Smith 527*8be712e4SBarry Smith Input Parameters: 528*8be712e4SBarry Smith + part - the partitioning context 529*8be712e4SBarry Smith - use_edge_weights - the flag indicateing whether or not to use edge weights. By default no edge weights will be used, 530*8be712e4SBarry Smith that is, use_edge_weights is set to FALSE. If set use_edge_weights to TRUE, users need to make sure legal 531*8be712e4SBarry Smith edge weights are stored in an ADJ matrix. 532*8be712e4SBarry Smith 533*8be712e4SBarry Smith Options Database Key: 534*8be712e4SBarry Smith . -mat_partitioning_use_edge_weights - (true or false) 535*8be712e4SBarry Smith 536*8be712e4SBarry Smith Level: beginner 537*8be712e4SBarry Smith 538*8be712e4SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatPartitioning`, `MatPartitioningCreate()`, `MatPartitioningSetType()`, `MatPartitioningSetVertexWeights()`, `MatPartitioningSetPartitionWeights()` 539*8be712e4SBarry Smith @*/ 540*8be712e4SBarry Smith PetscErrorCode MatPartitioningSetUseEdgeWeights(MatPartitioning part, PetscBool use_edge_weights) 541*8be712e4SBarry Smith { 542*8be712e4SBarry Smith PetscFunctionBegin; 543*8be712e4SBarry Smith PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1); 544*8be712e4SBarry Smith part->use_edge_weights = use_edge_weights; 545*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 546*8be712e4SBarry Smith } 547*8be712e4SBarry Smith 548*8be712e4SBarry Smith /*@ 549*8be712e4SBarry Smith MatPartitioningGetUseEdgeWeights - Get a flag that indicates whether or not to edge weights are used. 550*8be712e4SBarry Smith 551*8be712e4SBarry Smith Logically Collective 552*8be712e4SBarry Smith 553*8be712e4SBarry Smith Input Parameter: 554*8be712e4SBarry Smith . part - the partitioning context 555*8be712e4SBarry Smith 556*8be712e4SBarry Smith Output Parameter: 557*8be712e4SBarry Smith . use_edge_weights - the flag indicateing whether or not to edge weights are used. 558*8be712e4SBarry Smith 559*8be712e4SBarry Smith Level: beginner 560*8be712e4SBarry Smith 561*8be712e4SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatPartitioning`, `MatPartitioningCreate()`, `MatPartitioningSetType()`, `MatPartitioningSetVertexWeights()`, `MatPartitioningSetPartitionWeights()`, 562*8be712e4SBarry Smith `MatPartitioningSetUseEdgeWeights` 563*8be712e4SBarry Smith @*/ 564*8be712e4SBarry Smith PetscErrorCode MatPartitioningGetUseEdgeWeights(MatPartitioning part, PetscBool *use_edge_weights) 565*8be712e4SBarry Smith { 566*8be712e4SBarry Smith PetscFunctionBegin; 567*8be712e4SBarry Smith PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1); 568*8be712e4SBarry Smith PetscAssertPointer(use_edge_weights, 2); 569*8be712e4SBarry Smith *use_edge_weights = part->use_edge_weights; 570*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 571*8be712e4SBarry Smith } 572*8be712e4SBarry Smith 573*8be712e4SBarry Smith /*@ 574*8be712e4SBarry Smith MatPartitioningCreate - Creates a partitioning context. 575*8be712e4SBarry Smith 576*8be712e4SBarry Smith Collective 577*8be712e4SBarry Smith 578*8be712e4SBarry Smith Input Parameter: 579*8be712e4SBarry Smith . comm - MPI communicator 580*8be712e4SBarry Smith 581*8be712e4SBarry Smith Output Parameter: 582*8be712e4SBarry Smith . newp - location to put the context 583*8be712e4SBarry Smith 584*8be712e4SBarry Smith Level: beginner 585*8be712e4SBarry Smith 586*8be712e4SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatPartitioning`, `MatPartitioningSetType()`, `MatPartitioningApply()`, `MatPartitioningDestroy()`, 587*8be712e4SBarry Smith `MatPartitioningSetAdjacency()` 588*8be712e4SBarry Smith @*/ 589*8be712e4SBarry Smith PetscErrorCode MatPartitioningCreate(MPI_Comm comm, MatPartitioning *newp) 590*8be712e4SBarry Smith { 591*8be712e4SBarry Smith MatPartitioning part; 592*8be712e4SBarry Smith PetscMPIInt size; 593*8be712e4SBarry Smith 594*8be712e4SBarry Smith PetscFunctionBegin; 595*8be712e4SBarry Smith *newp = NULL; 596*8be712e4SBarry Smith 597*8be712e4SBarry Smith PetscCall(MatInitializePackage()); 598*8be712e4SBarry Smith PetscCall(PetscHeaderCreate(part, MAT_PARTITIONING_CLASSID, "MatPartitioning", "Matrix/graph partitioning", "MatGraphOperations", comm, MatPartitioningDestroy, MatPartitioningView)); 599*8be712e4SBarry Smith part->vertex_weights = NULL; 600*8be712e4SBarry Smith part->part_weights = NULL; 601*8be712e4SBarry Smith part->use_edge_weights = PETSC_FALSE; /* By default we don't use edge weights */ 602*8be712e4SBarry Smith 603*8be712e4SBarry Smith PetscCallMPI(MPI_Comm_size(comm, &size)); 604*8be712e4SBarry Smith part->n = (PetscInt)size; 605*8be712e4SBarry Smith part->ncon = 1; 606*8be712e4SBarry Smith 607*8be712e4SBarry Smith *newp = part; 608*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 609*8be712e4SBarry Smith } 610*8be712e4SBarry Smith 611*8be712e4SBarry Smith /*@C 612*8be712e4SBarry Smith MatPartitioningViewFromOptions - View a partitioning context from the options database 613*8be712e4SBarry Smith 614*8be712e4SBarry Smith Collective 615*8be712e4SBarry Smith 616*8be712e4SBarry Smith Input Parameters: 617*8be712e4SBarry Smith + A - the partitioning context 618*8be712e4SBarry Smith . obj - Optional object that provides the prefix used in the options database check 619*8be712e4SBarry Smith - name - command line option 620*8be712e4SBarry Smith 621*8be712e4SBarry Smith Options Database Key: 622*8be712e4SBarry Smith . -mat_partitioning_view [viewertype]:... - the viewer and its options 623*8be712e4SBarry Smith 624*8be712e4SBarry Smith Level: intermediate 625*8be712e4SBarry Smith 626*8be712e4SBarry Smith Note: 627*8be712e4SBarry Smith .vb 628*8be712e4SBarry Smith If no value is provided ascii:stdout is used 629*8be712e4SBarry Smith ascii[:[filename][:[format][:append]]] defaults to stdout - format can be one of ascii_info, ascii_info_detail, or ascii_matlab, 630*8be712e4SBarry Smith for example ascii::ascii_info prints just the information about the object not all details 631*8be712e4SBarry Smith unless :append is given filename opens in write mode, overwriting what was already there 632*8be712e4SBarry Smith binary[:[filename][:[format][:append]]] defaults to the file binaryoutput 633*8be712e4SBarry Smith draw[:drawtype[:filename]] for example, draw:tikz, draw:tikz:figure.tex or draw:x 634*8be712e4SBarry Smith socket[:port] defaults to the standard output port 635*8be712e4SBarry Smith saws[:communicatorname] publishes object to the Scientific Application Webserver (SAWs) 636*8be712e4SBarry Smith .ve 637*8be712e4SBarry Smith 638*8be712e4SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatPartitioning`, `MatPartitioningView()`, `PetscObjectViewFromOptions()`, `MatPartitioningCreate()` 639*8be712e4SBarry Smith @*/ 640*8be712e4SBarry Smith PetscErrorCode MatPartitioningViewFromOptions(MatPartitioning A, PetscObject obj, const char name[]) 641*8be712e4SBarry Smith { 642*8be712e4SBarry Smith PetscFunctionBegin; 643*8be712e4SBarry Smith PetscValidHeaderSpecific(A, MAT_PARTITIONING_CLASSID, 1); 644*8be712e4SBarry Smith PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name)); 645*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 646*8be712e4SBarry Smith } 647*8be712e4SBarry Smith 648*8be712e4SBarry Smith /*@C 649*8be712e4SBarry Smith MatPartitioningView - Prints the partitioning data structure. 650*8be712e4SBarry Smith 651*8be712e4SBarry Smith Collective 652*8be712e4SBarry Smith 653*8be712e4SBarry Smith Input Parameters: 654*8be712e4SBarry Smith + part - the partitioning context 655*8be712e4SBarry Smith - viewer - optional visualization context 656*8be712e4SBarry Smith 657*8be712e4SBarry Smith Level: intermediate 658*8be712e4SBarry Smith 659*8be712e4SBarry Smith Note: 660*8be712e4SBarry Smith The available visualization contexts include 661*8be712e4SBarry Smith + `PETSC_VIEWER_STDOUT_SELF` - standard output (default) 662*8be712e4SBarry Smith - `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard 663*8be712e4SBarry Smith output where only the first processor opens 664*8be712e4SBarry Smith the file. All other processors send their 665*8be712e4SBarry Smith data to the first processor to print. 666*8be712e4SBarry Smith 667*8be712e4SBarry Smith The user can open alternative visualization contexts with 668*8be712e4SBarry Smith . `PetscViewerASCIIOpen()` - output to a specified file 669*8be712e4SBarry Smith 670*8be712e4SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatPartitioning`, `PetscViewer`, `PetscViewerASCIIOpen()` 671*8be712e4SBarry Smith @*/ 672*8be712e4SBarry Smith PetscErrorCode MatPartitioningView(MatPartitioning part, PetscViewer viewer) 673*8be712e4SBarry Smith { 674*8be712e4SBarry Smith PetscBool iascii; 675*8be712e4SBarry Smith 676*8be712e4SBarry Smith PetscFunctionBegin; 677*8be712e4SBarry Smith PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1); 678*8be712e4SBarry Smith if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)part), &viewer)); 679*8be712e4SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 680*8be712e4SBarry Smith PetscCheckSameComm(part, 1, viewer, 2); 681*8be712e4SBarry Smith 682*8be712e4SBarry Smith PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 683*8be712e4SBarry Smith if (iascii) { 684*8be712e4SBarry Smith PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)part, viewer)); 685*8be712e4SBarry Smith if (part->vertex_weights) PetscCall(PetscViewerASCIIPrintf(viewer, " Using vertex weights\n")); 686*8be712e4SBarry Smith } 687*8be712e4SBarry Smith PetscCall(PetscViewerASCIIPushTab(viewer)); 688*8be712e4SBarry Smith PetscTryTypeMethod(part, view, viewer); 689*8be712e4SBarry Smith PetscCall(PetscViewerASCIIPopTab(viewer)); 690*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 691*8be712e4SBarry Smith } 692*8be712e4SBarry Smith 693*8be712e4SBarry Smith /*@C 694*8be712e4SBarry Smith MatPartitioningSetType - Sets the type of partitioner to use 695*8be712e4SBarry Smith 696*8be712e4SBarry Smith Collective 697*8be712e4SBarry Smith 698*8be712e4SBarry Smith Input Parameters: 699*8be712e4SBarry Smith + part - the partitioning context. 700*8be712e4SBarry Smith - type - a known method 701*8be712e4SBarry Smith 702*8be712e4SBarry Smith Options Database Key: 703*8be712e4SBarry Smith . -mat_partitioning_type <type> - (for instance, parmetis), use -help for a list of available methods or see `MatPartitioningType` 704*8be712e4SBarry Smith 705*8be712e4SBarry Smith Level: intermediate 706*8be712e4SBarry Smith 707*8be712e4SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatPartitioning`, `MatPartitioningCreate()`, `MatPartitioningApply()`, `MatPartitioningType` 708*8be712e4SBarry Smith @*/ 709*8be712e4SBarry Smith PetscErrorCode MatPartitioningSetType(MatPartitioning part, MatPartitioningType type) 710*8be712e4SBarry Smith { 711*8be712e4SBarry Smith PetscBool match; 712*8be712e4SBarry Smith PetscErrorCode (*r)(MatPartitioning); 713*8be712e4SBarry Smith 714*8be712e4SBarry Smith PetscFunctionBegin; 715*8be712e4SBarry Smith PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1); 716*8be712e4SBarry Smith PetscAssertPointer(type, 2); 717*8be712e4SBarry Smith 718*8be712e4SBarry Smith PetscCall(PetscObjectTypeCompare((PetscObject)part, type, &match)); 719*8be712e4SBarry Smith if (match) PetscFunctionReturn(PETSC_SUCCESS); 720*8be712e4SBarry Smith 721*8be712e4SBarry Smith PetscTryTypeMethod(part, destroy); 722*8be712e4SBarry Smith part->ops->destroy = NULL; 723*8be712e4SBarry Smith 724*8be712e4SBarry Smith part->setupcalled = 0; 725*8be712e4SBarry Smith part->data = NULL; 726*8be712e4SBarry Smith PetscCall(PetscMemzero(part->ops, sizeof(struct _MatPartitioningOps))); 727*8be712e4SBarry Smith 728*8be712e4SBarry Smith PetscCall(PetscFunctionListFind(MatPartitioningList, type, &r)); 729*8be712e4SBarry Smith PetscCheck(r, PetscObjectComm((PetscObject)part), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown partitioning type %s", type); 730*8be712e4SBarry Smith 731*8be712e4SBarry Smith PetscCall((*r)(part)); 732*8be712e4SBarry Smith 733*8be712e4SBarry Smith PetscCall(PetscFree(((PetscObject)part)->type_name)); 734*8be712e4SBarry Smith PetscCall(PetscStrallocpy(type, &((PetscObject)part)->type_name)); 735*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 736*8be712e4SBarry Smith } 737*8be712e4SBarry Smith 738*8be712e4SBarry Smith /*@ 739*8be712e4SBarry Smith MatPartitioningSetFromOptions - Sets various partitioning options from the 740*8be712e4SBarry Smith options database for the partitioning object 741*8be712e4SBarry Smith 742*8be712e4SBarry Smith Collective 743*8be712e4SBarry Smith 744*8be712e4SBarry Smith Input Parameter: 745*8be712e4SBarry Smith . part - the partitioning context. 746*8be712e4SBarry Smith 747*8be712e4SBarry Smith Options Database Keys: 748*8be712e4SBarry Smith + -mat_partitioning_type <type> - (for instance, parmetis), use -help for a list of available methods 749*8be712e4SBarry Smith - -mat_partitioning_nparts - number of subgraphs 750*8be712e4SBarry Smith 751*8be712e4SBarry Smith Level: beginner 752*8be712e4SBarry Smith 753*8be712e4SBarry Smith Note: 754*8be712e4SBarry Smith If the partitioner has not been set by the user it uses one of the installed partitioner such as ParMetis. If there are 755*8be712e4SBarry Smith no installed partitioners it does no repartioning. 756*8be712e4SBarry Smith 757*8be712e4SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatPartitioning` 758*8be712e4SBarry Smith @*/ 759*8be712e4SBarry Smith PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning part) 760*8be712e4SBarry Smith { 761*8be712e4SBarry Smith PetscBool flag; 762*8be712e4SBarry Smith char type[256]; 763*8be712e4SBarry Smith const char *def; 764*8be712e4SBarry Smith 765*8be712e4SBarry Smith PetscFunctionBegin; 766*8be712e4SBarry Smith PetscObjectOptionsBegin((PetscObject)part); 767*8be712e4SBarry Smith if (!((PetscObject)part)->type_name) { 768*8be712e4SBarry Smith #if defined(PETSC_HAVE_PARMETIS) 769*8be712e4SBarry Smith def = MATPARTITIONINGPARMETIS; 770*8be712e4SBarry Smith #elif defined(PETSC_HAVE_CHACO) 771*8be712e4SBarry Smith def = MATPARTITIONINGCHACO; 772*8be712e4SBarry Smith #elif defined(PETSC_HAVE_PARTY) 773*8be712e4SBarry Smith def = MATPARTITIONINGPARTY; 774*8be712e4SBarry Smith #elif defined(PETSC_HAVE_PTSCOTCH) 775*8be712e4SBarry Smith def = MATPARTITIONINGPTSCOTCH; 776*8be712e4SBarry Smith #else 777*8be712e4SBarry Smith def = MATPARTITIONINGCURRENT; 778*8be712e4SBarry Smith #endif 779*8be712e4SBarry Smith } else { 780*8be712e4SBarry Smith def = ((PetscObject)part)->type_name; 781*8be712e4SBarry Smith } 782*8be712e4SBarry Smith PetscCall(PetscOptionsFList("-mat_partitioning_type", "Type of partitioner", "MatPartitioningSetType", MatPartitioningList, def, type, 256, &flag)); 783*8be712e4SBarry Smith if (flag) PetscCall(MatPartitioningSetType(part, type)); 784*8be712e4SBarry Smith 785*8be712e4SBarry Smith PetscCall(PetscOptionsInt("-mat_partitioning_nparts", "number of fine parts", NULL, part->n, &part->n, &flag)); 786*8be712e4SBarry Smith 787*8be712e4SBarry Smith PetscCall(PetscOptionsBool("-mat_partitioning_use_edge_weights", "whether or not to use edge weights", NULL, part->use_edge_weights, &part->use_edge_weights, &flag)); 788*8be712e4SBarry Smith 789*8be712e4SBarry Smith /* 790*8be712e4SBarry Smith Set the type if it was never set. 791*8be712e4SBarry Smith */ 792*8be712e4SBarry Smith if (!((PetscObject)part)->type_name) PetscCall(MatPartitioningSetType(part, def)); 793*8be712e4SBarry Smith 794*8be712e4SBarry Smith PetscTryTypeMethod(part, setfromoptions, PetscOptionsObject); 795*8be712e4SBarry Smith PetscOptionsEnd(); 796*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 797*8be712e4SBarry Smith } 798*8be712e4SBarry Smith 799*8be712e4SBarry Smith /*@C 800*8be712e4SBarry Smith MatPartitioningSetNumberVertexWeights - Sets the number of weights per vertex 801*8be712e4SBarry Smith 802*8be712e4SBarry Smith Not Collective 803*8be712e4SBarry Smith 804*8be712e4SBarry Smith Input Parameters: 805*8be712e4SBarry Smith + partitioning - the partitioning context 806*8be712e4SBarry Smith - ncon - the number of weights 807*8be712e4SBarry Smith 808*8be712e4SBarry Smith Level: intermediate 809*8be712e4SBarry Smith 810*8be712e4SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatPartitioning`, `MatPartitioningSetVertexWeights()` 811*8be712e4SBarry Smith @*/ 812*8be712e4SBarry Smith PetscErrorCode MatPartitioningSetNumberVertexWeights(MatPartitioning partitioning, PetscInt ncon) 813*8be712e4SBarry Smith { 814*8be712e4SBarry Smith PetscFunctionBegin; 815*8be712e4SBarry Smith PetscValidHeaderSpecific(partitioning, MAT_PARTITIONING_CLASSID, 1); 816*8be712e4SBarry Smith partitioning->ncon = ncon; 817*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 818*8be712e4SBarry Smith } 819