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]);
148 template <
typename User>
162 : nids(nids_), dim(dim_)
166 typedef Kokkos::DualView<gno_t *> view_ids_t;
167 kokkos_gids = view_ids_t(
168 Kokkos::ViewAllocateWithoutInitializing(
"gids"), nids);
170 auto host_gids = kokkos_gids.h_view;
171 for(
size_t n = 0; n < nids; ++n) {
172 host_gids(n) = gids_[n];
175 kokkos_gids.template modify<typename view_ids_t::host_mirror_space>();
176 kokkos_gids.sync_host();
177 kokkos_gids.template sync<typename node_t::device_type>();
183 typedef Kokkos::DualView<scalar_t **> view_weights_t;
184 kokkos_weights = view_weights_t(
185 Kokkos::ViewAllocateWithoutInitializing(
"weights"), nids, 0);
186 auto host_kokkos_weights = kokkos_weights.h_view;
187 for(
size_t n = 0; n < nids; ++n) {
188 host_kokkos_weights(n,0) = weights_[n];
191 kokkos_weights.template modify<typename view_weights_t::host_mirror_space>();
192 kokkos_weights.sync_host();
193 kokkos_weights.template sync<typename node_t::device_type>();
198 typedef Kokkos::DualView<scalar_t **, Kokkos::LayoutLeft> kokkos_coords_t;
199 kokkos_coords = kokkos_coords_t(
200 Kokkos::ViewAllocateWithoutInitializing(
"coords"), nids, dim);
201 auto host_kokkos_coords = kokkos_coords.h_view;
203 for(
size_t n = 0; n < nids; ++n) {
204 for(
int idx = 0; idx < dim; ++idx) {
205 host_kokkos_coords(n,idx) = coords_[n+idx*nids];
209 kokkos_coords.template modify<typename kokkos_coords_t::host_mirror_space>();
210 kokkos_coords.sync_host();
211 kokkos_coords.template sync<typename node_t::device_type>();
218 auto kokkosIds = kokkos_gids.view_host();
219 ids = kokkosIds.data();
223 typename node_t::device_type> &ids)
const {
224 ids = kokkos_gids.template view<typename node_t::device_type>();;
230 int idx = 0)
const override
232 auto h_wgts_2d = kokkos_weights.view_host();
234 wgt = Kokkos::subview(h_wgts_2d, Kokkos::ALL, idx).data();
239 typename node_t::device_type> & wgt)
const {
240 wgt = kokkos_weights.template view<typename node_t::device_type>();;
246 int &stride,
int idx = 0)
const override {
247 elements = kokkos_coords.view_host().data();
252 Kokkos::LayoutLeft,
typename node_t::device_type> & coo)
const {
253 coo = kokkos_coords.template view<typename node_t::device_type>();
258 Kokkos::DualView<gno_t *> kokkos_gids;
260 Kokkos::DualView<scalar_t **, Kokkos::LayoutLeft> kokkos_coords;
261 Kokkos::DualView<scalar_t **> kokkos_weights;
269 const Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
271 int rank = comm->getRank();
272 int nprocs = comm->getSize();
274 typedef Tpetra::Map<> Map_t;
275 typedef Map_t::local_ordinal_type localId_t;
276 typedef Map_t::global_ordinal_type globalId_t;
277 typedef double scalar_t;
288 size_t localCount = 40;
292 scalar_t *cStrided =
new scalar_t [dim * localCount];
294 for (
size_t i = 0; i < localCount; i++)
295 for (
int d = 0; d < dim; d++)
296 cStrided[cnt++] = d*1000 + rank*100 + i;
299 scalar_t *cContig =
new scalar_t [dim * localCount];
301 for (
int d = 0; d < dim; d++)
302 for (
size_t i = 0; i < localCount; i++)
303 cContig[cnt++] = d*1000 + rank*100 + i;
306 globalId_t *globalIds =
new globalId_t [localCount];
307 globalId_t offset = rank * localCount;
308 for (
size_t i=0; i < localCount; i++) globalIds[i] = offset++;
313 Teuchos::ParameterList params(
"test params");
314 params.set(
"debug_level",
"basic_status");
315 params.set(
"error_check_level",
"debug_mode_assertions");
317 params.set(
"algorithm", algorithm);
318 params.set(
"num_global_parts", nprocs+1);
324 stridedAdapter_t *ia1 =
325 new stridedAdapter_t(localCount,globalIds,dim,cStrided);
340 std::cout <<
"no weights -- balance satisfied: " << imb << std::endl;
342 std::cout <<
"no weights -- balance failure: " << imb << std::endl;
345 std::cout << std::endl;
349 contigAdapter_t *ia2 =
new contigAdapter_t(localCount,globalIds,dim,cContig);
357 kokkosAdapter_t *ia3 =
new kokkosAdapter_t(localCount,globalIds,dim,cContig);
366 for (
size_t i = 0; i < localCount; i++) {
367 if((problem1->
getSolution().getPartListView()[i] !=
372 std::cout << rank <<
" Error: differing parts for index " << i
373 << problem1->
getSolution().getPartListView()[i] <<
" "
374 << problem2->
getSolution().getPartListView()[i] <<
" "
375 << problem3->
getSolution().getPartListView()[i] << std::endl;
380 if (ndiff > 0) nFail++;
381 else if (rank == 0) std::cout <<
"no weights -- comparisons OK " << std::endl;
383 delete metricObject1;
395 scalar_t *
weights =
new scalar_t [localCount];
396 for (
size_t i=0; i < localCount; i++)
weights[i] = 1 + scalar_t(rank);
399 ia1 =
new stridedAdapter_t(localCount, globalIds, dim, cStrided,
weights);
413 std::cout <<
"weighted -- balance satisfied " << imb << std::endl;
415 std::cout <<
"weighted -- balance failed " << imb << std::endl;
418 std::cout << std::endl;
422 ia2 =
new contigAdapter_t(localCount, globalIds, dim, cContig,
weights);
430 for (
size_t i = 0; i < localCount; i++) {
431 if (problem1->
getSolution().getPartListView()[i] !=
433 std::cout << rank <<
" Error: differing parts for index " << i
434 << problem1->
getSolution().getPartListView()[i] <<
" "
435 << problem2->
getSolution().getPartListView()[i] << std::endl;
440 if (ndiff > 0) nFail++;
441 else if (rank == 0) std::cout <<
"weighted -- comparisons OK " << std::endl;
443 delete metricObject1;
451 if (cStrided)
delete [] cStrided;
452 if (cContig)
delete [] cContig;
453 if (globalIds)
delete [] globalIds;
458 Teuchos::reduceAll(*comm, Teuchos::REDUCE_SUM, 1, &nFail, &gnFail);
463 int main(
int narg,
char *arg[])
465 Tpetra::ScopeGuard scope(&narg, &arg);
466 const Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
473 if (comm->getRank() == 0) {
474 if (err == 0) std::cout <<
"PASS" << std::endl;
475 else std::cout <<
"FAIL: " << err <<
" tests failed" << std::endl;
Defines the PartitioningProblem class.
Defines the PartitioningSolution class.
Defines the VectorAdapter interface.
int getNumEntriesPerID() const
Return the number of vectors.
void getIDsView(const gno_t *&ids) const override
Tpetra::Map ::node_type node_t
void getEntriesView(const scalar_t *&elements, int &stride, int idx=0) const override
size_t getLocalNumIDs() const
Returns the number of objects on this process.
int getNumWeightsPerID() const
Returns the number of weights per object. Number of weights per object should be zero or greater....
virtual void getWeightsKokkosView(Kokkos::View< scalar_t **, typename node_t::device_type > &wgt) const
KokkosVectorAdapter(const size_t nids_, const gno_t *gids_, const int dim_, const scalar_t *coords_, const scalar_t *weights_=NULL)
virtual void getEntriesKokkosView(Kokkos::View< scalar_t **, Kokkos::LayoutLeft, typename node_t::device_type > &coo) const
Zoltan2::InputTraits< User >::gno_t gno_t
virtual void getIDsKokkosView(Kokkos::View< const gno_t *, typename node_t::device_type > &ids) const
Zoltan2::InputTraits< User >::scalar_t scalar_t
void getWeightsView(const scalar_t *&wgt, int &stride, int idx=0) const override
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.
void printMetrics(std::ostream &os) const
Print all metrics.
scalar_t getWeightImbalance(int weightIndex) const
Return the imbalance for the requested weight.
scalar_t getObjectCountImbalance() const
Return the object count imbalance.
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.
map_t::global_ordinal_type gno_t
int run_test_strided_versus_contig(const std::string &algorithm)
int main(int narg, char *arg[])
Zoltan2::EvaluatePartition< matrixAdapter_t > quality_t