Zoltan2
Loading...
Searching...
No Matches
fix4785.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
20#include <Tpetra_Map.hpp>
21#include <vector>
22#include <cstdlib>
23
24typedef Tpetra::Map<> myMap_t;
25typedef myMap_t::local_ordinal_type myLocalId_t;
26typedef myMap_t::global_ordinal_type myGlobalId_t;
27typedef double myScalar_t;
28
30template <typename A, typename B>
31void copyArrays(size_t n, A *&a, B *b)
32{
33 a = new A[n];
34 for (size_t i = 0; i < n; i++) a[i] = b[i];
35}
36
37template <typename A, typename B>
38void copyArrays(size_t n, A *&a, A *b)
39{
40 a = b;
41}
42
44template <typename A, typename B>
45void freeArrays(A *&a, B *b)
46{
47 delete [] a;
48}
49
50template <typename A, typename B>
51void freeArrays(A *&a, A *b)
52{
53 // no delete needed since only copied pointer
54}
55
58 const Teuchos::RCP<const Teuchos::Comm<int> > &comm,
59 Teuchos::ParameterList &params,
60 size_t localCount,
61 myGlobalId_t *globalIds,
62 myScalar_t *coords,
63 int &nFail
64)
65{
66 typedef Tpetra::Map<>::node_type myNode_t;
67
69 myNode_t> myTypes;
70 typedef Zoltan2::BasicVectorAdapter<myTypes> inputAdapter_t;
72
73 myScalar_t *sx = coords;
74 myScalar_t *sy = sx + localCount;
75 myScalar_t *sz = sy + localCount;
76
77 inputAdapter_t *ia = new inputAdapter_t(localCount,globalIds,sx,sy,sz,1,1,1);
78
81
82 problem->solve();
83
84 quality_t *metricObject = new quality_t(ia, &params, comm,
85 &problem->getSolution());
86 if (comm->getRank() == 0){
87
88 metricObject->printMetrics(std::cout);
89
90 double imb = metricObject->getObjectCountImbalance();
91 if (imb <= 1.01) // Should get perfect balance
92 std::cout << "no weights -- balance satisfied: " << imb << std::endl;
93 else {
94 std::cout << "no weights -- balance failure: " << imb << std::endl;
95 nFail++;
96 }
97 std::cout << std::endl;
98 }
99
100 delete metricObject;
101 delete problem;
102 delete ia;
103}
104
107 const Teuchos::RCP<const Teuchos::Comm<int> > &comm,
108 Teuchos::ParameterList &params,
109 size_t localCount,
110 myGlobalId_t *globalIds,
111 myScalar_t *coords,
113 int &nFail
114)
115{
116 typedef Tpetra::Map<>::node_type myNode_t;
117
119 myNode_t> myTypes;
120 typedef Zoltan2::BasicVectorAdapter<myTypes> inputAdapter_t;
122
123 std::vector<const myScalar_t *> coordVec(3);
124 std::vector<int> coordStrides(3);
125
126 coordVec[0] = coords; coordStrides[0] = 1;
127 coordVec[1] = coords + localCount; coordStrides[1] = 1;
128 coordVec[2] = coords + localCount + localCount; coordStrides[2] = 1;
129
130 std::vector<const myScalar_t *> weightVec(1);
131 std::vector<int> weightStrides(1);
132
133 weightVec[0] = weights; weightStrides[0] = 1;
134
135 inputAdapter_t *ia=new inputAdapter_t(localCount, globalIds, coordVec,
136 coordStrides,weightVec,weightStrides);
137
140
141 problem->solve();
142
143 quality_t *metricObject = new quality_t(ia, &params, comm,
144 &problem->getSolution());
145 if (comm->getRank() == 0){
146
147 metricObject->printMetrics(std::cout);
148
149 double imb = metricObject->getWeightImbalance(0);
150 if (imb <= 1.01)
151 std::cout << "weighted -- balance satisfied " << imb << std::endl;
152 else {
153 std::cout << "weighted -- balance failed " << imb << std::endl;
154 nFail++;
155 }
156 std::cout << std::endl;
157 }
158
159 delete metricObject;
160 delete problem;
161 delete ia;
162}
163
165int main(int narg, char *arg[])
166{
167 Tpetra::ScopeGuard scope(&narg, &arg);
168 const Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
169 int me = comm->getRank();
170 int np = comm->getSize();
171 int nFail = 0;
172
174 // Create parameters for an MJ problem
175
176 Teuchos::ParameterList params("test params");
177 params.set("debug_level", "basic_status");
178 params.set("error_check_level", "debug_mode_assertions");
179
180 params.set("algorithm", "multijagged");
181 params.set("num_global_parts", 64);
182
184 // Create input data.
185
186 const size_t N = 9000000;
187 size_t maxLocalCount = N / np; // biggest test we'll run
188
189 // Create coordinates that range from 0 to 999
190 const int dim = 3;
191 myScalar_t *coords = new myScalar_t[dim * maxLocalCount];
192
193 srand(me);
194 for (size_t i=0; i < maxLocalCount*dim; i++)
195 coords[i] = myScalar_t(0);
196
197 // Create weights for the coordinates
198 myScalar_t *weights = new myScalar_t[maxLocalCount];
199 for (size_t i=0; i < maxLocalCount; i++)
200 weights[i] = myScalar_t(1+me);
201
202 // Allocate space for global IDs; they will be generated below
203 myGlobalId_t *globalIds = new myGlobalId_t[maxLocalCount];
204
205 size_t localCount = N / np;
206
207 // Create consecutive global ids for the coordinates
208 myGlobalId_t offset = me * localCount;
209 for (size_t i=0; i < localCount; i++)
210 globalIds[i] = offset++;
211
212 if (me == 0) {
213 std::cout << "---------------------------------------------" << std::endl;
214 std::cout << "myGlobalId_t = " << typeid(offset).name()
215 << " " << sizeof(offset)
216 << "; localCount = " << localCount
217 << "; globalCount = " << np * localCount << std::endl;
218 }
219
221 // Test one: No weights
222 if (me == 0) std::cout << "Test: no weights, scalar = double" << std::endl;
223 test_no_weights(comm, params, localCount, globalIds, coords, nFail);
224
226 // Test two: weighted
227 if (me == 0) std::cout << "Test: weights, scalar = double" << std::endl;
228 test_weights(comm, params, localCount, globalIds, coords, weights, nFail);
229
230 // Early exit when failure is detected
231 int gnFail;
232 Teuchos::reduceAll<int, int>(*comm, Teuchos::REDUCE_MAX, 1, &nFail, &gnFail);
233
234 delete [] weights;
235 delete [] coords;
236 delete [] globalIds;
237
238 if (me == 0) {
239 if (nFail == 0) std::cout << "PASS" << std::endl;
240 else std::cout << "FAIL: " << nFail << " tests failed" << std::endl;
241 }
242
243 return 0;
244}
245
Defines the BasicVectorAdapter class.
Traits for application input objects.
Defines the PartitioningProblem class.
Defines the PartitioningSolution class.
int main()
A simple class that can be the User template argument for an InputAdapter.
BasicVectorAdapter represents a vector (plus optional weights) supplied by the user as pointers to st...
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.
myMap_t::local_ordinal_type myLocalId_t
Definition fix4785.cpp:25
double myScalar_t
Definition fix4785.cpp:27
void freeArrays(A *&a, B *b)
Definition fix4785.cpp:45
void test_no_weights(const Teuchos::RCP< const Teuchos::Comm< int > > &comm, Teuchos::ParameterList &params, size_t localCount, myGlobalId_t *globalIds, myScalar_t *coords, int &nFail)
Definition fix4785.cpp:57
Tpetra::Map myMap_t
Definition fix4785.cpp:24
void copyArrays(size_t n, A *&a, B *b)
Definition fix4785.cpp:31
void test_weights(const Teuchos::RCP< const Teuchos::Comm< int > > &comm, Teuchos::ParameterList &params, size_t localCount, myGlobalId_t *globalIds, myScalar_t *coords, myScalar_t *weights, int &nFail)
Definition fix4785.cpp:106
myMap_t::global_ordinal_type myGlobalId_t
Definition fix4785.cpp:26
static ArrayRCP< ArrayRCP< zscalar_t > > weights
Zoltan2::EvaluatePartition< matrixAdapter_t > quality_t