45 #ifndef _ZOLTAN2_ALGPARMETIS_HPP_
46 #define _ZOLTAN2_ALGPARMETIS_HPP_
58 #ifndef HAVE_ZOLTAN2_PARMETIS
64 template <
typename Adapter>
69 const RCP<
const Comm<int> > &,
73 throw std::runtime_error(
74 "BUILD ERROR: ParMETIS requested but not compiled into Zoltan2.\n"
75 "Please set CMake flag Zoltan2_ENABLE_ParMETIS:BOOL=ON.");
85 #ifdef HAVE_ZOLTAN2_PARMETIS
87 #ifndef HAVE_ZOLTAN2_MPI
91 #error "TPL ParMETIS requires compilation with MPI. Configure with -DTPL_ENABLE_MPI:BOOL=ON or -DZoltan2_ENABLE_ParMETIS:BOOL=OFF"
99 #if (PARMETIS_MAJOR_VERSION < 4)
103 #error "Specified version of ParMETIS is not compatible with Zoltan2; upgrade to ParMETIS v4 or later, or build Zoltan2 without ParMETIS."
111 template <
typename Adapter>
112 class AlgParMETIS :
public Algorithm<Adapter>
116 typedef GraphModel<typename Adapter::base_adapter_t> graphModel_t;
117 typedef typename Adapter::lno_t
lno_t;
118 typedef typename Adapter::gno_t
gno_t;
119 typedef typename Adapter::offset_t offset_t;
120 typedef typename Adapter::scalar_t
scalar_t;
123 typedef idx_t pm_idx_t;
124 typedef real_t pm_real_t;
137 const RCP<
const Comm<int> > &problemComm__,
138 const RCP<graphModel_t> &model__) :
139 env(env__), problemComm(problemComm__),
143 void partition(
const RCP<PartitioningSolution<Adapter> > &solution);
147 const RCP<const Environment> env;
148 const RCP<const Comm<int> > problemComm;
149 const RCP<GraphModel<typename Adapter::base_adapter_t> > model;
151 void scale_weights(
size_t n, ArrayView<StridedData<lno_t, scalar_t> > &fwgts,
157 template <
typename Adapter>
159 const RCP<PartitioningSolution<Adapter> > &solution
164 size_t numGlobalParts = solution->getTargetGlobalNumberOfParts();
166 int me = problemComm->getRank();
167 int np = problemComm->getSize();
170 ArrayView<const gno_t> vtxgnos;
171 ArrayView<StridedData<lno_t, scalar_t> > vwgts;
172 int nVwgt = model->getNumWeightsPerVertex();
173 size_t nVtx = model->getVertexList(vtxgnos, vwgts);
177 pm_idx_t *pm_vwgts = NULL;
179 pm_vwgts =
new pm_idx_t[nVtx*nVwgt];
180 scale_weights(nVtx, vwgts, pm_vwgts);
184 ArrayView<const gno_t> adjgnos;
185 ArrayView<const offset_t> offsets;
186 ArrayView<StridedData<lno_t, scalar_t> > ewgts;
187 int nEwgt = model->getNumWeightsPerEdge();
188 size_t nEdge = model->getEdgeList(adjgnos, offsets, ewgts);
190 pm_idx_t *pm_ewgts = NULL;
192 pm_ewgts =
new pm_idx_t[nEdge*nEwgt];
193 scale_weights(nEdge, ewgts, pm_ewgts);
197 pm_idx_t *pm_offsets;
200 pm_idx_t pm_dummy_adj;
204 pm_adjs = &pm_dummy_adj;
208 pm_idx_t *pm_vtxdist;
209 ArrayView<size_t> vtxdist;
210 model->getVertexDist(vtxdist);
215 RCP<Comm<int> > subcomm;
220 Array<int> keepRanks(np);
221 for (
int i = 0; i < np; i++) {
222 if ((pm_vtxdist[i+1] - pm_vtxdist[i]) > 0) {
223 keepRanks[nKeep] = i;
224 pm_vtxdist[nKeep] = pm_vtxdist[i];
228 pm_vtxdist[nKeep] = pm_vtxdist[np];
230 subcomm = problemComm->createSubcommunicator(keepRanks.view(0,nKeep));
231 if (subcomm != Teuchos::null)
232 mpicomm = Teuchos::getRawMpiComm(*subcomm);
234 mpicomm = MPI_COMM_NULL;
237 mpicomm = Teuchos::getRawMpiComm(*problemComm);
241 mpicomm = Teuchos::getRawMpiComm(*problemComm);
245 pm_idx_t *pm_partList = NULL;
246 if (nVtx) pm_partList =
new pm_idx_t[nVtx];
247 for (
size_t i = 0; i < nVtx; i++) pm_partList[i] = 0;
248 int pm_return = METIS_OK;
250 if (mpicomm != MPI_COMM_NULL) {
254 pm_idx_t pm_nCon = (nVwgt == 0 ? 1 : pm_idx_t(nVwgt));
255 pm_real_t *pm_partsizes =
new pm_real_t[numGlobalParts*pm_nCon];
256 for (pm_idx_t dim = 0; dim < pm_nCon; dim++) {
257 if (!solution->criteriaHasUniformPartSizes(dim))
258 for (
size_t i=0; i<numGlobalParts; i++)
259 pm_partsizes[i*pm_nCon+dim] =
260 pm_real_t(solution->getCriteriaPartSize(dim,i));
262 for (
size_t i=0; i<numGlobalParts; i++)
263 pm_partsizes[i*pm_nCon+dim] = pm_real_t(1.)/pm_real_t(numGlobalParts);
267 double tolerance = 1.1;
268 const Teuchos::ParameterList &pl = env->getParameters();
269 const Teuchos::ParameterEntry *pe = pl.getEntryPtr(
"imbalance_tolerance");
270 if (pe) tolerance = pe->getValue<
double>(&tolerance);
274 if (tolerance <= 1.0) {
276 std::cerr <<
"Warning: ParMETIS requires imbalance_tolerance > 1.0; "
277 <<
"to comply, Zoltan2 reset imbalance_tolerance to 1.01."
282 pm_real_t *pm_imbTols =
new pm_real_t[pm_nCon];
283 for (pm_idx_t dim = 0; dim < pm_nCon; dim++)
284 pm_imbTols[dim] = pm_real_t(tolerance);
286 std::string parmetis_method(
"PARTKWAY");
287 pe = pl.getEntryPtr(
"partitioning_approach");
289 std::string approach;
290 approach = pe->getValue<std::string>(&approach);
291 if ((approach ==
"repartition") || (approach ==
"maximize_overlap")) {
294 parmetis_method =
"ADAPTIVE_REPART";
299 parmetis_method =
"PARTKWAY";
304 pm_idx_t pm_wgtflag = 2*(nVwgt > 0) + (nEwgt > 0);
305 pm_idx_t pm_numflag = 0;
306 pm_idx_t pm_edgecut = -1;
307 pm_idx_t pm_options[METIS_NOPTIONS];
309 for (
int i = 1; i < METIS_NOPTIONS; i++)
316 if (parmetis_method ==
"PARTKWAY") {
317 pm_return = ParMETIS_V3_PartKway(pm_vtxdist, pm_offsets, pm_adjs,
318 pm_vwgts, pm_ewgts, &pm_wgtflag,
319 &pm_numflag, &pm_nCon, &pm_nPart,
320 pm_partsizes, pm_imbTols, pm_options,
321 &pm_edgecut, pm_partList, &mpicomm);
323 else if (parmetis_method ==
"ADAPTIVE_REPART") {
327 pm_idx_t *pm_vsize =
new pm_idx_t[nVtx];
328 for (
size_t i = 0; i < nVtx; i++) pm_vsize[i] = 1;
330 pm_real_t itr = 100.;
331 pm_return = ParMETIS_V3_AdaptiveRepart(pm_vtxdist, pm_offsets, pm_adjs,
333 pm_vsize, pm_ewgts, &pm_wgtflag,
334 &pm_numflag, &pm_nCon, &pm_nPart,
335 pm_partsizes, pm_imbTols,
336 &itr, pm_options, &pm_edgecut,
337 pm_partList, &mpicomm);
350 throw std::logic_error(
"\nInvalid ParMETIS method requested.\n");
354 delete [] pm_partsizes;
355 delete [] pm_imbTols;
360 ArrayRCP<part_t> partList;
365 solution->setParts(partList);
367 env->memory(
"Zoltan2-ParMETIS: After creating solution");
375 if (nVwgt)
delete [] pm_vwgts;
376 if (nEwgt)
delete [] pm_ewgts;
378 if (pm_return != METIS_OK) {
379 throw std::runtime_error(
380 "\nParMETIS returned an error; no valid partition generated.\n"
381 "Look for 'PARMETIS ERROR' in your output for more details.\n");
395 template <
typename Adapter>
396 void AlgParMETIS<Adapter>::scale_weights(
398 ArrayView<StridedData<
typename Adapter::lno_t,
399 typename Adapter::scalar_t> > &fwgts,
403 const double INT_EPSILON = 1e-5;
404 const int nWgt = fwgts.size();
406 int *nonint_local =
new int[nWgt+nWgt];
407 int *nonint = nonint_local + nWgt;
409 double *sum_wgt_local =
new double[nWgt*4];
410 double *max_wgt_local = sum_wgt_local + nWgt;
411 double *sum_wgt = max_wgt_local + nWgt;
412 double *max_wgt = sum_wgt + nWgt;
414 for (
int i = 0; i < nWgt; i++) {
416 sum_wgt_local[i] = 0.;
417 max_wgt_local[i] = 0;
422 for (
int j = 0; j < nWgt; j++) {
423 for (
size_t i = 0; i < n; i++) {
424 double fw = double(fwgts[j][i]);
425 if (!nonint_local[j]) {
426 pm_idx_t tmp = (pm_idx_t) floor(fw + .5);
427 if (fabs((
double)tmp-fw) > INT_EPSILON) {
431 sum_wgt_local[j] += fw;
432 if (fw > max_wgt_local[j]) max_wgt_local[j] = fw;
436 Teuchos::reduceAll<int,int>(*problemComm, Teuchos::REDUCE_MAX, nWgt,
437 nonint_local, nonint);
438 Teuchos::reduceAll<int,double>(*problemComm, Teuchos::REDUCE_SUM, nWgt,
439 sum_wgt_local, sum_wgt);
440 Teuchos::reduceAll<int,double>(*problemComm, Teuchos::REDUCE_MAX, nWgt,
441 max_wgt_local, max_wgt);
443 const double max_wgt_sum = double(std::numeric_limits<pm_idx_t>::max()/8);
444 for (
int j = 0; j < nWgt; j++) {
449 if (nonint[j] || (max_wgt[j]<=INT_EPSILON) || (sum_wgt[j]>max_wgt_sum)) {
451 if (sum_wgt[j] != 0.) scale = max_wgt_sum/sum_wgt[j];
455 for (
size_t i = 0; i < n; i++)
456 iwgts[i*nWgt+j] = (pm_idx_t) ceil(
double(fwgts[j][i])*scale);
458 delete [] nonint_local;
459 delete [] sum_wgt_local;
464 #endif // PARMETIS VERSION 4 OR HIGHER CHECK
466 #endif // HAVE_ZOLTAN2_MPI
468 #endif // HAVE_ZOLTAN2_PARMETIS
virtual void partition(const RCP< PartitioningSolution< Adapter > > &)
Partitioning method.
Defines the PartitioningSolution class.
static void SAVE_ARRAYRCP(ArrayRCP< first_t > *a, second_t *b, size_t size)
SparseMatrixAdapter_t::part_t part_t
Adapter::scalar_t scalar_t
AlgParMETIS(const RCP< const Environment > &, const RCP< const Comm< int > > &, const RCP< GraphModel< typename Adapter::base_adapter_t > > &)
Algorithm defines the base class for all algorithms.
static void ASSIGN(first_t &a, second_t b)
GraphModel defines the interface required for graph models.
static void ASSIGN_ARRAY(first_t **a, ArrayView< second_t > &b)
Defines the GraphModel interface.
A gathering of useful namespace methods.
static void DELETE_ARRAY(first_t **a)