50 #ifndef ZOLTAN2_EVALUATEORDERING_HPP 51 #define ZOLTAN2_EVALUATEORDERING_HPP 61 template <
typename Adapter>
65 typedef typename Adapter::lno_t lno_t;
66 typedef typename Adapter::gno_t gno_t;
68 typedef typename Adapter::scalar_t scalar_t;
80 void sharedConstructor(
const Adapter *ia,
82 const RCP<
const Comm<int> > &problemComm,
101 RCP<const Comm<int> > problemComm = DefaultComm<int>::getComm();
102 sharedConstructor(ia, p, problemComm, localSoln, globalSoln);
117 const RCP<
const Comm<int> > &problemComm,
121 sharedConstructor(ia, p, problemComm, localSoln, globalSoln);
124 #ifdef HAVE_ZOLTAN2_MPI 141 RCP<Teuchos::OpaqueWrapper<MPI_Comm> > wrapper =
142 Teuchos::opaqueWrapper(comm);
143 RCP<const Comm<int> > problemComm =
144 rcp<const Comm<int> >(
new Teuchos::MpiComm<int>(wrapper));
145 sharedConstructor(ia, p, problemComm, localSoln, globalSoln);
159 os <<
"Ordering Metrics" << std::endl;
160 os << std::setw(20) <<
" " << std::setw(11) <<
"ordered" << std::endl;
161 os << std::setw(20) <<
"envelope" << std::setw(11) << std::setprecision(4)
162 << envelope << std::endl;
163 os << std::setw(20) <<
"bandwidth" << std::setw(11) << std::setprecision(4)
164 << bandwidth << std::endl;
168 const RCP<const Environment> &env,
169 const RCP<
const Comm<int> > &comm,
178 std::bitset<NUM_MODEL_FLAGS> modelFlags;
179 RCP<GraphModel<base_adapter_t> > graph;
180 const RCP<const base_adapter_t> bia =
181 rcp(dynamic_cast<const base_adapter_t *>(ia),
false);
183 ArrayView<const gno_t> Ids;
184 ArrayView<input_t> vwgts;
185 ArrayView<const gno_t> edgeIds;
186 ArrayView<const lno_t> offsets;
187 ArrayView<input_t> wgts;
188 ArrayView<input_t> vtx;
189 graph->getEdgeList(edgeIds, offsets, wgts);
190 lno_t numVertex = graph->getVertexList(Ids, vwgts);
197 for(
int checkRank = 0; checkRank < comm->getSize(); ++checkRank ) {
199 if( checkRank == comm->getRank() ) {
200 std::cout <<
"-----------------------------------------" << std::endl;
201 std::cout <<
"Inspect rank: " << checkRank << std::endl;
202 std::cout << std::endl;
204 Array<lno_t> oldMatrix(numVertex*numVertex);
205 Array<lno_t> newMatrix(numVertex*numVertex);
208 std::cout << std::endl <<
"perm: ";
209 for(lno_t n = 0; n < numVertex; ++n) {
210 std::cout <<
" " << perm[n] <<
" ";
214 std::cout << std::endl <<
"iperm: ";
215 for(lno_t n = 0; n < numVertex; ++n) {
216 std::cout <<
" " << iperm[n] <<
" ";
218 std::cout << std::endl;
220 for (lno_t y = 0; y < numVertex; y++) {
221 for (lno_t n = offsets[y]; n < offsets[y+1]; ++n) {
222 lno_t x =
static_cast<lno_t
>(edgeIds[n]);
223 if (x < numVertex && y < numVertex) {
224 oldMatrix[x + y*numVertex] = 1;
225 newMatrix[perm[x] + perm[y]*numVertex] = 1;
231 std::cout << std::endl <<
"original graph in matrix form:" << std::endl;
232 for(lno_t y = 0; y < numVertex; ++y) {
233 for(lno_t x = 0; x < numVertex; ++x) {
234 std::cout <<
" " << oldMatrix[x + y*numVertex];
236 std::cout << std::endl;
240 std::cout << std::endl <<
"reordered graph in matrix form:" << std::endl;
241 for(lno_t y = 0; y < numVertex; ++y) {
242 for(lno_t x = 0; x < numVertex; ++x) {
243 std::cout <<
" " << newMatrix[x + y*numVertex];
245 std::cout << std::endl;
247 std::cout << std::endl;
253 #endif // Ends temporary logging which can be deleted later 260 for (lno_t j = 0; j < numVertex; j++) {
262 for (
auto n = offsets[j]; n < offsets[j+1]; ++n) {
263 lno_t x =
static_cast<lno_t
>(edgeIds[n]);
269 lno_t delta_right = y2 - x2;
270 if (delta_right > bw_right) {
271 bw_right = delta_right;
273 lno_t delta_left = y2 - x2;
274 if (delta_left > bw_left) {
275 bw_left = delta_left;
279 if(delta_right > 0) {
280 envelope += delta_right;
283 envelope += delta_left;
290 bandwidth = (bw_left + bw_right + 1);
300 template <
typename Adapter>
304 const RCP<
const Comm<int> > &comm,
308 RCP<const Comm<int> > problemComm = (comm == Teuchos::null) ?
309 DefaultComm<int>::getComm() : comm;
311 RCP<Environment> env;
318 env->timerStart(
MACRO_TIMERS,
"Computing ordering metrics");
329 throw std::logic_error(
"EvaluateOrdering not set up for global ordering.");
333 env->timerStop(
MACRO_TIMERS,
"Computing ordering metrics");
338 template <
typename Adapter>
341 typedef typename Adapter::lno_t lno_t;
368 const RCP<
const Comm<int> > &problemComm,
372 #ifdef HAVE_ZOLTAN2_MPI 391 template <
typename Adapter>
394 typedef typename Adapter::gno_t gno_t;
421 const RCP<
const Comm<int> > &problemComm,
425 #ifdef HAVE_ZOLTAN2_MPI Zoltan2::BaseAdapter< userTypes_t > base_adapter_t
Time an algorithm (or other entity) as a whole.
EvaluateGlobalOrdering(const Adapter *ia, ParameterList *p, const RCP< const Comm< int > > &problemComm, const GlobalOrderingSolution< gno_t > *globalSoln)
Constructor where Teuchos communicator is specified.
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
EvaluateOrdering(const Adapter *ia, ParameterList *p, const LocalOrderingSolution< lno_t > *localSoln, const GlobalOrderingSolution< gno_t > *globalSoln)
Constructor where communicator is Teuchos default.
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.
Defines the OrderingSolution class.
A class that computes and returns quality metrics. base class for the local and global ordering vers...
sub-steps, each method's entry and exit
lno_t getSeparatorSize() const
EvaluateLocalOrdering(const Adapter *ia, ParameterList *p, const LocalOrderingSolution< lno_t > *localSoln)
Constructor where communicator is Teuchos default.
SparseMatrixAdapter_t::part_t part_t
The StridedData class manages lists of weights or coordinates.
void localOrderingMetrics(const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, const Adapter *ia, const LocalOrderingSolution< typename Adapter::lno_t > *localSoln)
The user parameters, debug, timing and memory profiling output objects, and error checking methods...
Base class for the EvaluatePartition and EvaluateOrdering classes.
lno_t getEnvelope() const
GraphModel defines the interface required for graph models.
virtual void printMetrics(std::ostream &os) const
Print all metrics of type metricType based on the metric object type Note that parent class currently...
lno_t * getPermutationView(bool inverse=false) const
Get pointer to permutation. If inverse = true, return inverse permutation. By default, perm[i] is where new index i can be found in the old ordering. When inverse==true, perm[i] is where old index i can be found in the new ordering.
EvaluateGlobalOrdering(const Adapter *ia, ParameterList *p, const GlobalOrderingSolution< gno_t > *globalSoln)
Constructor where communicator is Teuchos default.
lno_t getBandwidth() const
EvaluateLocalOrdering(const Adapter *ia, ParameterList *p, const RCP< const Comm< int > > &problemComm, const LocalOrderingSolution< lno_t > *localSoln)
Constructor where Teuchos communicator is specified.