18be712e4SBarry Smith #include <../src/mat/impls/adj/mpi/mpiadj.h> /*I "petscmat.h" I*/
28be712e4SBarry Smith
38be712e4SBarry Smith EXTERN_C_BEGIN
48be712e4SBarry Smith #include <ptscotch.h>
58be712e4SBarry Smith #if defined(PETSC_HAVE_SCOTCH_PARMETIS_V3_NODEND)
68be712e4SBarry Smith /* we define the prototype instead of include SCOTCH's parmetis.h */
78be712e4SBarry Smith void SCOTCH_ParMETIS_V3_NodeND(const SCOTCH_Num *const, SCOTCH_Num *const, SCOTCH_Num *const, const SCOTCH_Num *const, const SCOTCH_Num *const, SCOTCH_Num *const, SCOTCH_Num *const, MPI_Comm *);
88be712e4SBarry Smith #endif
98be712e4SBarry Smith EXTERN_C_END
108be712e4SBarry Smith
118be712e4SBarry Smith typedef struct {
128be712e4SBarry Smith double imbalance;
138be712e4SBarry Smith SCOTCH_Num strategy;
148be712e4SBarry Smith } MatPartitioning_PTScotch;
158be712e4SBarry Smith
168be712e4SBarry Smith /*@
178be712e4SBarry Smith MatPartitioningPTScotchSetImbalance - Sets the value of the load imbalance
188be712e4SBarry Smith ratio to be used during strategy selection.
198be712e4SBarry Smith
208be712e4SBarry Smith Collective
218be712e4SBarry Smith
228be712e4SBarry Smith Input Parameters:
238be712e4SBarry Smith + part - the partitioning context
248be712e4SBarry Smith - imb - the load imbalance ratio
258be712e4SBarry Smith
268be712e4SBarry Smith Options Database Key:
278be712e4SBarry Smith . -mat_partitioning_ptscotch_imbalance <imb> - set load imbalance ratio
288be712e4SBarry Smith
298be712e4SBarry Smith Note:
308be712e4SBarry Smith Must be in the range [0,1]. The default value is 0.01.
318be712e4SBarry Smith
328be712e4SBarry Smith Level: advanced
338be712e4SBarry Smith
348be712e4SBarry Smith .seealso: `MATPARTITIONINGSCOTCH`, `MatPartitioningPTScotchSetStrategy()`, `MatPartitioningPTScotchGetImbalance()`
358be712e4SBarry Smith @*/
MatPartitioningPTScotchSetImbalance(MatPartitioning part,PetscReal imb)368be712e4SBarry Smith PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning part, PetscReal imb)
378be712e4SBarry Smith {
388be712e4SBarry Smith PetscFunctionBegin;
398be712e4SBarry Smith PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
408be712e4SBarry Smith PetscValidLogicalCollectiveReal(part, imb, 2);
418be712e4SBarry Smith PetscTryMethod(part, "MatPartitioningPTScotchSetImbalance_C", (MatPartitioning, PetscReal), (part, imb));
428be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
438be712e4SBarry Smith }
448be712e4SBarry Smith
MatPartitioningPTScotchSetImbalance_PTScotch(MatPartitioning part,PetscReal imb)458be712e4SBarry Smith static PetscErrorCode MatPartitioningPTScotchSetImbalance_PTScotch(MatPartitioning part, PetscReal imb)
468be712e4SBarry Smith {
478be712e4SBarry Smith MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch *)part->data;
488be712e4SBarry Smith
498be712e4SBarry Smith PetscFunctionBegin;
508be712e4SBarry Smith if (imb == PETSC_DEFAULT) scotch->imbalance = 0.01;
518be712e4SBarry Smith else {
528be712e4SBarry Smith PetscCheck(imb >= 0.0 && imb <= 1.0, PetscObjectComm((PetscObject)part), PETSC_ERR_ARG_OUTOFRANGE, "Illegal value of imb. Must be in range [0,1]");
538be712e4SBarry Smith scotch->imbalance = (double)imb;
548be712e4SBarry Smith }
558be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
568be712e4SBarry Smith }
578be712e4SBarry Smith
588be712e4SBarry Smith /*@
598be712e4SBarry Smith MatPartitioningPTScotchGetImbalance - Gets the value of the load imbalance
608be712e4SBarry Smith ratio used during strategy selection.
618be712e4SBarry Smith
628be712e4SBarry Smith Not Collective
638be712e4SBarry Smith
648be712e4SBarry Smith Input Parameter:
658be712e4SBarry Smith . part - the partitioning context
668be712e4SBarry Smith
678be712e4SBarry Smith Output Parameter:
688be712e4SBarry Smith . imb - the load imbalance ratio
698be712e4SBarry Smith
708be712e4SBarry Smith Level: advanced
718be712e4SBarry Smith
728be712e4SBarry Smith .seealso: `MATPARTITIONINGSCOTCH`, `MatPartitioningPTScotchSetImbalance()`
738be712e4SBarry Smith @*/
MatPartitioningPTScotchGetImbalance(MatPartitioning part,PetscReal * imb)748be712e4SBarry Smith PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning part, PetscReal *imb)
758be712e4SBarry Smith {
768be712e4SBarry Smith PetscFunctionBegin;
778be712e4SBarry Smith PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
788be712e4SBarry Smith PetscAssertPointer(imb, 2);
798be712e4SBarry Smith PetscUseMethod(part, "MatPartitioningPTScotchGetImbalance_C", (MatPartitioning, PetscReal *), (part, imb));
808be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
818be712e4SBarry Smith }
828be712e4SBarry Smith
MatPartitioningPTScotchGetImbalance_PTScotch(MatPartitioning part,PetscReal * imb)838be712e4SBarry Smith static PetscErrorCode MatPartitioningPTScotchGetImbalance_PTScotch(MatPartitioning part, PetscReal *imb)
848be712e4SBarry Smith {
858be712e4SBarry Smith MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch *)part->data;
868be712e4SBarry Smith
878be712e4SBarry Smith PetscFunctionBegin;
888be712e4SBarry Smith *imb = scotch->imbalance;
898be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
908be712e4SBarry Smith }
918be712e4SBarry Smith
928be712e4SBarry Smith /*@
938be712e4SBarry Smith MatPartitioningPTScotchSetStrategy - Sets the strategy to be used in PTScotch.
948be712e4SBarry Smith
958be712e4SBarry Smith Collective
968be712e4SBarry Smith
978be712e4SBarry Smith Input Parameters:
988be712e4SBarry Smith + part - the partitioning context
998be712e4SBarry Smith - strategy - the strategy, one of
1008be712e4SBarry Smith .vb
1018be712e4SBarry Smith MP_PTSCOTCH_DEFAULT - Default behavior
1028be712e4SBarry Smith MP_PTSCOTCH_QUALITY - Prioritize quality over speed
1038be712e4SBarry Smith MP_PTSCOTCH_SPEED - Prioritize speed over quality
1048be712e4SBarry Smith MP_PTSCOTCH_BALANCE - Enforce load balance
1058be712e4SBarry Smith MP_PTSCOTCH_SAFETY - Avoid methods that may fail
1068be712e4SBarry Smith MP_PTSCOTCH_SCALABILITY - Favor scalability as much as possible
1078be712e4SBarry Smith .ve
1088be712e4SBarry Smith
1098be712e4SBarry Smith Options Database Key:
1108be712e4SBarry Smith . -mat_partitioning_ptscotch_strategy [quality,speed,balance,safety,scalability] - strategy
1118be712e4SBarry Smith
1128be712e4SBarry Smith Level: advanced
1138be712e4SBarry Smith
1148be712e4SBarry Smith Note:
1158be712e4SBarry Smith The default is `MP_SCOTCH_QUALITY`. See the PTScotch documentation for more information.
1168be712e4SBarry Smith
1178be712e4SBarry Smith .seealso: `MATPARTITIONINGSCOTCH`, `MatPartitioningPTScotchSetImbalance()`, `MatPartitioningPTScotchGetStrategy()`
1188be712e4SBarry Smith @*/
MatPartitioningPTScotchSetStrategy(MatPartitioning part,MPPTScotchStrategyType strategy)1198be712e4SBarry Smith PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning part, MPPTScotchStrategyType strategy)
1208be712e4SBarry Smith {
1218be712e4SBarry Smith PetscFunctionBegin;
1228be712e4SBarry Smith PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
1238be712e4SBarry Smith PetscValidLogicalCollectiveEnum(part, strategy, 2);
1248be712e4SBarry Smith PetscTryMethod(part, "MatPartitioningPTScotchSetStrategy_C", (MatPartitioning, MPPTScotchStrategyType), (part, strategy));
1258be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
1268be712e4SBarry Smith }
1278be712e4SBarry Smith
MatPartitioningPTScotchSetStrategy_PTScotch(MatPartitioning part,MPPTScotchStrategyType strategy)1288be712e4SBarry Smith static PetscErrorCode MatPartitioningPTScotchSetStrategy_PTScotch(MatPartitioning part, MPPTScotchStrategyType strategy)
1298be712e4SBarry Smith {
1308be712e4SBarry Smith MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch *)part->data;
1318be712e4SBarry Smith
1328be712e4SBarry Smith PetscFunctionBegin;
1338be712e4SBarry Smith switch (strategy) {
1348be712e4SBarry Smith case MP_PTSCOTCH_QUALITY:
1358be712e4SBarry Smith scotch->strategy = SCOTCH_STRATQUALITY;
1368be712e4SBarry Smith break;
1378be712e4SBarry Smith case MP_PTSCOTCH_SPEED:
1388be712e4SBarry Smith scotch->strategy = SCOTCH_STRATSPEED;
1398be712e4SBarry Smith break;
1408be712e4SBarry Smith case MP_PTSCOTCH_BALANCE:
1418be712e4SBarry Smith scotch->strategy = SCOTCH_STRATBALANCE;
1428be712e4SBarry Smith break;
1438be712e4SBarry Smith case MP_PTSCOTCH_SAFETY:
1448be712e4SBarry Smith scotch->strategy = SCOTCH_STRATSAFETY;
1458be712e4SBarry Smith break;
1468be712e4SBarry Smith case MP_PTSCOTCH_SCALABILITY:
1478be712e4SBarry Smith scotch->strategy = SCOTCH_STRATSCALABILITY;
1488be712e4SBarry Smith break;
1498be712e4SBarry Smith default:
1508be712e4SBarry Smith scotch->strategy = SCOTCH_STRATDEFAULT;
1518be712e4SBarry Smith break;
1528be712e4SBarry Smith }
1538be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
1548be712e4SBarry Smith }
1558be712e4SBarry Smith
1568be712e4SBarry Smith /*@
1578be712e4SBarry Smith MatPartitioningPTScotchGetStrategy - Gets the strategy used in PTScotch.
1588be712e4SBarry Smith
1598be712e4SBarry Smith Not Collective
1608be712e4SBarry Smith
1618be712e4SBarry Smith Input Parameter:
1628be712e4SBarry Smith . part - the partitioning context
1638be712e4SBarry Smith
1648be712e4SBarry Smith Output Parameter:
1658be712e4SBarry Smith . strategy - the strategy
1668be712e4SBarry Smith
1678be712e4SBarry Smith Level: advanced
1688be712e4SBarry Smith
1698be712e4SBarry Smith .seealso: `MATPARTITIONINGSCOTCH`, `MatPartitioningPTScotchSetStrategy()`
1708be712e4SBarry Smith @*/
MatPartitioningPTScotchGetStrategy(MatPartitioning part,MPPTScotchStrategyType * strategy)1718be712e4SBarry Smith PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning part, MPPTScotchStrategyType *strategy)
1728be712e4SBarry Smith {
1738be712e4SBarry Smith PetscFunctionBegin;
1748be712e4SBarry Smith PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
1758be712e4SBarry Smith PetscAssertPointer(strategy, 2);
1768be712e4SBarry Smith PetscUseMethod(part, "MatPartitioningPTScotchGetStrategy_C", (MatPartitioning, MPPTScotchStrategyType *), (part, strategy));
1778be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
1788be712e4SBarry Smith }
1798be712e4SBarry Smith
MatPartitioningPTScotchGetStrategy_PTScotch(MatPartitioning part,MPPTScotchStrategyType * strategy)1808be712e4SBarry Smith static PetscErrorCode MatPartitioningPTScotchGetStrategy_PTScotch(MatPartitioning part, MPPTScotchStrategyType *strategy)
1818be712e4SBarry Smith {
1828be712e4SBarry Smith MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch *)part->data;
1838be712e4SBarry Smith
1848be712e4SBarry Smith PetscFunctionBegin;
1858be712e4SBarry Smith switch (scotch->strategy) {
1868be712e4SBarry Smith case SCOTCH_STRATQUALITY:
1878be712e4SBarry Smith *strategy = MP_PTSCOTCH_QUALITY;
1888be712e4SBarry Smith break;
1898be712e4SBarry Smith case SCOTCH_STRATSPEED:
1908be712e4SBarry Smith *strategy = MP_PTSCOTCH_SPEED;
1918be712e4SBarry Smith break;
1928be712e4SBarry Smith case SCOTCH_STRATBALANCE:
1938be712e4SBarry Smith *strategy = MP_PTSCOTCH_BALANCE;
1948be712e4SBarry Smith break;
1958be712e4SBarry Smith case SCOTCH_STRATSAFETY:
1968be712e4SBarry Smith *strategy = MP_PTSCOTCH_SAFETY;
1978be712e4SBarry Smith break;
1988be712e4SBarry Smith case SCOTCH_STRATSCALABILITY:
1998be712e4SBarry Smith *strategy = MP_PTSCOTCH_SCALABILITY;
2008be712e4SBarry Smith break;
2018be712e4SBarry Smith default:
2028be712e4SBarry Smith *strategy = MP_PTSCOTCH_DEFAULT;
2038be712e4SBarry Smith break;
2048be712e4SBarry Smith }
2058be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
2068be712e4SBarry Smith }
2078be712e4SBarry Smith
MatPartitioningView_PTScotch(MatPartitioning part,PetscViewer viewer)2088be712e4SBarry Smith static PetscErrorCode MatPartitioningView_PTScotch(MatPartitioning part, PetscViewer viewer)
2098be712e4SBarry Smith {
2108be712e4SBarry Smith MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch *)part->data;
2118be712e4SBarry Smith PetscBool isascii;
2128be712e4SBarry Smith const char *str = NULL;
2138be712e4SBarry Smith
2148be712e4SBarry Smith PetscFunctionBegin;
2158be712e4SBarry Smith PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
2168be712e4SBarry Smith if (isascii) {
2178be712e4SBarry Smith switch (scotch->strategy) {
2188be712e4SBarry Smith case SCOTCH_STRATQUALITY:
2198be712e4SBarry Smith str = "Prioritize quality over speed";
2208be712e4SBarry Smith break;
2218be712e4SBarry Smith case SCOTCH_STRATSPEED:
2228be712e4SBarry Smith str = "Prioritize speed over quality";
2238be712e4SBarry Smith break;
2248be712e4SBarry Smith case SCOTCH_STRATBALANCE:
2258be712e4SBarry Smith str = "Enforce load balance";
2268be712e4SBarry Smith break;
2278be712e4SBarry Smith case SCOTCH_STRATSAFETY:
2288be712e4SBarry Smith str = "Avoid methods that may fail";
2298be712e4SBarry Smith break;
2308be712e4SBarry Smith case SCOTCH_STRATSCALABILITY:
2318be712e4SBarry Smith str = "Favor scalability as much as possible";
2328be712e4SBarry Smith break;
2338be712e4SBarry Smith default:
2348be712e4SBarry Smith str = "Default behavior";
2358be712e4SBarry Smith break;
2368be712e4SBarry Smith }
2378be712e4SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, " Strategy=%s\n", str));
2388be712e4SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, " Load imbalance ratio=%g\n", scotch->imbalance));
2398be712e4SBarry Smith }
2408be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
2418be712e4SBarry Smith }
2428be712e4SBarry Smith
MatPartitioningSetFromOptions_PTScotch(MatPartitioning part,PetscOptionItems PetscOptionsObject)243*ce78bad3SBarry Smith static PetscErrorCode MatPartitioningSetFromOptions_PTScotch(MatPartitioning part, PetscOptionItems PetscOptionsObject)
2448be712e4SBarry Smith {
2458be712e4SBarry Smith PetscBool flag;
2468be712e4SBarry Smith PetscReal r;
2478be712e4SBarry Smith MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch *)part->data;
2488be712e4SBarry Smith MPPTScotchStrategyType strat;
2498be712e4SBarry Smith
2508be712e4SBarry Smith PetscFunctionBegin;
2518be712e4SBarry Smith PetscCall(MatPartitioningPTScotchGetStrategy(part, &strat));
2528be712e4SBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "PTScotch partitioning options");
2538be712e4SBarry Smith PetscCall(PetscOptionsEnum("-mat_partitioning_ptscotch_strategy", "Strategy", "MatPartitioningPTScotchSetStrategy", MPPTScotchStrategyTypes, (PetscEnum)strat, (PetscEnum *)&strat, &flag));
2548be712e4SBarry Smith if (flag) PetscCall(MatPartitioningPTScotchSetStrategy(part, strat));
2558be712e4SBarry Smith PetscCall(PetscOptionsReal("-mat_partitioning_ptscotch_imbalance", "Load imbalance ratio", "MatPartitioningPTScotchSetImbalance", scotch->imbalance, &r, &flag));
2568be712e4SBarry Smith if (flag) PetscCall(MatPartitioningPTScotchSetImbalance(part, r));
2578be712e4SBarry Smith PetscOptionsHeadEnd();
2588be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
2598be712e4SBarry Smith }
2608be712e4SBarry Smith
MatPartitioningApply_PTScotch_Private(MatPartitioning part,PetscBool useND,IS * partitioning)2618be712e4SBarry Smith static PetscErrorCode MatPartitioningApply_PTScotch_Private(MatPartitioning part, PetscBool useND, IS *partitioning)
2628be712e4SBarry Smith {
2638be712e4SBarry Smith MPI_Comm pcomm, comm;
2648be712e4SBarry Smith MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch *)part->data;
2658be712e4SBarry Smith PetscMPIInt rank;
2668be712e4SBarry Smith Mat mat = part->adj;
2678be712e4SBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj *)mat->data;
2688be712e4SBarry Smith PetscBool flg, distributed;
2698be712e4SBarry Smith PetscBool proc_weight_flg;
2708be712e4SBarry Smith PetscInt i, j, p, bs = 1, nold;
2718be712e4SBarry Smith PetscInt *NDorder = NULL;
2728be712e4SBarry Smith PetscReal *vwgttab, deltval;
2738be712e4SBarry Smith SCOTCH_Num *locals, *velotab, *veloloctab, *edloloctab, vertlocnbr, edgelocnbr, nparts = part->n;
2748be712e4SBarry Smith
2758be712e4SBarry Smith PetscFunctionBegin;
2768be712e4SBarry Smith PetscCall(PetscObjectGetComm((PetscObject)part, &pcomm));
2778be712e4SBarry Smith /* Duplicate the communicator to be sure that PTSCOTCH attribute caching does not interfere with PETSc. */
2788be712e4SBarry Smith PetscCallMPI(MPI_Comm_dup(pcomm, &comm));
2798be712e4SBarry Smith PetscCallMPI(MPI_Comm_rank(comm, &rank));
2808be712e4SBarry Smith PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPIADJ, &flg));
2818be712e4SBarry Smith if (!flg) {
2828be712e4SBarry Smith /* bs indicates if the converted matrix is "reduced" from the original and hence the
2838be712e4SBarry Smith resulting partition results need to be stretched to match the original matrix */
2848be712e4SBarry Smith nold = mat->rmap->n;
2858be712e4SBarry Smith PetscCall(MatConvert(mat, MATMPIADJ, MAT_INITIAL_MATRIX, &mat));
2868be712e4SBarry Smith if (mat->rmap->n > 0) bs = nold / mat->rmap->n;
2878be712e4SBarry Smith adj = (Mat_MPIAdj *)mat->data;
2888be712e4SBarry Smith }
2898be712e4SBarry Smith
2908be712e4SBarry Smith proc_weight_flg = part->part_weights ? PETSC_TRUE : PETSC_FALSE;
2918be712e4SBarry Smith PetscCall(PetscOptionsGetBool(NULL, NULL, "-mat_partitioning_ptscotch_proc_weight", &proc_weight_flg, NULL));
2928be712e4SBarry Smith
2938be712e4SBarry Smith PetscCall(PetscMalloc1(mat->rmap->n + 1, &locals));
2948be712e4SBarry Smith
2958be712e4SBarry Smith if (useND) {
2968be712e4SBarry Smith #if defined(PETSC_HAVE_SCOTCH_PARMETIS_V3_NODEND)
2978be712e4SBarry Smith PetscInt *sizes, *seps, log2size, subd, *level, base = 0;
2988be712e4SBarry Smith PetscMPIInt size;
2998be712e4SBarry Smith
3008be712e4SBarry Smith PetscCallMPI(MPI_Comm_size(comm, &size));
3018be712e4SBarry Smith log2size = PetscLog2Real(size);
3028be712e4SBarry Smith subd = PetscPowInt(2, log2size);
3038be712e4SBarry Smith PetscCheck(subd == size, comm, PETSC_ERR_SUP, "Only power of 2 communicator sizes");
3048be712e4SBarry Smith PetscCall(PetscMalloc1(mat->rmap->n, &NDorder));
3058be712e4SBarry Smith PetscCall(PetscMalloc3(2 * size, &sizes, 4 * size, &seps, size, &level));
3068be712e4SBarry Smith SCOTCH_ParMETIS_V3_NodeND(mat->rmap->range, adj->i, adj->j, &base, NULL, NDorder, sizes, &comm);
3078be712e4SBarry Smith PetscCall(MatPartitioningSizesToSep_Private(subd, sizes, seps, level));
3088be712e4SBarry Smith for (i = 0; i < mat->rmap->n; i++) {
3098be712e4SBarry Smith PetscInt loc;
3108be712e4SBarry Smith
3118be712e4SBarry Smith PetscCall(PetscFindInt(NDorder[i], 2 * subd, seps, &loc));
3128be712e4SBarry Smith if (loc < 0) {
3138be712e4SBarry Smith loc = -(loc + 1);
3148be712e4SBarry Smith if (loc % 2) { /* part of subdomain */
3158be712e4SBarry Smith locals[i] = loc / 2;
3168be712e4SBarry Smith } else {
3178be712e4SBarry Smith PetscCall(PetscFindInt(NDorder[i], 2 * (subd - 1), seps + 2 * subd, &loc));
3188be712e4SBarry Smith loc = loc < 0 ? -(loc + 1) / 2 : loc / 2;
3198be712e4SBarry Smith locals[i] = level[loc];
3208be712e4SBarry Smith }
3218be712e4SBarry Smith } else locals[i] = loc / 2;
3228be712e4SBarry Smith }
3238be712e4SBarry Smith PetscCall(PetscFree3(sizes, seps, level));
3248be712e4SBarry Smith #else
3258be712e4SBarry Smith SETERRQ(pcomm, PETSC_ERR_SUP, "Need libptscotchparmetis.a compiled with -DSCOTCH_METIS_PREFIX");
3268be712e4SBarry Smith #endif
3278be712e4SBarry Smith } else {
3288be712e4SBarry Smith velotab = NULL;
3298be712e4SBarry Smith if (proc_weight_flg) {
3308be712e4SBarry Smith PetscCall(PetscMalloc1(nparts, &vwgttab));
3318be712e4SBarry Smith PetscCall(PetscMalloc1(nparts, &velotab));
3328be712e4SBarry Smith for (j = 0; j < nparts; j++) {
3338be712e4SBarry Smith if (part->part_weights) vwgttab[j] = part->part_weights[j] * nparts;
3348be712e4SBarry Smith else vwgttab[j] = 1.0;
3358be712e4SBarry Smith }
3368be712e4SBarry Smith for (i = 0; i < nparts; i++) {
3378be712e4SBarry Smith deltval = PetscAbsReal(vwgttab[i] - PetscFloorReal(vwgttab[i] + 0.5));
3388be712e4SBarry Smith if (deltval > 0.01) {
3398be712e4SBarry Smith for (j = 0; j < nparts; j++) vwgttab[j] /= deltval;
3408be712e4SBarry Smith }
3418be712e4SBarry Smith }
3428be712e4SBarry Smith for (i = 0; i < nparts; i++) velotab[i] = (SCOTCH_Num)(vwgttab[i] + 0.5);
3438be712e4SBarry Smith PetscCall(PetscFree(vwgttab));
3448be712e4SBarry Smith }
3458be712e4SBarry Smith
3468be712e4SBarry Smith vertlocnbr = mat->rmap->range[rank + 1] - mat->rmap->range[rank];
3478be712e4SBarry Smith edgelocnbr = adj->i[vertlocnbr];
3488be712e4SBarry Smith veloloctab = part->vertex_weights;
3498be712e4SBarry Smith edloloctab = part->use_edge_weights ? adj->values : NULL;
3508be712e4SBarry Smith
3518be712e4SBarry Smith /* detect whether all vertices are located at the same process in original graph */
352fbccb6d4SPierre Jolivet for (p = 0; !mat->rmap->range[p + 1] && p < nparts; ++p);
3538be712e4SBarry Smith distributed = (mat->rmap->range[p + 1] == mat->rmap->N) ? PETSC_FALSE : PETSC_TRUE;
3548be712e4SBarry Smith if (distributed) {
3558be712e4SBarry Smith SCOTCH_Arch archdat;
3568be712e4SBarry Smith SCOTCH_Dgraph grafdat;
3578be712e4SBarry Smith SCOTCH_Dmapping mappdat;
3588be712e4SBarry Smith SCOTCH_Strat stradat;
3598be712e4SBarry Smith
3608be712e4SBarry Smith PetscCallExternal(SCOTCH_dgraphInit, &grafdat, comm);
3618be712e4SBarry Smith PetscCallExternal(SCOTCH_dgraphBuild, &grafdat, 0, vertlocnbr, vertlocnbr, adj->i, adj->i + 1, veloloctab, NULL, edgelocnbr, edgelocnbr, adj->j, NULL, edloloctab);
3628be712e4SBarry Smith
3638be712e4SBarry Smith if (PetscDefined(USE_DEBUG)) PetscCallExternal(SCOTCH_dgraphCheck, &grafdat);
3648be712e4SBarry Smith
3658be712e4SBarry Smith PetscCallExternal(SCOTCH_archInit, &archdat);
3668be712e4SBarry Smith PetscCallExternal(SCOTCH_stratInit, &stradat);
3678be712e4SBarry Smith PetscCallExternal(SCOTCH_stratDgraphMapBuild, &stradat, scotch->strategy, nparts, nparts, scotch->imbalance);
3688be712e4SBarry Smith
3698be712e4SBarry Smith if (velotab) {
3708be712e4SBarry Smith PetscCallExternal(SCOTCH_archCmpltw, &archdat, nparts, velotab);
3718be712e4SBarry Smith } else {
3728be712e4SBarry Smith PetscCallExternal(SCOTCH_archCmplt, &archdat, nparts);
3738be712e4SBarry Smith }
3748be712e4SBarry Smith PetscCallExternal(SCOTCH_dgraphMapInit, &grafdat, &mappdat, &archdat, locals);
3758be712e4SBarry Smith PetscCallExternal(SCOTCH_dgraphMapCompute, &grafdat, &mappdat, &stradat);
3768be712e4SBarry Smith
3778be712e4SBarry Smith SCOTCH_dgraphMapExit(&grafdat, &mappdat);
3788be712e4SBarry Smith SCOTCH_archExit(&archdat);
3798be712e4SBarry Smith SCOTCH_stratExit(&stradat);
3808be712e4SBarry Smith SCOTCH_dgraphExit(&grafdat);
3818be712e4SBarry Smith
3828be712e4SBarry Smith } else if (rank == p) {
3838be712e4SBarry Smith SCOTCH_Graph grafdat;
3848be712e4SBarry Smith SCOTCH_Strat stradat;
3858be712e4SBarry Smith
3868be712e4SBarry Smith PetscCallExternal(SCOTCH_graphInit, &grafdat);
3878be712e4SBarry Smith PetscCallExternal(SCOTCH_graphBuild, &grafdat, 0, vertlocnbr, adj->i, adj->i + 1, veloloctab, NULL, edgelocnbr, adj->j, edloloctab);
3888be712e4SBarry Smith if (PetscDefined(USE_DEBUG)) PetscCallExternal(SCOTCH_graphCheck, &grafdat);
3898be712e4SBarry Smith PetscCallExternal(SCOTCH_stratInit, &stradat);
3908be712e4SBarry Smith PetscCallExternal(SCOTCH_stratGraphMapBuild, &stradat, scotch->strategy, nparts, scotch->imbalance);
3918be712e4SBarry Smith if (velotab) {
3928be712e4SBarry Smith SCOTCH_Arch archdat;
3938be712e4SBarry Smith PetscCallExternal(SCOTCH_archInit, &archdat);
3948be712e4SBarry Smith PetscCallExternal(SCOTCH_archCmpltw, &archdat, nparts, velotab);
3958be712e4SBarry Smith PetscCallExternal(SCOTCH_graphMap, &grafdat, &archdat, &stradat, locals);
3968be712e4SBarry Smith SCOTCH_archExit(&archdat);
3978be712e4SBarry Smith } else {
3988be712e4SBarry Smith PetscCallExternal(SCOTCH_graphPart, &grafdat, nparts, &stradat, locals);
3998be712e4SBarry Smith }
4008be712e4SBarry Smith SCOTCH_stratExit(&stradat);
4018be712e4SBarry Smith SCOTCH_graphExit(&grafdat);
4028be712e4SBarry Smith }
4038be712e4SBarry Smith
4048be712e4SBarry Smith PetscCall(PetscFree(velotab));
4058be712e4SBarry Smith }
4068be712e4SBarry Smith PetscCallMPI(MPI_Comm_free(&comm));
4078be712e4SBarry Smith
4088be712e4SBarry Smith if (bs > 1) {
4098be712e4SBarry Smith PetscInt *newlocals;
4108be712e4SBarry Smith PetscCall(PetscMalloc1(bs * mat->rmap->n, &newlocals));
4118be712e4SBarry Smith for (i = 0; i < mat->rmap->n; i++) {
4128be712e4SBarry Smith for (j = 0; j < bs; j++) newlocals[bs * i + j] = locals[i];
4138be712e4SBarry Smith }
4148be712e4SBarry Smith PetscCall(PetscFree(locals));
4158be712e4SBarry Smith PetscCall(ISCreateGeneral(pcomm, bs * mat->rmap->n, newlocals, PETSC_OWN_POINTER, partitioning));
4168be712e4SBarry Smith } else {
4178be712e4SBarry Smith PetscCall(ISCreateGeneral(pcomm, mat->rmap->n, locals, PETSC_OWN_POINTER, partitioning));
4188be712e4SBarry Smith }
4198be712e4SBarry Smith if (useND) {
4208be712e4SBarry Smith IS ndis;
4218be712e4SBarry Smith
4228be712e4SBarry Smith if (bs > 1) {
4238be712e4SBarry Smith PetscCall(ISCreateBlock(pcomm, bs, mat->rmap->n, NDorder, PETSC_OWN_POINTER, &ndis));
4248be712e4SBarry Smith } else {
4258be712e4SBarry Smith PetscCall(ISCreateGeneral(pcomm, mat->rmap->n, NDorder, PETSC_OWN_POINTER, &ndis));
4268be712e4SBarry Smith }
4278be712e4SBarry Smith PetscCall(ISSetPermutation(ndis));
428f4f49eeaSPierre Jolivet PetscCall(PetscObjectCompose((PetscObject)*partitioning, "_petsc_matpartitioning_ndorder", (PetscObject)ndis));
4298be712e4SBarry Smith PetscCall(ISDestroy(&ndis));
4308be712e4SBarry Smith }
4318be712e4SBarry Smith
4328be712e4SBarry Smith if (!flg) PetscCall(MatDestroy(&mat));
4338be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
4348be712e4SBarry Smith }
4358be712e4SBarry Smith
MatPartitioningApply_PTScotch(MatPartitioning part,IS * partitioning)4368be712e4SBarry Smith static PetscErrorCode MatPartitioningApply_PTScotch(MatPartitioning part, IS *partitioning)
4378be712e4SBarry Smith {
4388be712e4SBarry Smith PetscFunctionBegin;
4398be712e4SBarry Smith PetscCall(MatPartitioningApply_PTScotch_Private(part, PETSC_FALSE, partitioning));
4408be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
4418be712e4SBarry Smith }
4428be712e4SBarry Smith
MatPartitioningApplyND_PTScotch(MatPartitioning part,IS * partitioning)4438be712e4SBarry Smith static PetscErrorCode MatPartitioningApplyND_PTScotch(MatPartitioning part, IS *partitioning)
4448be712e4SBarry Smith {
4458be712e4SBarry Smith PetscFunctionBegin;
4468be712e4SBarry Smith PetscCall(MatPartitioningApply_PTScotch_Private(part, PETSC_TRUE, partitioning));
4478be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
4488be712e4SBarry Smith }
4498be712e4SBarry Smith
MatPartitioningDestroy_PTScotch(MatPartitioning part)4508be712e4SBarry Smith static PetscErrorCode MatPartitioningDestroy_PTScotch(MatPartitioning part)
4518be712e4SBarry Smith {
4528be712e4SBarry Smith MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch *)part->data;
4538be712e4SBarry Smith
4548be712e4SBarry Smith PetscFunctionBegin;
4558be712e4SBarry Smith PetscCall(PetscFree(scotch));
4568be712e4SBarry Smith /* clear composed functions */
4578be712e4SBarry Smith PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPTScotchSetImbalance_C", NULL));
4588be712e4SBarry Smith PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPTScotchGetImbalance_C", NULL));
4598be712e4SBarry Smith PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPTScotchSetStrategy_C", NULL));
4608be712e4SBarry Smith PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPTScotchGetStrategy_C", NULL));
4618be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
4628be712e4SBarry Smith }
4638be712e4SBarry Smith
4648be712e4SBarry Smith /*MC
465613ce9feSSatish Balay MATPARTITIONINGPTSCOTCH - Creates a partitioning context that uses the external package SCOTCH <http://www.labri.fr/perso/pelegrin/scotch/>.
4668be712e4SBarry Smith
4678be712e4SBarry Smith Level: beginner
4688be712e4SBarry Smith
4698be712e4SBarry Smith .seealso: `MatPartitioningSetType()`, `MatPartitioningType`, `MatPartitioningPTScotchSetImbalance()`, `MatPartitioningPTScotchGetImbalance()`,
470613ce9feSSatish Balay `MatPartitioningPTScotchSetStrategy()`, `MatPartitioningPTScotchGetStrategy()`, `MatPartitioning`
4718be712e4SBarry Smith M*/
4728be712e4SBarry Smith
MatPartitioningCreate_PTScotch(MatPartitioning part)4738be712e4SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningCreate_PTScotch(MatPartitioning part)
4748be712e4SBarry Smith {
4758be712e4SBarry Smith MatPartitioning_PTScotch *scotch;
4768be712e4SBarry Smith
4778be712e4SBarry Smith PetscFunctionBegin;
4788be712e4SBarry Smith PetscCall(PetscNew(&scotch));
4798be712e4SBarry Smith part->data = (void *)scotch;
4808be712e4SBarry Smith
4818be712e4SBarry Smith scotch->imbalance = 0.01;
4828be712e4SBarry Smith scotch->strategy = SCOTCH_STRATDEFAULT;
4838be712e4SBarry Smith
4848be712e4SBarry Smith part->ops->apply = MatPartitioningApply_PTScotch;
4858be712e4SBarry Smith part->ops->applynd = MatPartitioningApplyND_PTScotch;
4868be712e4SBarry Smith part->ops->view = MatPartitioningView_PTScotch;
4878be712e4SBarry Smith part->ops->setfromoptions = MatPartitioningSetFromOptions_PTScotch;
4888be712e4SBarry Smith part->ops->destroy = MatPartitioningDestroy_PTScotch;
4898be712e4SBarry Smith
4908be712e4SBarry Smith PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPTScotchSetImbalance_C", MatPartitioningPTScotchSetImbalance_PTScotch));
4918be712e4SBarry Smith PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPTScotchGetImbalance_C", MatPartitioningPTScotchGetImbalance_PTScotch));
4928be712e4SBarry Smith PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPTScotchSetStrategy_C", MatPartitioningPTScotchSetStrategy_PTScotch));
4938be712e4SBarry Smith PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPTScotchGetStrategy_C", MatPartitioningPTScotchGetStrategy_PTScotch));
4948be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
4958be712e4SBarry Smith }
496