Zoltan2
Loading...
Searching...
No Matches
Zoltan2_EvaluateOrdering.hpp
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#ifndef ZOLTAN2_EVALUATEORDERING_HPP
15#define ZOLTAN2_EVALUATEORDERING_HPP
16
19
20namespace Zoltan2{
21
25template <typename Adapter>
26class EvaluateOrdering : public EvaluateBaseClass<Adapter> {
27
28private:
29 typedef typename Adapter::lno_t lno_t;
30 typedef typename Adapter::gno_t gno_t;
31 typedef typename Adapter::offset_t offset_t;
32 typedef typename Adapter::part_t part_t;
33 typedef typename Adapter::scalar_t scalar_t;
34 typedef typename Adapter::base_adapter_t base_adapter_t;
35
36 // To do - this is only appropriate for the local ordering
37 // Need to decide how to organize these classes
38 // Do we potentially want local + global in this class
39 // Do we want to eliminate the EvaluateLocalOrdering and EvaluateGlobalOrdering
40 // derived classes? Or perhaps make them completely independent of each other
41 lno_t bandwidth;
42 lno_t envelope;
43 lno_t separatorSize;
44
45 void sharedConstructor(const Adapter *ia,
46 ParameterList *p,
47 const RCP<const Comm<int> > &problemComm,
48 const LocalOrderingSolution<lno_t> *localSoln,
49 const GlobalOrderingSolution<gno_t> *globalSoln);
50public:
51
61 const Adapter *ia,
62 ParameterList *p,
63 const LocalOrderingSolution<lno_t> *localSoln,
64 const GlobalOrderingSolution<gno_t> *globalSoln)
65 {
66 Teuchos::RCP<const Comm<int> > problemComm = Tpetra::getDefaultComm();
67 sharedConstructor(ia, p, problemComm, localSoln, globalSoln);
68 }
69
80 const Adapter *ia,
81 ParameterList *p,
82 const RCP<const Comm<int> > &problemComm,
83 const LocalOrderingSolution<lno_t> *localSoln,
84 const GlobalOrderingSolution<gno_t> *globalSoln)
85 {
86 sharedConstructor(ia, p, problemComm, localSoln, globalSoln);
87 }
88
89#ifdef HAVE_ZOLTAN2_MPI
100 const Adapter *ia,
101 ParameterList *p,
102 MPI_Comm comm,
103 const LocalOrderingSolution<lno_t> *localSoln,
104 const GlobalOrderingSolution<gno_t> *globalSoln)
105 {
106 RCP<Teuchos::OpaqueWrapper<MPI_Comm> > wrapper =
107 Teuchos::opaqueWrapper(comm);
108 RCP<const Comm<int> > problemComm =
109 rcp<const Comm<int> >(new Teuchos::MpiComm<int>(wrapper));
110 sharedConstructor(ia, p, problemComm, localSoln, globalSoln);
111 }
112#endif
113
114 lno_t getBandwidth() const { return bandwidth; }
115 lno_t getEnvelope() const { return envelope; }
116 lno_t getSeparatorSize() const { return separatorSize; }
117
121 virtual void printMetrics(std::ostream &os) const {
122
123 // To Do - complete this formatting
124 os << "Ordering Metrics" << std::endl;
125 os << std::setw(20) << " " << std::setw(11) << "ordered" << std::endl;
126 os << std::setw(20) << "envelope" << std::setw(11) << std::setprecision(4)
127 << envelope << std::endl;
128 os << std::setw(20) << "bandwidth" << std::setw(11) << std::setprecision(4)
129 << bandwidth << std::endl;
130 }
131
133 const RCP<const Environment> &env,
134 const RCP<const Comm<int> > &comm,
135 const Adapter *ia,
137 {
138 env->debug(DETAILED_STATUS, "Entering orderingMetrics"); // begin
139
140 typedef StridedData<lno_t, scalar_t> input_t;
141
142 // get graph
143 std::bitset<NUM_MODEL_FLAGS> modelFlags;
144 RCP<GraphModel<base_adapter_t> > graph;
145 const RCP<const base_adapter_t> bia =
146 rcp(dynamic_cast<const base_adapter_t *>(ia), false);
147 graph = rcp(new GraphModel<base_adapter_t>(bia,env,comm,modelFlags));
148 ArrayView<const gno_t> Ids;
149 ArrayView<input_t> vwgts;
150 ArrayView<const gno_t> edgeIds;
151 ArrayView<const offset_t> offsets;
152 ArrayView<input_t> wgts;
153 ArrayView<input_t> vtx;
154 graph->getEdgeList(edgeIds, offsets, wgts);
155 lno_t numVertex = graph->getVertexList(Ids, vwgts);
156
157 lno_t * perm = localSoln->getPermutationView();
158
159 // print as matrix - this was debugging code which can be deleted later
160 #define MDM
161 #ifdef MDM
162 for( int checkRank = 0; checkRank < comm->getSize(); ++checkRank ) {
163 comm->barrier();
164 if( checkRank == comm->getRank() ) {
165 std::cout << "-----------------------------------------" << std::endl;
166 std::cout << "Inspect rank: " << checkRank << std::endl;
167 std::cout << std::endl;
168 if(numVertex < 30) { // don't spam if it's too many...
169 Array<lno_t> oldMatrix(numVertex*numVertex);
170 Array<lno_t> newMatrix(numVertex*numVertex);
171
172 // print the solution permutation
173 std::cout << std::endl << "perm: ";
174 for(lno_t n = 0; n < numVertex; ++n) {
175 std::cout << " " << perm[n] << " ";
176 }
177
178 lno_t * iperm = localSoln->getPermutationView(true);
179 std::cout << std::endl << "iperm: ";
180 for(lno_t n = 0; n < numVertex; ++n) {
181 std::cout << " " << iperm[n] << " ";
182 }
183 std::cout << std::endl;
184 // write 1's to old matrix (original form) and new matrix (using solution)
185 for (lno_t y = 0; y < numVertex; y++) {
186 for (offset_t n = offsets[y]; n < offsets[y+1]; ++n) {
187 lno_t x = static_cast<lno_t>(edgeIds[n]); // to resolve
188 if (x < numVertex && y < numVertex) { // to develop - for MPI this may not be local
189 oldMatrix[x + y*numVertex] = 1;
190 newMatrix[perm[x] + perm[y]*numVertex] = 1;
191 }
192 }
193 }
194
195 // print oldMatrix
196 std::cout << std::endl << "original graph in matrix form:" << std::endl;
197 for(lno_t y = 0; y < numVertex; ++y) {
198 for(lno_t x = 0; x < numVertex; ++x) {
199 std::cout << " " << oldMatrix[x + y*numVertex];
200 }
201 std::cout << std::endl;
202 }
203
204 // print newMatrix
205 std::cout << std::endl << "reordered graph in matrix form:" << std::endl;
206 for(lno_t y = 0; y < numVertex; ++y) {
207 for(lno_t x = 0; x < numVertex; ++x) {
208 std::cout << " " << newMatrix[x + y*numVertex];
209 }
210 std::cout << std::endl;
211 }
212 std::cout << std::endl;
213 }
214 }
215
216 comm->barrier();
217 }
218 #endif // Ends temporary logging which can be deleted later
219
220 // calculate bandwidth and envelope for unsolved and solved case
221 lno_t bw_right = 0;
222 lno_t bw_left = 0;
223 envelope = 0;
224
225 for (lno_t j = 0; j < numVertex; j++) {
226 lno_t y = Ids[j];
227 for (offset_t n = offsets[j]; n < offsets[j+1]; ++n) {
228 lno_t x = static_cast<lno_t>(edgeIds[n]); // to resolve
229 if(x < numVertex) {
230 lno_t x2 = perm[x];
231 lno_t y2 = perm[y];
232
233 // solved bandwidth calculation
234 lno_t delta_right = y2 - x2;
235 if (delta_right > bw_right) {
236 bw_right = delta_right;
237 }
238 lno_t delta_left = y2 - x2;
239 if (delta_left > bw_left) {
240 bw_left = delta_left;
241 }
242
243 // solved envelope calculation
244 if(delta_right > 0) {
245 envelope += delta_right;
246 }
247 if(delta_left > 0) {
248 envelope += delta_left;
249 }
250 envelope += 1; // need to check this - do we count center?
251 }
252 }
253 }
254
255 bandwidth = (bw_left + bw_right + 1);
256
257 // TO DO - No implementation yet for this metric
258 separatorSize = 0;
259
260 env->debug(DETAILED_STATUS, "Exiting orderingMetrics"); // end
261 }
262};
263
264// sharedConstructor
265template <typename Adapter>
266void EvaluateOrdering<Adapter>::sharedConstructor(
267 const Adapter *ia,
268 ParameterList *p,
269 const RCP<const Comm<int> > &comm,
270 const LocalOrderingSolution<lno_t> *localSoln,
271 const GlobalOrderingSolution<gno_t> *globalSoln)
272{
273 RCP<const Comm<int> > problemComm = (comm == Teuchos::null) ?
274 Tpetra::getDefaultComm() : comm;
275
276 RCP<Environment> env;
277 try{
278 env = rcp(new Environment(*p, problemComm));
279 }
281
282 env->debug(DETAILED_STATUS, std::string("Entering EvaluateOrdering"));
283 env->timerStart(MACRO_TIMERS, "Computing ordering metrics");
284
285 try{
286 // May want to move these into the specific derived classes
287 // But it depends on whether we eventually may have both types and perhaps
288 // want to combine the metrics
289 if(localSoln) {
290 localOrderingMetrics(env, problemComm, ia, localSoln);
291 }
292
293 if(globalSoln) {
294 throw std::logic_error("EvaluateOrdering not set up for global ordering.");
295 }
296 }
298 env->timerStop(MACRO_TIMERS, "Computing ordering metrics");
299
300 env->debug(DETAILED_STATUS, std::string("Exiting EvaluateOrdering"));
301}
302
303template <typename Adapter>
305private:
306 typedef typename Adapter::lno_t lno_t;
307
308public:
317 const Adapter *ia,
318 ParameterList *p,
319 const LocalOrderingSolution<lno_t> *localSoln) :
320 EvaluateOrdering<Adapter>(ia, p, localSoln, nullptr) {}
321
331 const Adapter *ia,
332 ParameterList *p,
333 const RCP<const Comm<int> > &problemComm,
334 const LocalOrderingSolution<lno_t> *localSoln) :
335 EvaluateOrdering<Adapter>(ia, p, problemComm, localSoln, nullptr) {}
336
337#ifdef HAVE_ZOLTAN2_MPI
348 const Adapter *ia,
349 ParameterList *p,
350 MPI_Comm comm,
351 const LocalOrderingSolution<lno_t> *localSoln) :
352 EvaluateOrdering<Adapter>(ia, p, comm, localSoln, nullptr) {}
353#endif
354};
355
356template <typename Adapter>
358private:
359 typedef typename Adapter::gno_t gno_t;
360
361public:
370 const Adapter *ia,
371 ParameterList *p,
372 const GlobalOrderingSolution<gno_t> *globalSoln) :
373 EvaluateOrdering<Adapter>(ia, p, nullptr, globalSoln) {}
374
384 const Adapter *ia,
385 ParameterList *p,
386 const RCP<const Comm<int> > &problemComm,
387 const GlobalOrderingSolution<gno_t> *globalSoln) :
388 EvaluateOrdering<Adapter>(ia, p, problemComm, nullptr, globalSoln) {}
389
390#ifdef HAVE_ZOLTAN2_MPI
401 const Adapter *ia,
402 ParameterList *p,
403 MPI_Comm comm,
404 const GlobalOrderingSolution<gno_t> *globalSoln) :
405 EvaluateOrdering<Adapter>(ia, p, comm, nullptr, globalSoln) {}
406#endif
407};
408
409} // namespace Zoltan2
410
411#endif
Base class for the EvaluatePartition and EvaluateOrdering classes.
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
Defines the OrderingSolution class.
EvaluateGlobalOrdering(const Adapter *ia, ParameterList *p, const RCP< const Comm< int > > &problemComm, const GlobalOrderingSolution< gno_t > *globalSoln)
Constructor where Teuchos communicator is specified.
EvaluateGlobalOrdering(const Adapter *ia, ParameterList *p, const GlobalOrderingSolution< gno_t > *globalSoln)
Constructor where communicator is Teuchos default.
EvaluateLocalOrdering(const Adapter *ia, ParameterList *p, const LocalOrderingSolution< lno_t > *localSoln)
Constructor where communicator is Teuchos default.
EvaluateLocalOrdering(const Adapter *ia, ParameterList *p, const RCP< const Comm< int > > &problemComm, const LocalOrderingSolution< lno_t > *localSoln)
Constructor where Teuchos communicator is specified.
A class that computes and returns quality metrics. \A base class for the local and global ordering ve...
EvaluateOrdering(const Adapter *ia, ParameterList *p, const LocalOrderingSolution< lno_t > *localSoln, const GlobalOrderingSolution< gno_t > *globalSoln)
Constructor where communicator is Teuchos default.
void localOrderingMetrics(const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, const Adapter *ia, const LocalOrderingSolution< typename Adapter::lno_t > *localSoln)
EvaluateOrdering(const Adapter *ia, ParameterList *p, const RCP< const Comm< int > > &problemComm, const LocalOrderingSolution< lno_t > *localSoln, const GlobalOrderingSolution< gno_t > *globalSoln)
Constructor where Teuchos communicator is specified.
virtual void printMetrics(std::ostream &os) const
Print all metrics of type metricType based on the metric object type Note that parent class currently...
GraphModel defines the interface required for graph models.
index_t * getPermutationView(bool inverse=false) const
Get pointer to permutation. If inverse = true, return inverse permutation. By default,...
The StridedData class manages lists of weights or coordinates.
Created by mbenlioglu on Aug 31, 2020.
@ MACRO_TIMERS
Time an algorithm (or other entity) as a whole.
@ DETAILED_STATUS
sub-steps, each method's entry and exit