Zoltan2
Loading...
Searching...
No Matches
blockTest.cpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Zoltan2: A package of combinatorial algorithms for scientific computing
4//
5// Copyright 2012 NTESS and the Zoltan2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
14
15#include <Teuchos_DefaultComm.hpp>
16#include <Teuchos_CommandLineProcessor.hpp>
17#include <Teuchos_ParameterList.hpp>
18
19
20int main(int narg, char **arg)
21{
22 Tpetra::ScopeGuard tscope(&narg, &arg);
23 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
24
25 int fail=0, gfail=0;
26
27 int rank = comm->getRank();
28 int nprocs = comm->getSize();
29
30 // Construct the user input
31 int numGlobalIdentifiers = 100;
32 int numMyIdentifiers = numGlobalIdentifiers / nprocs;
33 if (rank < numGlobalIdentifiers % nprocs)
34 numMyIdentifiers += 1;
35
36 zgno_t myBaseId = zgno_t(numGlobalIdentifiers * rank);
37
38 zgno_t *myIds = new zgno_t[numMyIdentifiers];
39 zscalar_t *myWeights = new zscalar_t[numMyIdentifiers];
40
41 if (!myIds || !myWeights){
42 fail = 1;
43 }
44
45 gfail = globalFail(*comm, fail);
46
47 if (gfail){
48 if (rank==0){
49 std::cout << "Memory allocation failure" << std::endl;
50 std::cout << "FAIL" << std::endl;
51 }
52 return 1;
53 }
54
55 zscalar_t origsumwgts = 0;
56 for (int i=0; i < numMyIdentifiers; i++){
57 myIds[i] = myBaseId+i;
58 myWeights[i] = rank%3 + 1;
59 origsumwgts += myWeights[i];
60 }
61
62 // Some output
63 int *origcnt = new int[nprocs];
64 zscalar_t *origwgts = new zscalar_t[nprocs];
65 Teuchos::gather<int, int>(&numMyIdentifiers, 1, origcnt, 1, 0, *comm);
66 Teuchos::gather<int, zscalar_t>(&origsumwgts, 1, origwgts, 1, 0, *comm);
67 if (rank == 0) {
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;
76 }
77 delete [] origcnt;
78 delete [] origwgts;
79
80 // Building Zoltan2 adapters
81 std::vector<const zscalar_t *> weightValues;
82 std::vector<int> weightStrides; // default is one
83 weightValues.push_back(const_cast<const zscalar_t *>(myWeights));
84
88 typedef adapter_t::part_t part_t;
89
90 adapter_t *adapter =
91 new adapter_t(zlno_t(numMyIdentifiers),myIds,weightValues,weightStrides);
92
93 // Set up the parameters and problem
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);
99
100 Teuchos::ParameterList params("test parameters");
101 //params.set("compute_metrics", true); // bool parameter
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");
106
107 Zoltan2::PartitioningProblem<adapter_t> problem(adapter, &params);
108
109 problem.solve();
110
112
113 // create metric object
114
115 quality_t *metricObject = new quality_t(adapter, &params, comm, &solution);
116
117 // Some output
118 zscalar_t *totalWeight = new zscalar_t [nprocs];
119 zscalar_t *sumWeight = new zscalar_t [nprocs];
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));
124
125 const part_t *partList = solution.getPartListView();
126 zscalar_t libImbalance = metricObject->getWeightImbalance(0);
127 delete metricObject;
128 delete adapter;
129
130 for (int i=0; i < numMyIdentifiers; i++){
131 totalCnt[partList[i]]++;
132 totalWeight[partList[i]] += myWeights[i];
133 }
134
135 Teuchos::reduceAll<int, int>(*comm, Teuchos::REDUCE_SUM, nprocs,
136 totalCnt, sumCnt);
137 Teuchos::reduceAll<int, zscalar_t>(*comm, Teuchos::REDUCE_SUM, nprocs,
138 totalWeight, sumWeight);
139
140 double epsilon = 10e-6;
141
142 if (rank == 0){
143 std::cout << "AFTER PART CNTS: ";
144 for (int i=0; i < nprocs; i++)
145 std::cout << sumCnt[i] << " ";
146 std::cout << std::endl;
147
148 zscalar_t total = 0;
149 std::cout << "AFTER PART WGTS: ";
150 for (int i=0; i < nprocs; i++){
151 std::cout << sumWeight[i] << " ";
152 total += sumWeight[i];
153 }
154 std::cout << std::endl;
155
156 zscalar_t avg = total / zscalar_t(nprocs);
157
158 zscalar_t imbalance = -1.0;
159
160 for (int i=0; i < nprocs; i++){
161 zscalar_t imb = 0;
162 if (sumWeight[i] > avg)
163 imb = (sumWeight[i] - avg) / avg;
164 else
165 imb = (avg - sumWeight[i]) / avg;
166
167 if (imb > imbalance)
168 imbalance = imb;
169 }
170 imbalance += 1.0;
171
172 std::cout << "Computed imbalance: " << imbalance << std::endl;
173 std::cout << "Library's imbalance: " << libImbalance << std::endl;
174
175 double err;
176 if (imbalance > libImbalance)
177 err = imbalance - libImbalance;
178 else
179 err = libImbalance - imbalance;
180
181 if (err > epsilon)
182 fail = 1;
183 }
184 else{
185 fail = 0;
186 }
187
188 gfail = globalFail(*comm, fail);
189
190 if (gfail){
191 if (rank==0){
192 std::cout << "failure in solution's imbalance data" << std::endl;
193 std::cout << "FAIL" << std::endl;
194 }
195 return 1;
196 }
197
198 if (rank==0)
199 std::cout << "PASS" << std::endl;
200
201 delete [] myWeights;
202 delete [] myIds;
203 delete [] sumCnt;
204 delete [] totalCnt;
205 delete [] sumWeight;
206 delete [] totalWeight;
207
208 return 0;
209}
int globalFail(const Comm< int > &comm, int fail)
Defines the BasicIdentifierAdapter class.
Defines the PartitioningProblem class.
Defines the PartitioningSolution class.
common code used by tests
float zscalar_t
Tpetra::Map ::local_ordinal_type zlno_t
Tpetra::Map ::global_ordinal_type zgno_t
int main()
This class represents a collection of global Identifiers and their associated weights,...
A simple class that can be the User template argument for an InputAdapter.
A class that computes and returns quality metrics.
scalar_t getWeightImbalance(int weightIndex) const
Return the imbalance for the requested weight.
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.
A PartitioningSolution is a solution to a partitioning problem.
const part_t * getPartListView() const
Returns the part list corresponding to the global ID list.
static const std::string fail
#define epsilon
Definition nd.cpp:47
SparseMatrixAdapter_t::part_t part_t
Zoltan2::EvaluatePartition< matrixAdapter_t > quality_t