15 #include <Teuchos_DefaultComm.hpp>
16 #include <Teuchos_CommandLineProcessor.hpp>
17 #include <Teuchos_ParameterList.hpp>
20 int main(
int narg,
char **arg)
22 Tpetra::ScopeGuard tscope(&narg, &arg);
23 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
27 int rank = comm->getRank();
28 int nprocs = comm->getSize();
31 int numGlobalIdentifiers = 100;
32 int numMyIdentifiers = numGlobalIdentifiers / nprocs;
33 if (rank < numGlobalIdentifiers % nprocs)
34 numMyIdentifiers += 1;
41 if (!myIds || !myWeights){
49 std::cout <<
"Memory allocation failure" << std::endl;
50 std::cout <<
"FAIL" << std::endl;
56 for (
int i=0; i < numMyIdentifiers; i++){
57 myIds[i] = myBaseId+i;
58 myWeights[i] = rank%3 + 1;
59 origsumwgts += myWeights[i];
63 int *origcnt =
new int[nprocs];
65 Teuchos::gather<int, int>(&numMyIdentifiers, 1, origcnt, 1, 0, *comm);
66 Teuchos::gather<int, zscalar_t>(&origsumwgts, 1, origwgts, 1, 0, *comm);
68 std::cout <<
"BEFORE PART CNTS: ";
69 for (
int i = 0; i < nprocs; i++)
70 std::cout << origcnt[i] <<
" ";
71 std::cout << std::endl;
72 std::cout <<
"BEFORE PART WGTS: ";
73 for (
int i = 0; i < nprocs; i++)
74 std::cout << origwgts[i] <<
" ";
75 std::cout << std::endl;
81 std::vector<const zscalar_t *> weightValues;
82 std::vector<int> weightStrides;
83 weightValues.push_back(const_cast<const zscalar_t *>(myWeights));
91 new adapter_t(
zlno_t(numMyIdentifiers),myIds,weightValues,weightStrides);
94 bool useWeights =
true;
95 Teuchos::CommandLineProcessor cmdp (
false,
false);
96 cmdp.setOption(
"weights",
"no-weights", &useWeights,
97 "Indicated whether to use identifier weights in partitioning");
98 cmdp.parse(narg, arg);
100 Teuchos::ParameterList params(
"test parameters");
102 params.set(
"num_global_parts", nprocs);
103 params.set(
"algorithm",
"block");
104 params.set(
"partitioning_approach",
"partition");
105 if (!useWeights) params.set(
"partitioning_objective",
"balance_object_count");
115 quality_t *metricObject =
new quality_t(adapter, ¶ms, comm, &solution);
120 memset(totalWeight, 0, nprocs *
sizeof(
zscalar_t));
121 int *totalCnt =
new int [nprocs];
122 int *sumCnt =
new int [nprocs];
123 memset(totalCnt, 0, nprocs *
sizeof(
int));
130 for (
int i=0; i < numMyIdentifiers; i++){
131 totalCnt[partList[i]]++;
132 totalWeight[partList[i]] += myWeights[i];
135 Teuchos::reduceAll<int, int>(*comm, Teuchos::REDUCE_SUM, nprocs,
137 Teuchos::reduceAll<int, zscalar_t>(*comm, Teuchos::REDUCE_SUM, nprocs,
138 totalWeight, sumWeight);
143 std::cout <<
"AFTER PART CNTS: ";
144 for (
int i=0; i < nprocs; i++)
145 std::cout << sumCnt[i] <<
" ";
146 std::cout << std::endl;
149 std::cout <<
"AFTER PART WGTS: ";
150 for (
int i=0; i < nprocs; i++){
151 std::cout << sumWeight[i] <<
" ";
152 total += sumWeight[i];
154 std::cout << std::endl;
160 for (
int i=0; i < nprocs; i++){
162 if (sumWeight[i] > avg)
163 imb = (sumWeight[i] - avg) / avg;
165 imb = (avg - sumWeight[i]) / avg;
172 std::cout <<
"Computed imbalance: " << imbalance << std::endl;
173 std::cout <<
"Library's imbalance: " << libImbalance << std::endl;
176 if (imbalance > libImbalance)
177 err = imbalance - libImbalance;
179 err = libImbalance - imbalance;
192 std::cout <<
"failure in solution's imbalance data" << std::endl;
193 std::cout <<
"FAIL" << std::endl;
199 std::cout <<
"PASS" << std::endl;
206 delete [] totalWeight;
Zoltan2::EvaluatePartition< matrixAdapter_t > quality_t
int main(int narg, char **arg)
Defines the PartitioningSolution class.
common code used by tests
scalar_t getWeightImbalance(int weightIndex) const
Return the imbalance for the requested weight.
SparseMatrixAdapter_t::part_t part_t
This class represents a collection of global Identifiers and their associated weights, if any.
A PartitioningSolution is a solution to a partitioning problem.
const part_t * getPartListView() const
Returns the part list corresponding to the global ID list.
Tpetra::Map::local_ordinal_type zlno_t
static const std::string fail
Defines the BasicIdentifierAdapter class.
int globalFail(const Comm< int > &comm, int fail)
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.
Tpetra::Map::global_ordinal_type zgno_t
void solve(bool updateInputData=true)
Direct the problem to create a solution.