43 #include "Ifpack_ConfigDefs.h"
44 #include "Ifpack_Reordering.h"
45 #include "Ifpack_METISReordering.h"
46 #include "Ifpack_Graph.h"
47 #include "Ifpack_Graph_Epetra_CrsGraph.h"
48 #include "Ifpack_Graph_Epetra_RowMatrix.h"
49 #include "Epetra_Comm.h"
50 #include "Epetra_MultiVector.h"
51 #include "Epetra_CrsGraph.h"
52 #include "Epetra_Map.h"
53 #include "Teuchos_ParameterList.hpp"
56 #ifdef HAVE_IFPACK_METIS
58 void METIS_NodeND(
int *n, idxtype *xadj, idxtype *adjncy,
59 int *numflag,
int *options,
int *perm,
int *iperm);
65 UseSymmetricGraph_(false),
81 NumMyRows_ = Graph.NumMyRows();
82 Reorder_.resize(NumMyRows_);
83 InvReorder_.resize(NumMyRows_);
87 Teuchos::RefCountPtr<Epetra_CrsGraph> SymGraph;
88 Teuchos::RefCountPtr<Epetra_Map> SymMap;
89 Teuchos::RefCountPtr<Ifpack_Graph_Epetra_CrsGraph> SymIFPACKGraph;
90 Teuchos::RefCountPtr<Ifpack_Graph> IFPACKGraph = Teuchos::rcp( (
Ifpack_Graph*)&Graph,
false );
92 int Length = 2 * Graph.MaxMyNumEntries();
94 std::vector<int> Indices;
95 Indices.resize(Length);
97 std::vector<int> options;
101 #ifdef HAVE_IFPACK_METIS
105 if (UseSymmetricGraph_) {
107 #if !defined(EPETRA_NO_32BIT_GLOBAL_INDICES) || !defined(EPETRA_NO_64BIT_GLOBAL_INDICES)
112 SymMap = Teuchos::rcp(
new Epetra_Map(NumMyRows_,0,Graph.Comm()) );
116 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
117 if(SymGraph->RowMap().GlobalIndicesInt()) {
118 for (
int i = 0; i < NumMyRows_ ; ++i) {
120 ierr = Graph.ExtractMyRowCopy(i, Length, NumIndices,
122 IFPACK_CHK_ERR(ierr);
124 for (
int j = 0 ; j < NumIndices ; ++j) {
128 SymGraph->InsertGlobalIndices(i,1,&jj);
129 SymGraph->InsertGlobalIndices(jj,1,&i);
136 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
137 if(SymGraph->RowMap().GlobalIndicesLongLong()) {
138 for (
int i = 0; i < NumMyRows_ ; ++i) {
141 ierr = Graph.ExtractMyRowCopy(i, Length, NumIndices,
143 IFPACK_CHK_ERR(ierr);
145 for (
int j = 0 ; j < NumIndices ; ++j) {
146 long long jj = Indices[j];
149 SymGraph->InsertGlobalIndices(i_LL,1,&jj);
150 SymGraph->InsertGlobalIndices(jj,1,&i_LL);
157 throw "Ifpack_METISReordering::Compute: GlobalIndices type unknown";
159 IFPACK_CHK_ERR(SymGraph->OptimizeStorage());
160 IFPACK_CHK_ERR(SymGraph->FillComplete());
162 IFPACKGraph = SymIFPACKGraph;
166 std::vector<idxtype> xadj;
167 xadj.resize(NumMyRows_ + 1);
169 std::vector<idxtype> adjncy;
170 adjncy.resize(Graph.NumMyNonzeros());
176 for (
int i = 0; i < NumMyRows_ ; ++i) {
178 xadj[count2+1] = xadj[count2];
180 ierr = IFPACKGraph->ExtractMyRowCopy(i, Length, NumIndices, &Indices[0]);
181 IFPACK_CHK_ERR(ierr);
183 for (
int j = 0 ; j < NumIndices ; ++j) {
186 adjncy[count++] = jj;
193 #ifdef HAVE_IFPACK_METIS
200 METIS_NodeND(&NumMyRows_, &xadj[0], &adjncy[0],
201 &numflag, &options[0],
202 &InvReorder_[0], &Reorder_[0]);
207 cerr <<
"Please configure with --enable-ifpack-metis" << endl;
208 cerr <<
"to use METIS Reordering." << endl;
220 IFPACK_CHK_ERR(
Compute(Graph));
231 if ((i < 0) || (i >= NumMyRows_))
244 if ((i < 0) || (i >= NumMyRows_))
248 return(InvReorder_[i]);
256 for (
int j = 0 ; j < NumVectors ; ++j) {
257 for (
int i = 0 ; i < NumMyRows_ ; ++i) {
258 int np = Reorder_[i];
259 X[j][np] = Xorig[j][i];
272 for (
int j = 0 ; j < NumVectors ; ++j) {
273 for (
int i = 0 ; i < NumMyRows_ ; ++i) {
274 int np = Reorder_[i];
275 X[j][i] = Xorig[j][np];
287 os <<
"*** Ifpack_METISReordering" << endl << endl;
289 os <<
"*** Reordering not yet computed." << endl;
291 os <<
"*** Number of local rows = " << NumMyRows_ << endl;
292 os <<
"Local Row\tReorder[i]\tInvReorder[i]" << endl;
293 for (
int i = 0 ; i < NumMyRows_ ; ++i) {
294 os <<
'\t' << i <<
"\t\t" << Reorder_[i] <<
"\t\t" << InvReorder_[i] << endl;