14 #include <Kokkos_Core.hpp>
15 #include <Teuchos_DefaultComm.hpp>
29 int main(
int narg,
char *arg[]) {
30 Tpetra::ScopeGuard tscope(&narg, &arg);
31 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
33 int rank = comm->getRank();
34 int nprocs = comm->getSize();
38 typedef Tpetra::Map<> Map_t;
39 typedef Map_t::local_ordinal_type localId_t;
40 typedef Map_t::global_ordinal_type globalId_t;
41 typedef Tpetra::Details::DefaultTypes::scalar_type scalar_t;
42 typedef Tpetra::Map<>::node_type
node_t;
47 int localCount = 40 * (rank + 1);
48 int totalCount = 20 * nprocs * (nprocs + 1);
49 int targetCount = totalCount / nprocs;
51 Kokkos::View<globalId_t*, typename node_t::device_type>
52 globalIds(Kokkos::ViewAllocateWithoutInitializing(
"globalIds"), localCount);
53 auto host_globalIds = Kokkos::create_mirror_view(globalIds);
56 for (
int i = 0, num = 40; i < nprocs ; i++, num += 40) {
57 std::cout <<
"Rank " << i <<
" generates " << num <<
" ids." << std::endl;
61 globalId_t offset = 0;
62 for (
int i = 1; i <= rank; i++) {
66 for (
int i = 0; i < localCount; i++) {
67 host_globalIds(i) = offset++;
69 Kokkos::deep_copy(globalIds, host_globalIds);
77 const int nWeights = 1;
78 Kokkos::View<scalar_t **, typename node_t::device_type>
79 weights(
"weights", localCount, nWeights);
80 auto host_weights = Kokkos::create_mirror_view(weights);
82 for (
int index = 0; index < localCount; index++) {
83 host_weights(index, 0) = 1;
86 Kokkos::deep_copy(weights, host_weights);
88 inputAdapter_t ia(globalIds, weights);
93 Teuchos::ParameterList params(
"test params");
94 params.set(
"debug_level",
"basic_status");
95 params.set(
"debug_procs",
"0");
96 params.set(
"error_check_level",
"debug_mode_assertions");
98 params.set(
"algorithm",
"block");
99 params.set(
"imbalance_tolerance", 1.1);
100 params.set(
"num_global_parts", nprocs);
117 Kokkos::View<const globalId_t *, typename node_t::device_type> ids;
118 ia.getIDsKokkosView(ids);
120 auto host_ids = Kokkos::create_mirror_view(ids);
121 Kokkos::deep_copy(host_ids, ids);
123 Kokkos::View<int*, Kokkos::HostSpace> partCounts(
"partCounts", nprocs);
125 Kokkos::View<int*, Kokkos::HostSpace>
126 globalPartCounts(
"globalPartCounts", nprocs);
128 for (
size_t i = 0; i < ia.getLocalNumIDs(); i++) {
129 int pp = problem->
getSolution().getPartListView()[i];
130 std::cout << rank <<
" LID " << i <<
" GID " << host_ids(i)
131 <<
" PART " << pp << std::endl;
135 Teuchos::reduceAll<int, int>(*comm, Teuchos::REDUCE_SUM, nprocs,
136 partCounts.data(), globalPartCounts.data());
140 for (
int i = 0; i < nprocs; i++) {
141 if (globalPartCounts(i) != targetCount) {
142 std::cout <<
"FAIL: part " << i <<
" has " << globalPartCounts(i)
143 <<
" != " << targetCount <<
"; " << ++ierr <<
" errors"
148 std::cout <<
"PASS" << std::endl;
A simple class that can be the User template argument for an InputAdapter.
int main(int narg, char **arg)
Defines the PartitioningSolution class.
This class represents a collection of global Identifiers and their associated weights, if any.
Defines the BasicKokkosIdentifierAdapter class.
const PartitioningSolution< Adapter > & getSolution()
Get the solution to the problem.
PartitioningProblem sets up partitioning problems for the user.
Defines the PartitioningProblem class.
void solve(bool updateInputData=true)
Direct the problem to create a solution.