54 #include <Tpetra_Map.hpp>
60 template <
typename User>
73 : nids(nids_), gids(gids_), dim(dim_), coords(coords_), weights(weights_)
83 { wgt = weights; stride = 1; }
89 coo = &(coords[idx*nids]);
103 template <
typename User>
116 : nids(nids_), gids(gids_), dim(dim_), coords(coords_), weights(weights_)
126 { wgt = weights; stride = 1; }
132 coo = &(coords[idx]);
145 int main(
int narg,
char *arg[])
147 Tpetra::ScopeGuard scope(&narg, &arg);
148 const Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
149 int rank = comm->getRank();
150 int nprocs = comm->getSize();
153 typedef Tpetra::Map<> Map_t;
154 typedef Map_t::local_ordinal_type localId_t;
155 typedef Map_t::global_ordinal_type globalId_t;
156 typedef double scalar_t;
166 size_t localCount = 40;
170 scalar_t *cStrided =
new scalar_t [dim * localCount];
172 for (
size_t i = 0; i < localCount; i++)
173 for (
int d = 0; d < dim; d++)
174 cStrided[cnt++] = d*1000 + rank*100 + i;
177 scalar_t *cContig =
new scalar_t [dim * localCount];
179 for (
int d = 0; d < dim; d++)
180 for (
size_t i = 0; i < localCount; i++)
181 cContig[cnt++] = d*1000 + rank*100 + i;
184 globalId_t *globalIds =
new globalId_t [localCount];
185 globalId_t offset = rank * localCount;
186 for (
size_t i=0; i < localCount; i++) globalIds[i] = offset++;
191 Teuchos::ParameterList params(
"test params");
192 params.set(
"debug_level",
"basic_status");
193 params.set(
"error_check_level",
"debug_mode_assertions");
195 params.set(
"algorithm",
"multijagged");
196 params.set(
"num_global_parts", nprocs+1);
202 stridedAdapter_t *ia1 =
203 new stridedAdapter_t(localCount,globalIds,dim,cStrided);
210 quality_t *metricObject1 =
new quality_t(ia1, ¶ms, comm,
214 metricObject1->printMetrics(std::cout);
216 double imb = metricObject1->getObjectCountImbalance();
218 std::cout <<
"no weights -- balance satisfied: " << imb << std::endl;
220 std::cout <<
"no weights -- balance failure: " << imb << std::endl;
223 std::cout << std::endl;
227 contigAdapter_t *ia2 =
new contigAdapter_t(localCount,globalIds,dim,cContig);
236 for (
size_t i = 0; i < localCount; i++) {
237 if (problem1->
getSolution().getPartListView()[i] !=
239 std::cout << rank <<
" Error: differing parts for index " << i
240 << problem1->
getSolution().getPartListView()[i] <<
" "
241 << problem2->
getSolution().getPartListView()[i] << std::endl;
246 if (ndiff > 0) nFail++;
247 else if (rank == 0) std::cout <<
"no weights -- comparisons OK " << std::endl;
249 delete metricObject1;
259 scalar_t *
weights =
new scalar_t [localCount];
260 for (
size_t i=0; i < localCount; i++)
weights[i] = 1 + scalar_t(rank);
263 ia1 =
new stridedAdapter_t(localCount, globalIds, dim, cStrided,
weights);
269 metricObject1 =
new quality_t(ia1, ¶ms, comm, &problem1->
getSolution());
273 metricObject1->printMetrics(std::cout);
275 double imb = metricObject1->getWeightImbalance(0);
277 std::cout <<
"weighted -- balance satisfied " << imb << std::endl;
279 std::cout <<
"weighted -- balance failed " << imb << std::endl;
282 std::cout << std::endl;
286 ia2 =
new contigAdapter_t(localCount, globalIds, dim, cContig,
weights);
294 for (
size_t i = 0; i < localCount; i++) {
295 if (problem1->
getSolution().getPartListView()[i] !=
297 std::cout << rank <<
" Error: differing parts for index " << i
298 << problem1->
getSolution().getPartListView()[i] <<
" "
299 << problem2->
getSolution().getPartListView()[i] << std::endl;
304 if (ndiff > 0) nFail++;
305 else if (rank == 0) std::cout <<
"weighted -- comparisons OK " << std::endl;
307 delete metricObject1;
315 if (cStrided)
delete [] cStrided;
316 if (cContig)
delete [] cContig;
317 if (globalIds)
delete [] globalIds;
322 Teuchos::reduceAll(*comm, Teuchos::REDUCE_SUM, 1, &nFail, &gnFail);
325 if (gnFail == 0) std::cout <<
"PASS" << std::endl;
326 else std::cout <<
"FAIL: " << gnFail <<
" tests failed" << std::endl;
Defines the PartitioningProblem class.
Defines the PartitioningSolution class.
Defines the VectorAdapter interface.
Zoltan2::InputTraits< User >::scalar_t scalar_t
int getNumEntriesPerID() const
Return the number of vectors.
int getNumWeightsPerID() const
Returns the number of weights per object. Number of weights per object should be zero or greater....
OldSchoolVectorAdapterContig(const size_t nids_, const gno_t *gids_, const int dim_, const scalar_t *coords_, const scalar_t *weights_=NULL)
void getIDsView(const gno_t *&ids) const
size_t getLocalNumIDs() const
Returns the number of objects on this process.
Zoltan2::InputTraits< User >::gno_t gno_t
void getWeightsView(const scalar_t *&wgt, int &stride, int idx=0) const
void getEntriesView(const scalar_t *&coo, int &stride, int idx=0) const
Zoltan2::InputTraits< User >::scalar_t scalar_t
int getNumWeightsPerID() const
Returns the number of weights per object. Number of weights per object should be zero or greater....
size_t getLocalNumIDs() const
Returns the number of objects on this process.
int getNumEntriesPerID() const
Return the number of vectors.
void getIDsView(const gno_t *&ids) const
void getEntriesView(const scalar_t *&coo, int &stride, int idx=0) const
void getWeightsView(const scalar_t *&wgt, int &stride, int idx=0) const
OldSchoolVectorAdapterStrided(const size_t nids_, const gno_t *gids_, const int dim_, const scalar_t *coords_, const scalar_t *weights_=NULL)
Zoltan2::InputTraits< User >::gno_t gno_t
A simple class that can be the User template argument for an InputAdapter.
A class that computes and returns quality metrics.
PartitioningProblem sets up partitioning problems for the user.
const PartitioningSolution< Adapter > & getSolution()
Get the solution to the problem.
void solve(bool updateInputData=true)
Direct the problem to create a solution.
VectorAdapter defines the interface for vector input.
int main(int narg, char *arg[])