xref: /petsc/src/mat/graphops/partition/impls/scotch/scotch.c (revision 8be712e46db5d855f641c6bd97b4543e0efe65bd)
1*8be712e4SBarry Smith #include <../src/mat/impls/adj/mpi/mpiadj.h> /*I "petscmat.h" I*/
2*8be712e4SBarry Smith 
3*8be712e4SBarry Smith EXTERN_C_BEGIN
4*8be712e4SBarry Smith #include <ptscotch.h>
5*8be712e4SBarry Smith #if defined(PETSC_HAVE_SCOTCH_PARMETIS_V3_NODEND)
6*8be712e4SBarry Smith /* we define the prototype instead of include SCOTCH's parmetis.h */
7*8be712e4SBarry 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 *);
8*8be712e4SBarry Smith #endif
9*8be712e4SBarry Smith EXTERN_C_END
10*8be712e4SBarry Smith 
11*8be712e4SBarry Smith typedef struct {
12*8be712e4SBarry Smith   double     imbalance;
13*8be712e4SBarry Smith   SCOTCH_Num strategy;
14*8be712e4SBarry Smith } MatPartitioning_PTScotch;
15*8be712e4SBarry Smith 
16*8be712e4SBarry Smith /*@
17*8be712e4SBarry Smith   MatPartitioningPTScotchSetImbalance - Sets the value of the load imbalance
18*8be712e4SBarry Smith   ratio to be used during strategy selection.
19*8be712e4SBarry Smith 
20*8be712e4SBarry Smith   Collective
21*8be712e4SBarry Smith 
22*8be712e4SBarry Smith   Input Parameters:
23*8be712e4SBarry Smith + part - the partitioning context
24*8be712e4SBarry Smith - imb  - the load imbalance ratio
25*8be712e4SBarry Smith 
26*8be712e4SBarry Smith   Options Database Key:
27*8be712e4SBarry Smith . -mat_partitioning_ptscotch_imbalance <imb> - set load imbalance ratio
28*8be712e4SBarry Smith 
29*8be712e4SBarry Smith   Note:
30*8be712e4SBarry Smith   Must be in the range [0,1]. The default value is 0.01.
31*8be712e4SBarry Smith 
32*8be712e4SBarry Smith   Level: advanced
33*8be712e4SBarry Smith 
34*8be712e4SBarry Smith .seealso: `MATPARTITIONINGSCOTCH`, `MatPartitioningPTScotchSetStrategy()`, `MatPartitioningPTScotchGetImbalance()`
35*8be712e4SBarry Smith @*/
36*8be712e4SBarry Smith PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning part, PetscReal imb)
37*8be712e4SBarry Smith {
38*8be712e4SBarry Smith   PetscFunctionBegin;
39*8be712e4SBarry Smith   PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
40*8be712e4SBarry Smith   PetscValidLogicalCollectiveReal(part, imb, 2);
41*8be712e4SBarry Smith   PetscTryMethod(part, "MatPartitioningPTScotchSetImbalance_C", (MatPartitioning, PetscReal), (part, imb));
42*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
43*8be712e4SBarry Smith }
44*8be712e4SBarry Smith 
45*8be712e4SBarry Smith static PetscErrorCode MatPartitioningPTScotchSetImbalance_PTScotch(MatPartitioning part, PetscReal imb)
46*8be712e4SBarry Smith {
47*8be712e4SBarry Smith   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch *)part->data;
48*8be712e4SBarry Smith 
49*8be712e4SBarry Smith   PetscFunctionBegin;
50*8be712e4SBarry Smith   if (imb == PETSC_DEFAULT) scotch->imbalance = 0.01;
51*8be712e4SBarry Smith   else {
52*8be712e4SBarry 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]");
53*8be712e4SBarry Smith     scotch->imbalance = (double)imb;
54*8be712e4SBarry Smith   }
55*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
56*8be712e4SBarry Smith }
57*8be712e4SBarry Smith 
58*8be712e4SBarry Smith /*@
59*8be712e4SBarry Smith   MatPartitioningPTScotchGetImbalance - Gets the value of the load imbalance
60*8be712e4SBarry Smith   ratio used during strategy selection.
61*8be712e4SBarry Smith 
62*8be712e4SBarry Smith   Not Collective
63*8be712e4SBarry Smith 
64*8be712e4SBarry Smith   Input Parameter:
65*8be712e4SBarry Smith . part - the partitioning context
66*8be712e4SBarry Smith 
67*8be712e4SBarry Smith   Output Parameter:
68*8be712e4SBarry Smith . imb - the load imbalance ratio
69*8be712e4SBarry Smith 
70*8be712e4SBarry Smith   Level: advanced
71*8be712e4SBarry Smith 
72*8be712e4SBarry Smith .seealso: `MATPARTITIONINGSCOTCH`, `MatPartitioningPTScotchSetImbalance()`
73*8be712e4SBarry Smith @*/
74*8be712e4SBarry Smith PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning part, PetscReal *imb)
75*8be712e4SBarry Smith {
76*8be712e4SBarry Smith   PetscFunctionBegin;
77*8be712e4SBarry Smith   PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
78*8be712e4SBarry Smith   PetscAssertPointer(imb, 2);
79*8be712e4SBarry Smith   PetscUseMethod(part, "MatPartitioningPTScotchGetImbalance_C", (MatPartitioning, PetscReal *), (part, imb));
80*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
81*8be712e4SBarry Smith }
82*8be712e4SBarry Smith 
83*8be712e4SBarry Smith static PetscErrorCode MatPartitioningPTScotchGetImbalance_PTScotch(MatPartitioning part, PetscReal *imb)
84*8be712e4SBarry Smith {
85*8be712e4SBarry Smith   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch *)part->data;
86*8be712e4SBarry Smith 
87*8be712e4SBarry Smith   PetscFunctionBegin;
88*8be712e4SBarry Smith   *imb = scotch->imbalance;
89*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
90*8be712e4SBarry Smith }
91*8be712e4SBarry Smith 
92*8be712e4SBarry Smith /*@
93*8be712e4SBarry Smith   MatPartitioningPTScotchSetStrategy - Sets the strategy to be used in PTScotch.
94*8be712e4SBarry Smith 
95*8be712e4SBarry Smith   Collective
96*8be712e4SBarry Smith 
97*8be712e4SBarry Smith   Input Parameters:
98*8be712e4SBarry Smith + part     - the partitioning context
99*8be712e4SBarry Smith - strategy - the strategy, one of
100*8be712e4SBarry Smith .vb
101*8be712e4SBarry Smith      MP_PTSCOTCH_DEFAULT     - Default behavior
102*8be712e4SBarry Smith      MP_PTSCOTCH_QUALITY     - Prioritize quality over speed
103*8be712e4SBarry Smith      MP_PTSCOTCH_SPEED       - Prioritize speed over quality
104*8be712e4SBarry Smith      MP_PTSCOTCH_BALANCE     - Enforce load balance
105*8be712e4SBarry Smith      MP_PTSCOTCH_SAFETY      - Avoid methods that may fail
106*8be712e4SBarry Smith      MP_PTSCOTCH_SCALABILITY - Favor scalability as much as possible
107*8be712e4SBarry Smith .ve
108*8be712e4SBarry Smith 
109*8be712e4SBarry Smith   Options Database Key:
110*8be712e4SBarry Smith . -mat_partitioning_ptscotch_strategy [quality,speed,balance,safety,scalability] - strategy
111*8be712e4SBarry Smith 
112*8be712e4SBarry Smith   Level: advanced
113*8be712e4SBarry Smith 
114*8be712e4SBarry Smith   Note:
115*8be712e4SBarry Smith   The default is `MP_SCOTCH_QUALITY`. See the PTScotch documentation for more information.
116*8be712e4SBarry Smith 
117*8be712e4SBarry Smith .seealso: `MATPARTITIONINGSCOTCH`, `MatPartitioningPTScotchSetImbalance()`, `MatPartitioningPTScotchGetStrategy()`
118*8be712e4SBarry Smith @*/
119*8be712e4SBarry Smith PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning part, MPPTScotchStrategyType strategy)
120*8be712e4SBarry Smith {
121*8be712e4SBarry Smith   PetscFunctionBegin;
122*8be712e4SBarry Smith   PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
123*8be712e4SBarry Smith   PetscValidLogicalCollectiveEnum(part, strategy, 2);
124*8be712e4SBarry Smith   PetscTryMethod(part, "MatPartitioningPTScotchSetStrategy_C", (MatPartitioning, MPPTScotchStrategyType), (part, strategy));
125*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
126*8be712e4SBarry Smith }
127*8be712e4SBarry Smith 
128*8be712e4SBarry Smith static PetscErrorCode MatPartitioningPTScotchSetStrategy_PTScotch(MatPartitioning part, MPPTScotchStrategyType strategy)
129*8be712e4SBarry Smith {
130*8be712e4SBarry Smith   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch *)part->data;
131*8be712e4SBarry Smith 
132*8be712e4SBarry Smith   PetscFunctionBegin;
133*8be712e4SBarry Smith   switch (strategy) {
134*8be712e4SBarry Smith   case MP_PTSCOTCH_QUALITY:
135*8be712e4SBarry Smith     scotch->strategy = SCOTCH_STRATQUALITY;
136*8be712e4SBarry Smith     break;
137*8be712e4SBarry Smith   case MP_PTSCOTCH_SPEED:
138*8be712e4SBarry Smith     scotch->strategy = SCOTCH_STRATSPEED;
139*8be712e4SBarry Smith     break;
140*8be712e4SBarry Smith   case MP_PTSCOTCH_BALANCE:
141*8be712e4SBarry Smith     scotch->strategy = SCOTCH_STRATBALANCE;
142*8be712e4SBarry Smith     break;
143*8be712e4SBarry Smith   case MP_PTSCOTCH_SAFETY:
144*8be712e4SBarry Smith     scotch->strategy = SCOTCH_STRATSAFETY;
145*8be712e4SBarry Smith     break;
146*8be712e4SBarry Smith   case MP_PTSCOTCH_SCALABILITY:
147*8be712e4SBarry Smith     scotch->strategy = SCOTCH_STRATSCALABILITY;
148*8be712e4SBarry Smith     break;
149*8be712e4SBarry Smith   default:
150*8be712e4SBarry Smith     scotch->strategy = SCOTCH_STRATDEFAULT;
151*8be712e4SBarry Smith     break;
152*8be712e4SBarry Smith   }
153*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
154*8be712e4SBarry Smith }
155*8be712e4SBarry Smith 
156*8be712e4SBarry Smith /*@
157*8be712e4SBarry Smith   MatPartitioningPTScotchGetStrategy - Gets the strategy used in PTScotch.
158*8be712e4SBarry Smith 
159*8be712e4SBarry Smith   Not Collective
160*8be712e4SBarry Smith 
161*8be712e4SBarry Smith   Input Parameter:
162*8be712e4SBarry Smith . part - the partitioning context
163*8be712e4SBarry Smith 
164*8be712e4SBarry Smith   Output Parameter:
165*8be712e4SBarry Smith . strategy - the strategy
166*8be712e4SBarry Smith 
167*8be712e4SBarry Smith   Level: advanced
168*8be712e4SBarry Smith 
169*8be712e4SBarry Smith .seealso: `MATPARTITIONINGSCOTCH`, `MatPartitioningPTScotchSetStrategy()`
170*8be712e4SBarry Smith @*/
171*8be712e4SBarry Smith PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning part, MPPTScotchStrategyType *strategy)
172*8be712e4SBarry Smith {
173*8be712e4SBarry Smith   PetscFunctionBegin;
174*8be712e4SBarry Smith   PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
175*8be712e4SBarry Smith   PetscAssertPointer(strategy, 2);
176*8be712e4SBarry Smith   PetscUseMethod(part, "MatPartitioningPTScotchGetStrategy_C", (MatPartitioning, MPPTScotchStrategyType *), (part, strategy));
177*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
178*8be712e4SBarry Smith }
179*8be712e4SBarry Smith 
180*8be712e4SBarry Smith static PetscErrorCode MatPartitioningPTScotchGetStrategy_PTScotch(MatPartitioning part, MPPTScotchStrategyType *strategy)
181*8be712e4SBarry Smith {
182*8be712e4SBarry Smith   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch *)part->data;
183*8be712e4SBarry Smith 
184*8be712e4SBarry Smith   PetscFunctionBegin;
185*8be712e4SBarry Smith   switch (scotch->strategy) {
186*8be712e4SBarry Smith   case SCOTCH_STRATQUALITY:
187*8be712e4SBarry Smith     *strategy = MP_PTSCOTCH_QUALITY;
188*8be712e4SBarry Smith     break;
189*8be712e4SBarry Smith   case SCOTCH_STRATSPEED:
190*8be712e4SBarry Smith     *strategy = MP_PTSCOTCH_SPEED;
191*8be712e4SBarry Smith     break;
192*8be712e4SBarry Smith   case SCOTCH_STRATBALANCE:
193*8be712e4SBarry Smith     *strategy = MP_PTSCOTCH_BALANCE;
194*8be712e4SBarry Smith     break;
195*8be712e4SBarry Smith   case SCOTCH_STRATSAFETY:
196*8be712e4SBarry Smith     *strategy = MP_PTSCOTCH_SAFETY;
197*8be712e4SBarry Smith     break;
198*8be712e4SBarry Smith   case SCOTCH_STRATSCALABILITY:
199*8be712e4SBarry Smith     *strategy = MP_PTSCOTCH_SCALABILITY;
200*8be712e4SBarry Smith     break;
201*8be712e4SBarry Smith   default:
202*8be712e4SBarry Smith     *strategy = MP_PTSCOTCH_DEFAULT;
203*8be712e4SBarry Smith     break;
204*8be712e4SBarry Smith   }
205*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
206*8be712e4SBarry Smith }
207*8be712e4SBarry Smith 
208*8be712e4SBarry Smith static PetscErrorCode MatPartitioningView_PTScotch(MatPartitioning part, PetscViewer viewer)
209*8be712e4SBarry Smith {
210*8be712e4SBarry Smith   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch *)part->data;
211*8be712e4SBarry Smith   PetscBool                 isascii;
212*8be712e4SBarry Smith   const char               *str = NULL;
213*8be712e4SBarry Smith 
214*8be712e4SBarry Smith   PetscFunctionBegin;
215*8be712e4SBarry Smith   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
216*8be712e4SBarry Smith   if (isascii) {
217*8be712e4SBarry Smith     switch (scotch->strategy) {
218*8be712e4SBarry Smith     case SCOTCH_STRATQUALITY:
219*8be712e4SBarry Smith       str = "Prioritize quality over speed";
220*8be712e4SBarry Smith       break;
221*8be712e4SBarry Smith     case SCOTCH_STRATSPEED:
222*8be712e4SBarry Smith       str = "Prioritize speed over quality";
223*8be712e4SBarry Smith       break;
224*8be712e4SBarry Smith     case SCOTCH_STRATBALANCE:
225*8be712e4SBarry Smith       str = "Enforce load balance";
226*8be712e4SBarry Smith       break;
227*8be712e4SBarry Smith     case SCOTCH_STRATSAFETY:
228*8be712e4SBarry Smith       str = "Avoid methods that may fail";
229*8be712e4SBarry Smith       break;
230*8be712e4SBarry Smith     case SCOTCH_STRATSCALABILITY:
231*8be712e4SBarry Smith       str = "Favor scalability as much as possible";
232*8be712e4SBarry Smith       break;
233*8be712e4SBarry Smith     default:
234*8be712e4SBarry Smith       str = "Default behavior";
235*8be712e4SBarry Smith       break;
236*8be712e4SBarry Smith     }
237*8be712e4SBarry Smith     PetscCall(PetscViewerASCIIPrintf(viewer, "  Strategy=%s\n", str));
238*8be712e4SBarry Smith     PetscCall(PetscViewerASCIIPrintf(viewer, "  Load imbalance ratio=%g\n", scotch->imbalance));
239*8be712e4SBarry Smith   }
240*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
241*8be712e4SBarry Smith }
242*8be712e4SBarry Smith 
243*8be712e4SBarry Smith static PetscErrorCode MatPartitioningSetFromOptions_PTScotch(MatPartitioning part, PetscOptionItems *PetscOptionsObject)
244*8be712e4SBarry Smith {
245*8be712e4SBarry Smith   PetscBool                 flag;
246*8be712e4SBarry Smith   PetscReal                 r;
247*8be712e4SBarry Smith   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch *)part->data;
248*8be712e4SBarry Smith   MPPTScotchStrategyType    strat;
249*8be712e4SBarry Smith 
250*8be712e4SBarry Smith   PetscFunctionBegin;
251*8be712e4SBarry Smith   PetscCall(MatPartitioningPTScotchGetStrategy(part, &strat));
252*8be712e4SBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "PTScotch partitioning options");
253*8be712e4SBarry Smith   PetscCall(PetscOptionsEnum("-mat_partitioning_ptscotch_strategy", "Strategy", "MatPartitioningPTScotchSetStrategy", MPPTScotchStrategyTypes, (PetscEnum)strat, (PetscEnum *)&strat, &flag));
254*8be712e4SBarry Smith   if (flag) PetscCall(MatPartitioningPTScotchSetStrategy(part, strat));
255*8be712e4SBarry Smith   PetscCall(PetscOptionsReal("-mat_partitioning_ptscotch_imbalance", "Load imbalance ratio", "MatPartitioningPTScotchSetImbalance", scotch->imbalance, &r, &flag));
256*8be712e4SBarry Smith   if (flag) PetscCall(MatPartitioningPTScotchSetImbalance(part, r));
257*8be712e4SBarry Smith   PetscOptionsHeadEnd();
258*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
259*8be712e4SBarry Smith }
260*8be712e4SBarry Smith 
261*8be712e4SBarry Smith static PetscErrorCode MatPartitioningApply_PTScotch_Private(MatPartitioning part, PetscBool useND, IS *partitioning)
262*8be712e4SBarry Smith {
263*8be712e4SBarry Smith   MPI_Comm                  pcomm, comm;
264*8be712e4SBarry Smith   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch *)part->data;
265*8be712e4SBarry Smith   PetscMPIInt               rank;
266*8be712e4SBarry Smith   Mat                       mat = part->adj;
267*8be712e4SBarry Smith   Mat_MPIAdj               *adj = (Mat_MPIAdj *)mat->data;
268*8be712e4SBarry Smith   PetscBool                 flg, distributed;
269*8be712e4SBarry Smith   PetscBool                 proc_weight_flg;
270*8be712e4SBarry Smith   PetscInt                  i, j, p, bs = 1, nold;
271*8be712e4SBarry Smith   PetscInt                 *NDorder = NULL;
272*8be712e4SBarry Smith   PetscReal                *vwgttab, deltval;
273*8be712e4SBarry Smith   SCOTCH_Num               *locals, *velotab, *veloloctab, *edloloctab, vertlocnbr, edgelocnbr, nparts = part->n;
274*8be712e4SBarry Smith 
275*8be712e4SBarry Smith   PetscFunctionBegin;
276*8be712e4SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)part, &pcomm));
277*8be712e4SBarry Smith   /* Duplicate the communicator to be sure that PTSCOTCH attribute caching does not interfere with PETSc. */
278*8be712e4SBarry Smith   PetscCallMPI(MPI_Comm_dup(pcomm, &comm));
279*8be712e4SBarry Smith   PetscCallMPI(MPI_Comm_rank(comm, &rank));
280*8be712e4SBarry Smith   PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPIADJ, &flg));
281*8be712e4SBarry Smith   if (!flg) {
282*8be712e4SBarry Smith     /* bs indicates if the converted matrix is "reduced" from the original and hence the
283*8be712e4SBarry Smith        resulting partition results need to be stretched to match the original matrix */
284*8be712e4SBarry Smith     nold = mat->rmap->n;
285*8be712e4SBarry Smith     PetscCall(MatConvert(mat, MATMPIADJ, MAT_INITIAL_MATRIX, &mat));
286*8be712e4SBarry Smith     if (mat->rmap->n > 0) bs = nold / mat->rmap->n;
287*8be712e4SBarry Smith     adj = (Mat_MPIAdj *)mat->data;
288*8be712e4SBarry Smith   }
289*8be712e4SBarry Smith 
290*8be712e4SBarry Smith   proc_weight_flg = part->part_weights ? PETSC_TRUE : PETSC_FALSE;
291*8be712e4SBarry Smith   PetscCall(PetscOptionsGetBool(NULL, NULL, "-mat_partitioning_ptscotch_proc_weight", &proc_weight_flg, NULL));
292*8be712e4SBarry Smith 
293*8be712e4SBarry Smith   PetscCall(PetscMalloc1(mat->rmap->n + 1, &locals));
294*8be712e4SBarry Smith 
295*8be712e4SBarry Smith   if (useND) {
296*8be712e4SBarry Smith #if defined(PETSC_HAVE_SCOTCH_PARMETIS_V3_NODEND)
297*8be712e4SBarry Smith     PetscInt   *sizes, *seps, log2size, subd, *level, base = 0;
298*8be712e4SBarry Smith     PetscMPIInt size;
299*8be712e4SBarry Smith 
300*8be712e4SBarry Smith     PetscCallMPI(MPI_Comm_size(comm, &size));
301*8be712e4SBarry Smith     log2size = PetscLog2Real(size);
302*8be712e4SBarry Smith     subd     = PetscPowInt(2, log2size);
303*8be712e4SBarry Smith     PetscCheck(subd == size, comm, PETSC_ERR_SUP, "Only power of 2 communicator sizes");
304*8be712e4SBarry Smith     PetscCall(PetscMalloc1(mat->rmap->n, &NDorder));
305*8be712e4SBarry Smith     PetscCall(PetscMalloc3(2 * size, &sizes, 4 * size, &seps, size, &level));
306*8be712e4SBarry Smith     SCOTCH_ParMETIS_V3_NodeND(mat->rmap->range, adj->i, adj->j, &base, NULL, NDorder, sizes, &comm);
307*8be712e4SBarry Smith     PetscCall(MatPartitioningSizesToSep_Private(subd, sizes, seps, level));
308*8be712e4SBarry Smith     for (i = 0; i < mat->rmap->n; i++) {
309*8be712e4SBarry Smith       PetscInt loc;
310*8be712e4SBarry Smith 
311*8be712e4SBarry Smith       PetscCall(PetscFindInt(NDorder[i], 2 * subd, seps, &loc));
312*8be712e4SBarry Smith       if (loc < 0) {
313*8be712e4SBarry Smith         loc = -(loc + 1);
314*8be712e4SBarry Smith         if (loc % 2) { /* part of subdomain */
315*8be712e4SBarry Smith           locals[i] = loc / 2;
316*8be712e4SBarry Smith         } else {
317*8be712e4SBarry Smith           PetscCall(PetscFindInt(NDorder[i], 2 * (subd - 1), seps + 2 * subd, &loc));
318*8be712e4SBarry Smith           loc       = loc < 0 ? -(loc + 1) / 2 : loc / 2;
319*8be712e4SBarry Smith           locals[i] = level[loc];
320*8be712e4SBarry Smith         }
321*8be712e4SBarry Smith       } else locals[i] = loc / 2;
322*8be712e4SBarry Smith     }
323*8be712e4SBarry Smith     PetscCall(PetscFree3(sizes, seps, level));
324*8be712e4SBarry Smith #else
325*8be712e4SBarry Smith     SETERRQ(pcomm, PETSC_ERR_SUP, "Need libptscotchparmetis.a compiled with -DSCOTCH_METIS_PREFIX");
326*8be712e4SBarry Smith #endif
327*8be712e4SBarry Smith   } else {
328*8be712e4SBarry Smith     velotab = NULL;
329*8be712e4SBarry Smith     if (proc_weight_flg) {
330*8be712e4SBarry Smith       PetscCall(PetscMalloc1(nparts, &vwgttab));
331*8be712e4SBarry Smith       PetscCall(PetscMalloc1(nparts, &velotab));
332*8be712e4SBarry Smith       for (j = 0; j < nparts; j++) {
333*8be712e4SBarry Smith         if (part->part_weights) vwgttab[j] = part->part_weights[j] * nparts;
334*8be712e4SBarry Smith         else vwgttab[j] = 1.0;
335*8be712e4SBarry Smith       }
336*8be712e4SBarry Smith       for (i = 0; i < nparts; i++) {
337*8be712e4SBarry Smith         deltval = PetscAbsReal(vwgttab[i] - PetscFloorReal(vwgttab[i] + 0.5));
338*8be712e4SBarry Smith         if (deltval > 0.01) {
339*8be712e4SBarry Smith           for (j = 0; j < nparts; j++) vwgttab[j] /= deltval;
340*8be712e4SBarry Smith         }
341*8be712e4SBarry Smith       }
342*8be712e4SBarry Smith       for (i = 0; i < nparts; i++) velotab[i] = (SCOTCH_Num)(vwgttab[i] + 0.5);
343*8be712e4SBarry Smith       PetscCall(PetscFree(vwgttab));
344*8be712e4SBarry Smith     }
345*8be712e4SBarry Smith 
346*8be712e4SBarry Smith     vertlocnbr = mat->rmap->range[rank + 1] - mat->rmap->range[rank];
347*8be712e4SBarry Smith     edgelocnbr = adj->i[vertlocnbr];
348*8be712e4SBarry Smith     veloloctab = part->vertex_weights;
349*8be712e4SBarry Smith     edloloctab = part->use_edge_weights ? adj->values : NULL;
350*8be712e4SBarry Smith 
351*8be712e4SBarry Smith     /* detect whether all vertices are located at the same process in original graph */
352*8be712e4SBarry Smith     for (p = 0; !mat->rmap->range[p + 1] && p < nparts; ++p)
353*8be712e4SBarry Smith       ;
354*8be712e4SBarry Smith     distributed = (mat->rmap->range[p + 1] == mat->rmap->N) ? PETSC_FALSE : PETSC_TRUE;
355*8be712e4SBarry Smith     if (distributed) {
356*8be712e4SBarry Smith       SCOTCH_Arch     archdat;
357*8be712e4SBarry Smith       SCOTCH_Dgraph   grafdat;
358*8be712e4SBarry Smith       SCOTCH_Dmapping mappdat;
359*8be712e4SBarry Smith       SCOTCH_Strat    stradat;
360*8be712e4SBarry Smith 
361*8be712e4SBarry Smith       PetscCallExternal(SCOTCH_dgraphInit, &grafdat, comm);
362*8be712e4SBarry Smith       PetscCallExternal(SCOTCH_dgraphBuild, &grafdat, 0, vertlocnbr, vertlocnbr, adj->i, adj->i + 1, veloloctab, NULL, edgelocnbr, edgelocnbr, adj->j, NULL, edloloctab);
363*8be712e4SBarry Smith 
364*8be712e4SBarry Smith       if (PetscDefined(USE_DEBUG)) PetscCallExternal(SCOTCH_dgraphCheck, &grafdat);
365*8be712e4SBarry Smith 
366*8be712e4SBarry Smith       PetscCallExternal(SCOTCH_archInit, &archdat);
367*8be712e4SBarry Smith       PetscCallExternal(SCOTCH_stratInit, &stradat);
368*8be712e4SBarry Smith       PetscCallExternal(SCOTCH_stratDgraphMapBuild, &stradat, scotch->strategy, nparts, nparts, scotch->imbalance);
369*8be712e4SBarry Smith 
370*8be712e4SBarry Smith       if (velotab) {
371*8be712e4SBarry Smith         PetscCallExternal(SCOTCH_archCmpltw, &archdat, nparts, velotab);
372*8be712e4SBarry Smith       } else {
373*8be712e4SBarry Smith         PetscCallExternal(SCOTCH_archCmplt, &archdat, nparts);
374*8be712e4SBarry Smith       }
375*8be712e4SBarry Smith       PetscCallExternal(SCOTCH_dgraphMapInit, &grafdat, &mappdat, &archdat, locals);
376*8be712e4SBarry Smith       PetscCallExternal(SCOTCH_dgraphMapCompute, &grafdat, &mappdat, &stradat);
377*8be712e4SBarry Smith 
378*8be712e4SBarry Smith       SCOTCH_dgraphMapExit(&grafdat, &mappdat);
379*8be712e4SBarry Smith       SCOTCH_archExit(&archdat);
380*8be712e4SBarry Smith       SCOTCH_stratExit(&stradat);
381*8be712e4SBarry Smith       SCOTCH_dgraphExit(&grafdat);
382*8be712e4SBarry Smith 
383*8be712e4SBarry Smith     } else if (rank == p) {
384*8be712e4SBarry Smith       SCOTCH_Graph grafdat;
385*8be712e4SBarry Smith       SCOTCH_Strat stradat;
386*8be712e4SBarry Smith 
387*8be712e4SBarry Smith       PetscCallExternal(SCOTCH_graphInit, &grafdat);
388*8be712e4SBarry Smith       PetscCallExternal(SCOTCH_graphBuild, &grafdat, 0, vertlocnbr, adj->i, adj->i + 1, veloloctab, NULL, edgelocnbr, adj->j, edloloctab);
389*8be712e4SBarry Smith       if (PetscDefined(USE_DEBUG)) PetscCallExternal(SCOTCH_graphCheck, &grafdat);
390*8be712e4SBarry Smith       PetscCallExternal(SCOTCH_stratInit, &stradat);
391*8be712e4SBarry Smith       PetscCallExternal(SCOTCH_stratGraphMapBuild, &stradat, scotch->strategy, nparts, scotch->imbalance);
392*8be712e4SBarry Smith       if (velotab) {
393*8be712e4SBarry Smith         SCOTCH_Arch archdat;
394*8be712e4SBarry Smith         PetscCallExternal(SCOTCH_archInit, &archdat);
395*8be712e4SBarry Smith         PetscCallExternal(SCOTCH_archCmpltw, &archdat, nparts, velotab);
396*8be712e4SBarry Smith         PetscCallExternal(SCOTCH_graphMap, &grafdat, &archdat, &stradat, locals);
397*8be712e4SBarry Smith         SCOTCH_archExit(&archdat);
398*8be712e4SBarry Smith       } else {
399*8be712e4SBarry Smith         PetscCallExternal(SCOTCH_graphPart, &grafdat, nparts, &stradat, locals);
400*8be712e4SBarry Smith       }
401*8be712e4SBarry Smith       SCOTCH_stratExit(&stradat);
402*8be712e4SBarry Smith       SCOTCH_graphExit(&grafdat);
403*8be712e4SBarry Smith     }
404*8be712e4SBarry Smith 
405*8be712e4SBarry Smith     PetscCall(PetscFree(velotab));
406*8be712e4SBarry Smith   }
407*8be712e4SBarry Smith   PetscCallMPI(MPI_Comm_free(&comm));
408*8be712e4SBarry Smith 
409*8be712e4SBarry Smith   if (bs > 1) {
410*8be712e4SBarry Smith     PetscInt *newlocals;
411*8be712e4SBarry Smith     PetscCall(PetscMalloc1(bs * mat->rmap->n, &newlocals));
412*8be712e4SBarry Smith     for (i = 0; i < mat->rmap->n; i++) {
413*8be712e4SBarry Smith       for (j = 0; j < bs; j++) newlocals[bs * i + j] = locals[i];
414*8be712e4SBarry Smith     }
415*8be712e4SBarry Smith     PetscCall(PetscFree(locals));
416*8be712e4SBarry Smith     PetscCall(ISCreateGeneral(pcomm, bs * mat->rmap->n, newlocals, PETSC_OWN_POINTER, partitioning));
417*8be712e4SBarry Smith   } else {
418*8be712e4SBarry Smith     PetscCall(ISCreateGeneral(pcomm, mat->rmap->n, locals, PETSC_OWN_POINTER, partitioning));
419*8be712e4SBarry Smith   }
420*8be712e4SBarry Smith   if (useND) {
421*8be712e4SBarry Smith     IS ndis;
422*8be712e4SBarry Smith 
423*8be712e4SBarry Smith     if (bs > 1) {
424*8be712e4SBarry Smith       PetscCall(ISCreateBlock(pcomm, bs, mat->rmap->n, NDorder, PETSC_OWN_POINTER, &ndis));
425*8be712e4SBarry Smith     } else {
426*8be712e4SBarry Smith       PetscCall(ISCreateGeneral(pcomm, mat->rmap->n, NDorder, PETSC_OWN_POINTER, &ndis));
427*8be712e4SBarry Smith     }
428*8be712e4SBarry Smith     PetscCall(ISSetPermutation(ndis));
429*8be712e4SBarry Smith     PetscCall(PetscObjectCompose((PetscObject)(*partitioning), "_petsc_matpartitioning_ndorder", (PetscObject)ndis));
430*8be712e4SBarry Smith     PetscCall(ISDestroy(&ndis));
431*8be712e4SBarry Smith   }
432*8be712e4SBarry Smith 
433*8be712e4SBarry Smith   if (!flg) PetscCall(MatDestroy(&mat));
434*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
435*8be712e4SBarry Smith }
436*8be712e4SBarry Smith 
437*8be712e4SBarry Smith static PetscErrorCode MatPartitioningApply_PTScotch(MatPartitioning part, IS *partitioning)
438*8be712e4SBarry Smith {
439*8be712e4SBarry Smith   PetscFunctionBegin;
440*8be712e4SBarry Smith   PetscCall(MatPartitioningApply_PTScotch_Private(part, PETSC_FALSE, partitioning));
441*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
442*8be712e4SBarry Smith }
443*8be712e4SBarry Smith 
444*8be712e4SBarry Smith static PetscErrorCode MatPartitioningApplyND_PTScotch(MatPartitioning part, IS *partitioning)
445*8be712e4SBarry Smith {
446*8be712e4SBarry Smith   PetscFunctionBegin;
447*8be712e4SBarry Smith   PetscCall(MatPartitioningApply_PTScotch_Private(part, PETSC_TRUE, partitioning));
448*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
449*8be712e4SBarry Smith }
450*8be712e4SBarry Smith 
451*8be712e4SBarry Smith static PetscErrorCode MatPartitioningDestroy_PTScotch(MatPartitioning part)
452*8be712e4SBarry Smith {
453*8be712e4SBarry Smith   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch *)part->data;
454*8be712e4SBarry Smith 
455*8be712e4SBarry Smith   PetscFunctionBegin;
456*8be712e4SBarry Smith   PetscCall(PetscFree(scotch));
457*8be712e4SBarry Smith   /* clear composed functions */
458*8be712e4SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPTScotchSetImbalance_C", NULL));
459*8be712e4SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPTScotchGetImbalance_C", NULL));
460*8be712e4SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPTScotchSetStrategy_C", NULL));
461*8be712e4SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPTScotchGetStrategy_C", NULL));
462*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
463*8be712e4SBarry Smith }
464*8be712e4SBarry Smith 
465*8be712e4SBarry Smith /*MC
466*8be712e4SBarry Smith    MATPARTITIONINGPTSCOTCH - Creates a partitioning context that uses the external package SCOTCH.
467*8be712e4SBarry Smith 
468*8be712e4SBarry Smith    Level: beginner
469*8be712e4SBarry Smith 
470*8be712e4SBarry Smith    Note:
471*8be712e4SBarry Smith     See http://www.labri.fr/perso/pelegrin/scotch/
472*8be712e4SBarry Smith 
473*8be712e4SBarry Smith .seealso: `MatPartitioningSetType()`, `MatPartitioningType`, `MatPartitioningPTScotchSetImbalance()`, `MatPartitioningPTScotchGetImbalance()`,
474*8be712e4SBarry Smith           `MatPartitioningPTScotchSetStrategy()`, `MatPartitioningPTScotchGetStrategy()`
475*8be712e4SBarry Smith M*/
476*8be712e4SBarry Smith 
477*8be712e4SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningCreate_PTScotch(MatPartitioning part)
478*8be712e4SBarry Smith {
479*8be712e4SBarry Smith   MatPartitioning_PTScotch *scotch;
480*8be712e4SBarry Smith 
481*8be712e4SBarry Smith   PetscFunctionBegin;
482*8be712e4SBarry Smith   PetscCall(PetscNew(&scotch));
483*8be712e4SBarry Smith   part->data = (void *)scotch;
484*8be712e4SBarry Smith 
485*8be712e4SBarry Smith   scotch->imbalance = 0.01;
486*8be712e4SBarry Smith   scotch->strategy  = SCOTCH_STRATDEFAULT;
487*8be712e4SBarry Smith 
488*8be712e4SBarry Smith   part->ops->apply          = MatPartitioningApply_PTScotch;
489*8be712e4SBarry Smith   part->ops->applynd        = MatPartitioningApplyND_PTScotch;
490*8be712e4SBarry Smith   part->ops->view           = MatPartitioningView_PTScotch;
491*8be712e4SBarry Smith   part->ops->setfromoptions = MatPartitioningSetFromOptions_PTScotch;
492*8be712e4SBarry Smith   part->ops->destroy        = MatPartitioningDestroy_PTScotch;
493*8be712e4SBarry Smith 
494*8be712e4SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPTScotchSetImbalance_C", MatPartitioningPTScotchSetImbalance_PTScotch));
495*8be712e4SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPTScotchGetImbalance_C", MatPartitioningPTScotchGetImbalance_PTScotch));
496*8be712e4SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPTScotchSetStrategy_C", MatPartitioningPTScotchSetStrategy_PTScotch));
497*8be712e4SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPTScotchGetStrategy_C", MatPartitioningPTScotchGetStrategy_PTScotch));
498*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
499*8be712e4SBarry Smith }
500