19 #include <Tpetra_Map.hpp>
23 int main(
int narg,
char *arg[])
25 Tpetra::ScopeGuard scope(&narg, &arg);
26 const Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
27 int rank = comm->getRank();
28 int nprocs = comm->getSize();
31 typedef Tpetra::Map<> Map_t;
32 typedef Map_t::local_ordinal_type localId_t;
33 typedef Map_t::global_ordinal_type globalId_t;
34 typedef int int_scalar_t;
44 size_t localCount = 40;
48 int_scalar_t *coords =
new int_scalar_t [dim * localCount];
49 int_scalar_t *x = coords;
50 int_scalar_t *y = x + localCount;
51 int_scalar_t *z = y + localCount;
54 for (
size_t i=0; i < localCount*dim; i++)
55 coords[i] = int_scalar_t(rand() % 10);
58 globalId_t *globalIds =
new globalId_t [localCount];
59 globalId_t offset = rank * localCount;
60 for (
size_t i=0; i < localCount; i++) globalIds[i] = offset++;
65 Teuchos::ParameterList params(
"test params");
66 params.set(
"debug_level",
"basic_status");
67 params.set(
"error_check_level",
"debug_mode_assertions");
69 params.set(
"algorithm",
"multijagged");
70 params.set(
"num_global_parts", nprocs+1);
75 inputAdapter_t *ia1 =
new inputAdapter_t(localCount,globalIds,x,y,z,1,1,1);
82 quality_t *metricObject1 =
new quality_t(ia1, ¶ms, comm,
90 std::cout <<
"no weights -- balance satisfied: " << imb << std::endl;
92 std::cout <<
"no weights -- balance failure: " << imb << std::endl;
95 std::cout << std::endl;
105 int_scalar_t *
weights =
new int_scalar_t [localCount];
106 for (
size_t i=0; i < localCount; i++) weights[i] = 1 + int_scalar_t(rank);
108 std::vector<const int_scalar_t *>coordVec(3);
109 std::vector<int> coordStrides(3);
111 coordVec[0] = x; coordStrides[0] = 1;
112 coordVec[1] = y; coordStrides[1] = 1;
113 coordVec[2] = z; coordStrides[2] = 1;
115 std::vector<const int_scalar_t *>weightVec(1);
116 std::vector<int> weightStrides(1);
118 weightVec[0] =
weights; weightStrides[0] = 1;
120 inputAdapter_t *ia2=
new inputAdapter_t(localCount, globalIds, coordVec,
121 coordStrides,weightVec,weightStrides);
128 quality_t *metricObject2 =
new quality_t(ia2, ¶ms, comm,
136 std::cout <<
"weighted -- balance satisfied " << imb << std::endl;
138 std::cout <<
"weighted -- balance failed " << imb << std::endl;
141 std::cout << std::endl;
143 delete metricObject2;
145 if (weights)
delete []
weights;
146 if (coords)
delete [] coords;
147 if (globalIds)
delete [] globalIds;
152 if (nFail == 0) std::cout <<
"PASS" << std::endl;
153 else std::cout <<
"FAIL: " << nFail <<
" tests failed" << std::endl;
void printMetrics(std::ostream &os) const
Print all metrics.
A simple class that can be the User template argument for an InputAdapter.
Zoltan2::EvaluatePartition< matrixAdapter_t > quality_t
int main(int narg, char **arg)
Defines the PartitioningSolution class.
scalar_t getWeightImbalance(int weightIndex) const
Return the imbalance for the requested weight.
BasicVectorAdapter represents a vector (plus optional weights) supplied by the user as pointers to st...
const PartitioningSolution< Adapter > & getSolution()
Get the solution to the problem.
PartitioningProblem sets up partitioning problems for the user.
Defines the PartitioningProblem class.
A class that computes and returns quality metrics.
Defines the BasicVectorAdapter class.
void solve(bool updateInputData=true)
Direct the problem to create a solution.
scalar_t getObjectCountImbalance() const
Return the object count imbalance.