59 using Teuchos::ArrayRCP;
66 void doTest(RCP<
const Comm<int> > comm,
int numLocalObj,
67 int nWeights,
int numLocalParts,
bool givePartSizes);
69 int main(
int narg,
char *arg[])
71 Tpetra::ScopeGuard tscope(&narg, &arg);
72 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
74 int rank = comm->getRank();
76 doTest(comm, 10, 0, -1,
false);
77 doTest(comm, 10, 0, 1,
false);
78 doTest(comm, 10, 0, 1,
true);
79 doTest(comm, 10, 1, 1,
false);
80 doTest(comm, 10, 1, 1,
true);
81 doTest(comm, 10, 2, 1,
false);
82 doTest(comm, 10, 2, 1,
true);
83 doTest(comm, 10, 1, 2,
true);
84 doTest(comm, 10, 1, 2,
false);
85 doTest(comm, 10, 1, -1,
false);
86 doTest(comm, 10, 1, -1,
true);
87 doTest(comm, 10, 2, -1,
false);
90 std::cout <<
"PASS" << std::endl;
95 void doTest(RCP<
const Comm<int> > comm,
int numLocalObj,
96 int nWeights,
int numLocalParts,
bool givePartSizes)
103 int rank = comm->getRank();
104 int nprocs = comm->getSize();
107 bool testEmptyParts = (numLocalParts < 1);
108 int numGlobalParts = 0;
111 numGlobalParts = nprocs / 2;
112 if (numGlobalParts >= 1)
113 numLocalParts = (rank < numGlobalParts ? 1 : 0);
116 testEmptyParts =
false;
120 numGlobalParts = nprocs * numLocalParts;
124 std::cout << std::endl;
125 std::cout <<
"Test: number of weights " << nWeights;
126 std::cout <<
", desired number of parts " << numGlobalParts;
128 std::cout <<
", with differing part sizes." << std::endl;
130 std::cout <<
", with uniform part sizes." << std::endl;
131 std::cout <<
"Number of procs " << nprocs;
132 std::cout <<
", each with " << numLocalObj <<
" objects, part = rank." << std::endl;
137 Teuchos::ParameterList pl(
"test list");
138 pl.set(
"num_local_parts", numLocalParts);
140 RCP<const Zoltan2::Environment> env =
146 for (
int i=0, x=rank*numLocalObj; i < numLocalObj; i++, x++){
153 int partSizeDim = (givePartSizes ? (nWeights ? nWeights : 1) : 0);
154 ArrayRCP<ArrayRCP<part_t> > ids(partSizeDim);
155 ArrayRCP<ArrayRCP<zscalar_t> > sizes(partSizeDim);
157 if (givePartSizes && numLocalParts > 0){
159 myParts[0] = rank * numLocalParts;
160 for (
int i=1; i < numLocalParts; i++)
161 myParts[i] = myParts[i-1] + 1;
162 ArrayRCP<part_t> partNums(myParts, 0, numLocalParts,
true);
165 if (sizeFactor < 0) sizeFactor *= -1;
168 for (
int dim=0; dim < partSizeDim; dim++){
170 for (
int i=0; i < numLocalParts; i++)
171 psizes[i] = sizeFactor;
172 sizes[dim] = arcp(psizes, 0, numLocalParts,
true);
180 std::vector<const zscalar_t *>
weights;
181 std::vector<int> strides;
183 int len = numLocalObj*nWeights;
184 ArrayRCP<zscalar_t> wgtBuf;
189 wgtBuf = arcp(wgts, 0, len,
true);
190 for (
int i=0; i < len; i++)
194 for (
int i=0; i < nWeights; i++, wgts+=numLocalObj)
202 catch (std::exception &e){
210 RCP<Zoltan2::PartitioningSolution<idInput_t> > solution;
216 ids.view(0,partSizeDim), sizes.view(0,partSizeDim)));
219 env, comm, nWeights));
221 catch (std::exception &e){
230 ArrayRCP<part_t> partAssignment(partNum, 0, numLocalObj,
true);
231 for (
int i=0; i < numLocalObj; i++)
234 solution->setParts(partAssignment);
238 RCP<quality_t> metricObject;
241 metricObject = rcp(
new quality_t(ia, &pl, comm, solution.getRawPtr()));
243 catch (std::exception &e){
253 zscalar_t imb = metricObject->getObjectCountImbalance();
254 std::cout <<
"Object imbalance: " << imb << std::endl;
256 catch (std::exception &e){
263 if (rank==0 && nWeights > 0){
265 for (
int i=0; i < nWeights; i++){
266 zscalar_t imb = metricObject->getWeightImbalance(i);
267 std::cout <<
"Weight " << i <<
" imbalance: " << imb << std::endl;
270 catch (std::exception &e){
273 if (!
fail && nWeights > 1){
275 zscalar_t imb = metricObject->getNormedImbalance();
276 std::cout <<
"Normed weight imbalance: " << imb << std::endl;
278 catch (std::exception &e){
288 metricObject->printMetrics(std::cout);
290 catch (std::exception &e){