Commit 3c199079 authored by Kirill Terekhov's avatar Kirill Terekhov
Browse files

Reverse-Cuthill-McKee reordering in MPTILUC preconditioner

Added Reverse-Cuthill-McKee reordering into INNER_MPTILUC preconditioner
as default reordering when Metis_NodeND is not provided. On the tests
RCM shows 10% more fill-in in inverse factors then with Metis_NodeND,
but the RCM reordering procedure is 4 times faster then Metis_NodeND and
overall computational time decreases with RCM.

Fix type mismatch warnings in autodiff.cpp
parent 9c538219
...@@ -255,8 +255,8 @@ namespace INMOST ...@@ -255,8 +255,8 @@ namespace INMOST
Storage::integer_array indarr = jt->IntegerArray(it->indices); Storage::integer_array indarr = jt->IntegerArray(it->indices);
for (Storage::integer_array::iterator qt = indarr.begin(); qt != indarr.end(); ++qt) for (Storage::integer_array::iterator qt = indarr.begin(); qt != indarr.end(); ++qt)
{ {
if( *qt < first_num ) Pre.insert(*qt); if( static_cast<INMOST_DATA_ENUM_TYPE>(*qt) < first_num ) Pre.insert(*qt);
else if( *qt >= last_num ) Post.insert(*qt); else if( static_cast<INMOST_DATA_ENUM_TYPE>(*qt) >= last_num ) Post.insert(*qt);
} }
} }
} }
...@@ -271,8 +271,8 @@ namespace INMOST ...@@ -271,8 +271,8 @@ namespace INMOST
indarr.resize(jt->RealArray(it->d.t).size()); indarr.resize(jt->RealArray(it->d.t).size());
for (Storage::integer_array::iterator qt = indarr.begin(); qt != indarr.end(); ++qt) for (Storage::integer_array::iterator qt = indarr.begin(); qt != indarr.end(); ++qt)
{ {
if( *qt < first_num ) Pre.insert(*qt); if( static_cast<INMOST_DATA_ENUM_TYPE>(*qt) < first_num ) Pre.insert(*qt);
else if( *qt >= last_num ) Post.insert(*qt); else if( static_cast<INMOST_DATA_ENUM_TYPE>(*qt) >= last_num ) Post.insert(*qt);
} }
} }
} }
...@@ -289,8 +289,8 @@ namespace INMOST ...@@ -289,8 +289,8 @@ namespace INMOST
Storage::integer_array indarr = jt->IntegerArray(it->indices); Storage::integer_array indarr = jt->IntegerArray(it->indices);
for (Storage::integer_array::iterator qt = indarr.begin(); qt != indarr.end(); ++qt) for (Storage::integer_array::iterator qt = indarr.begin(); qt != indarr.end(); ++qt)
{ {
if( *qt < first_num ) Pre.insert(*qt); if( static_cast<INMOST_DATA_ENUM_TYPE>(*qt) < first_num ) Pre.insert(*qt);
else if( *qt >= last_num ) Post.insert(*qt); else if( static_cast<INMOST_DATA_ENUM_TYPE>(*qt) >= last_num ) Post.insert(*qt);
} }
} }
} }
...@@ -304,8 +304,8 @@ namespace INMOST ...@@ -304,8 +304,8 @@ namespace INMOST
Storage::integer_array indarr = jt->IntegerArray(it->indices); Storage::integer_array indarr = jt->IntegerArray(it->indices);
for (Storage::integer_array::iterator qt = indarr.begin(); qt != indarr.end(); ++qt) for (Storage::integer_array::iterator qt = indarr.begin(); qt != indarr.end(); ++qt)
{ {
if( *qt < first_num ) Pre.insert(*qt); if( static_cast<INMOST_DATA_ENUM_TYPE>(*qt) < first_num ) Pre.insert(*qt);
else if( *qt >= last_num ) Post.insert(*qt); else if( static_cast<INMOST_DATA_ENUM_TYPE>(*qt) >= last_num ) Post.insert(*qt);
} }
} }
} }
......
...@@ -4,6 +4,7 @@ ...@@ -4,6 +4,7 @@
#if defined(USE_SOLVER) #if defined(USE_SOLVER)
#include "solver_mtiluc2.hpp" #include "solver_mtiluc2.hpp"
#include <sstream> #include <sstream>
#include <deque>
//#define REPORT_ILU //#define REPORT_ILU
//#undef REPORT_ILU //#undef REPORT_ILU
//#define REPORT_ILU_PROGRESS //#define REPORT_ILU_PROGRESS
...@@ -20,7 +21,7 @@ using namespace INMOST; ...@@ -20,7 +21,7 @@ using namespace INMOST;
#endif #endif
//#define REORDER_RCM #define REORDER_RCM
//#define REORDER_NNZ //#define REORDER_NNZ
#if defined(USE_SOLVER_METIS) #if defined(USE_SOLVER_METIS)
#define REORDER_METIS_ND #define REORDER_METIS_ND
...@@ -77,6 +78,20 @@ __INLINE static bool compare(INMOST_DATA_REAL_TYPE * a, INMOST_DATA_REAL_TYPE * ...@@ -77,6 +78,20 @@ __INLINE static bool compare(INMOST_DATA_REAL_TYPE * a, INMOST_DATA_REAL_TYPE *
(*reinterpret_cast< make_integer<INMOST_DATA_REAL_TYPE>::type * >(b)); (*reinterpret_cast< make_integer<INMOST_DATA_REAL_TYPE>::type * >(b));
} }
//this structure is used to provide sorting routing for Reverse-Cuthill-McKee reordering
struct RCM_Comparator
{
INMOST_DATA_ENUM_TYPE wbeg;
std::vector<INMOST_DATA_ENUM_TYPE> & xadj;
public:
RCM_Comparator(INMOST_DATA_ENUM_TYPE _wbeg, std::vector<INMOST_DATA_ENUM_TYPE> & _xadj)
: wbeg(_wbeg), xadj(_xadj) {}
RCM_Comparator(const RCM_Comparator & b) : wbeg(b.wbeg), xadj(b.xadj) {}
RCM_Comparator & operator = (RCM_Comparator const & b) {wbeg = b.wbeg; xadj = b.xadj; return *this;}
bool operator () (INMOST_DATA_ENUM_TYPE i, INMOST_DATA_ENUM_TYPE j)
{return xadj[i+1-wbeg] -xadj[i-wbeg] < xadj[j+1-wbeg] - xadj[j-wbeg];}
};
class BinaryHeap class BinaryHeap
{ {
...@@ -507,6 +522,9 @@ public: ...@@ -507,6 +522,9 @@ public:
double tfactor = 0.0, trescale = 0.0, treorder = 0.0, treassamble = 0.0, ttotal, tt; double tfactor = 0.0, trescale = 0.0, treorder = 0.0, treassamble = 0.0, ttotal, tt;
#if defined(REORDER_METIS_ND) #if defined(REORDER_METIS_ND)
double tmetisgraph = 0, tmetisnd = 0; double tmetisgraph = 0, tmetisnd = 0;
#endif
#if defined(REORDER_RCM)
double trcmgraph = 0, trcmorder = 0;
#endif #endif
double tlfactor, tlrescale, tlreorder, tlreassamble; double tlfactor, tlrescale, tlreorder, tlreassamble;
ttotal = Timer(); ttotal = Timer();
...@@ -936,6 +954,7 @@ public: ...@@ -936,6 +954,7 @@ public:
/////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
/////////// END MAXIMUM TRANSVERSE REORDERING ///////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////// /////////// END MAXIMUM TRANSVERSE REORDERING /////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
/////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
tt = Timer();
#if defined(REORDER_DDPQ) #if defined(REORDER_DDPQ)
#if defined(REPORT_ILU) #if defined(REPORT_ILU)
...@@ -1152,7 +1171,6 @@ public: ...@@ -1152,7 +1171,6 @@ public:
cend = wend; cend = wend;
i = wend; i = wend;
treorder += Timer() - tt;
#elif defined(REORDER_NNZ) #elif defined(REORDER_NNZ)
//for (k = cbeg; k < cend; ++k) V[k] = DR[localQ[k]]; //for (k = cbeg; k < cend; ++k) V[k] = DR[localQ[k]];
//for (k = cbeg; k < cend; ++k) DR[k] = V[k]; //for (k = cbeg; k < cend; ++k) DR[k] = V[k];
...@@ -1202,31 +1220,153 @@ public: ...@@ -1202,31 +1220,153 @@ public:
} }
cend = i; cend = i;
#elif defined(REORDER_RCM) #elif defined(REORDER_RCM)
//compute order of each entry {
std::fill(Ulist.begin() + wbeg - mobeg, Ulist.begin() + wend - mobeg, 0); //create a symmetric graph of the matrix A + A^T
std::fill(Llist.begin() + wbeg - mobeg, Llist.begin() + wend - mobeg, 0); std::vector<INMOST_DATA_ENUM_TYPE> xadj(wend-wbeg+1), adjncy;
for (k = wbeg; k < wend; ++k) adjncy.reserve(nzA*2);
{
for (INMOST_DATA_ENUM_TYPE it = A_Address[k].first; it != A_Address[k].last; ++it) for (k = cbeg; k < cend; ++k) if( localQ[k] == ENUMUNDEF ) printf("%s:%d No column permutation for row %d. Matrix is structurally singular\n",__FILE__,__LINE__,k);
{
i = A_Entries[it].first; for (k = cbeg; k < cend; ++k) localP[k] = k;
//nonzeros for (k = cbeg; k < cend; ++k) V[localQ[k]] = DR[k]; //if you expirience a bug here then the matrix is structurally singular
Ulist[k]++; //row nnz for (k = cbeg; k < cend; ++k) DR[k] = V[k];
Llist[i]++; //col nnz
} for(k = wbeg; k < wend; ++k)
} {
//find node with the lowest order for (INMOST_DATA_ENUM_TYPE jt = A_Address[k].first; jt != A_Address[k].last; ++jt)
INMOST_DATA_ENUM_TYPE start = wbeg; {
INMOST_DATA_ENUM_TYPE index = wbeg; A_Entries[jt].first = localQ[A_Entries[jt].first];
for(k = wbeg; k < wend; ++k) //A_Entries[jt].second *= DL[k]*DR[A_Entries[jt].first];
if( Ulist[k] + Llist[Blist[k]] < Ulist[start] + Llist[Blist[start]] ) }
start = k; std::sort(A_Entries.begin()+A_Address[k].first,A_Entries.begin()+A_Address[k].last);
localP[start] = localQ[Blist[start]] = index++; }
inversePQ(wbeg,wend,localP,localQ, invP,invQ);
applyPQ(wbeg, wend, localP, localQ, invP, invQ);
trcmgraph = Timer();
for (k = wend; k > wbeg; --k)
{
//vwgt[k-1] = 1;
if (A_Address[k-1].Size() > 0)
{
i = A_Entries[A_Address[k-1].first].first;
Blist[k-1] = Bbeg[i];
Bbeg[i] = k-1;
}
Ulist[k-1] = A_Address[k-1].first;
}
xadj[0] = 0;
for(i = wbeg; i < wend; ++i)
{
for (INMOST_DATA_ENUM_TYPE jt = A_Address[i].first; jt != A_Address[i].last; ++jt)
{
if( A_Entries[jt].first != i )
adjncy.push_back(static_cast<INMOST_DATA_ENUM_TYPE>(A_Entries[jt].first));
}
Li = Bbeg[i];
while (Li != EOL)
{
if( Li != i )
adjncy.push_back(static_cast<INMOST_DATA_ENUM_TYPE>(Li));
Li = Blist[Li];
}
Li = Bbeg[i];
while (Li != EOL)
{
Ui = Blist[Li];
Ulist[Li]++;
if (A_Address[Li].last - Ulist[Li] > 0)
{
k = A_Entries[Ulist[Li]].first;
Blist[Li] = Bbeg[k];
Bbeg[k] = Li;
}
Li = Ui;
}
std::sort(adjncy.begin()+xadj[i-wbeg],adjncy.end());
adjncy.resize(std::unique(adjncy.begin()+xadj[i-wbeg],adjncy.end())-adjncy.begin());
xadj[i-wbeg+1] = static_cast<INMOST_DATA_ENUM_TYPE>(adjncy.size());
}
std::fill(Bbeg.begin(),Bbeg.end(),EOL);
trcmgraph = Timer()-trcmgraph;
trcmorder = Timer();
std::fill(Ulist.begin() + wbeg - mobeg, Ulist.begin() + wend - mobeg, ENUMUNDEF);
//find node with the lowest order
INMOST_DATA_ENUM_TYPE start = wbeg;
INMOST_DATA_ENUM_TYPE index = wbeg;
INMOST_DATA_ENUM_TYPE cur = ENUMUNDEF;
std::deque<INMOST_DATA_ENUM_TYPE> q;
std::vector<INMOST_DATA_ENUM_TYPE> conns;
//for(k = wbeg+1; k < wend; ++k)
// if( RCM_Comparator(wbeg,xadj)(k,start) )
// start = k;
//Ulist[start] = index++;
do
{
cur = ENUMUNDEF;
for(k = wbeg; k < wend && cur == ENUMUNDEF; ++k)
{
if( Ulist[k] == ENUMUNDEF )
cur = k;
}
assert(cur != ENUMUNDEF);
for(k = cur+1; k < wend; ++k) if( Ulist[k] == ENUMUNDEF )
{
if( RCM_Comparator(wbeg,xadj)(k,cur) )
cur = k;
}
q.push_back(cur);
Ulist[cur] = index++;
while(!q.empty())
{
cur = q.front();
q.pop_front();
for (INMOST_DATA_ENUM_TYPE it = xadj[cur-wbeg]; it < xadj[cur-wbeg+1]; ++it)
if( Ulist[adjncy[it]] == ENUMUNDEF )
conns.push_back(adjncy[it]);
std::sort(conns.begin(),conns.end(),RCM_Comparator(wbeg,xadj));
for (k = 0; k < static_cast<INMOST_DATA_ENUM_TYPE>(conns.size()); ++k)
{
Ulist[conns[k]] = index++;
q.push_back(conns[k]);
}
conns.clear();
}
}
while( index < wend );
for(k = wbeg; k < wend; ++k)
Ulist[k] = wend-(Ulist[k]-wbeg)-1;
for(k = wbeg; k < wend; ++k)
{
localP[k] = Ulist[k];//wend - Ulist[k] + 1;
localQ[k] = Ulist[k];//wend - Ulist[k] + 1;
//localP[Ulist[k]] = k;
//localQ[Ulist[k]] = k;
}
cend = wend;
i = wend;
trcmorder = Timer() - trcmorder;
}
#else #else
cend = wend; cend = wend;
i = wbeg; i = wbeg;
#endif #endif
tt = Timer() - tt;
treorder += tt;
if (cbeg == cend && cbeg != wend) if (cbeg == cend && cbeg != wend)
{ {
//std::cout << __FILE__ << ":" << __LINE__ << " singular matrix, factored " << mobeg << ".." << cend << " out of " << mobeg << ".." << moend << std::endl; //std::cout << __FILE__ << ":" << __LINE__ << " singular matrix, factored " << mobeg << ".." << cend << " out of " << mobeg << ".." << moend << std::endl;
...@@ -1299,7 +1439,21 @@ public: ...@@ -1299,7 +1439,21 @@ public:
nzA += A_Address[k].Size(); nzA += A_Address[k].Size();
} }
//DumpMatrix(B_Address,B_Entries,cbeg,cend,"mat_reorder.mtx"); /*
{
int rank = 0;
#if defined(USE_MPI)
MPI_Comm_rank(MPI_COMM_WORLD,&rank);
#endif
std::stringstream str;
str << "mat_reorder" << rank << ".mtx";
DumpMatrix(B_Address,B_Entries,cbeg,cend,str.str());
}
#if defined(USE_MPI)
MPI_Barrier(MPI_COMM_WORLD);
#endif
exit(0);
*/
//Setup column addressing for B,F, in descending order to keep list ordered //Setup column addressing for B,F, in descending order to keep list ordered
for (k = cend; k > cbeg; --k) for (k = cend; k > cbeg; --k)
{ {
...@@ -1552,7 +1706,7 @@ public: ...@@ -1552,7 +1706,7 @@ public:
//initialize diagonal values //initialize diagonal values
for (k = cbeg; k < cend; k++) for (k = cbeg; k < cend; k++)
{ {
LU_Diag[k] = 0.0; LU_Diag[k] = 0.0;
localQ[k] = B_Address[k].first; //remember B block starts permutation for pivoting localQ[k] = B_Address[k].first; //remember B block starts permutation for pivoting
for (i = B_Address[k].first; i < B_Address[k].last; i++) for (i = B_Address[k].first; i < B_Address[k].last; i++)
{ {
...@@ -1606,44 +1760,44 @@ public: ...@@ -1606,44 +1760,44 @@ public:
LU_Beg = static_cast<INMOST_DATA_ENUM_TYPE>(LU_Entries.size()); LU_Beg = static_cast<INMOST_DATA_ENUM_TYPE>(LU_Entries.size());
U_Address[cbeg].first = static_cast<INMOST_DATA_ENUM_TYPE>(LU_Entries.size()); U_Address[cbeg].first = static_cast<INMOST_DATA_ENUM_TYPE>(LU_Entries.size());
#if defined(ILUC2) #if defined(ILUC2)
U2_Address[cbeg].first = 0; U2_Address[cbeg].first = 0;
#endif #endif
for (INMOST_DATA_ENUM_TYPE r = B_Address[cbeg].first + (B_Entries[B_Address[cbeg].first].first == cbeg ? 1 : 0); r != B_Address[cbeg].last; ++r) for (INMOST_DATA_ENUM_TYPE r = B_Address[cbeg].first + (B_Entries[B_Address[cbeg].first].first == cbeg ? 1 : 0); r != B_Address[cbeg].last; ++r)
{ {
u = B_Entries[r].second / LU_Diag[cbeg]; u = B_Entries[r].second / LU_Diag[cbeg];
if( fabs(u) > tau ) if( fabs(u) > tau )
LU_Entries.push_back(Sparse::Row::make_entry(B_Entries[r].first, u)); LU_Entries.push_back(Sparse::Row::make_entry(B_Entries[r].first, u));
#if defined(ILUC2) #if defined(ILUC2)
else if( fabs(u) > tau2 ) else if( fabs(u) > tau2 )
LU2_Entries.push_back(Sparse::Row::make_entry(B_Entries[r].first, u)); LU2_Entries.push_back(Sparse::Row::make_entry(B_Entries[r].first, u));
#endif #endif
} }
U_Address[cbeg].last = static_cast<INMOST_DATA_ENUM_TYPE>(LU_Entries.size()); U_Address[cbeg].last = static_cast<INMOST_DATA_ENUM_TYPE>(LU_Entries.size());
#if defined(ILUC2) #if defined(ILUC2)
U2_Address[cbeg].last = static_cast<INMOST_DATA_ENUM_TYPE>(LU2_Entries.size()); U2_Address[cbeg].last = static_cast<INMOST_DATA_ENUM_TYPE>(LU2_Entries.size());
#endif #endif
/////////////////////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////////////////////
// form first line of L and second-order L // // form first line of L and second-order L //
/////////////////////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////////////////////
L_Address[cbeg].first = static_cast<INMOST_DATA_ENUM_TYPE>(LU_Entries.size()); L_Address[cbeg].first = static_cast<INMOST_DATA_ENUM_TYPE>(LU_Entries.size());
#if defined(ILUC2) #if defined(ILUC2)
L2_Address[cbeg].first = static_cast<INMOST_DATA_ENUM_TYPE>(LU2_Entries.size()); L2_Address[cbeg].first = static_cast<INMOST_DATA_ENUM_TYPE>(LU2_Entries.size());
#endif #endif
Li = Bbeg[cbeg]; Li = Bbeg[cbeg];
while (Li != EOL) while (Li != EOL)
{ {
l = B_Entries[B_Address[Li].first].second / LU_Diag[cbeg]; l = B_Entries[B_Address[Li].first].second / LU_Diag[cbeg];
if( fabs(l) > tau ) if( fabs(l) > tau )
LU_Entries.push_back(Sparse::Row::make_entry(Li, l)); LU_Entries.push_back(Sparse::Row::make_entry(Li, l));
#if defined(ILUC2) #if defined(ILUC2)
else if( fabs(l) > tau2 ) else if( fabs(l) > tau2 )
LU2_Entries.push_back(Sparse::Row::make_entry(Li, l)); LU2_Entries.push_back(Sparse::Row::make_entry(Li, l));
#endif #endif
Li = Blist[Li]; Li = Blist[Li];
} }
L_Address[cbeg].last = static_cast<INMOST_DATA_ENUM_TYPE>(LU_Entries.size()); L_Address[cbeg].last = static_cast<INMOST_DATA_ENUM_TYPE>(LU_Entries.size());
#if defined(ILUC2) #if defined(ILUC2)
L2_Address[cbeg].last = static_cast<INMOST_DATA_ENUM_TYPE>(LU2_Entries.size()); L2_Address[cbeg].last = static_cast<INMOST_DATA_ENUM_TYPE>(LU2_Entries.size());
#endif #endif
//update diagonal by factors //update diagonal by factors
/////////////////////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////////////////////
...@@ -1734,31 +1888,31 @@ public: ...@@ -1734,31 +1888,31 @@ public:
// initialize condition estimator // // initialize condition estimator //
/////////////////////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////////////////////
#if defined(ESTIMATOR) #if defined(ESTIMATOR)
if( estimator ) if( estimator )
{ {
NuU = NuL = 1.0; NuU = NuL = 1.0;
if( estimator ) if( estimator )
{ {
EstU1[cbeg] = EstL1[cbeg] = 0.0; EstU1[cbeg] = EstL1[cbeg] = 0.0;
#if defined(ESTIMATOR_REFINE) #if defined(ESTIMATOR_REFINE)
EstU2[cbeg] = EstL2[cbeg] = 0.0; EstU2[cbeg] = EstL2[cbeg] = 0.0;
#endif #endif
for (INMOST_DATA_ENUM_TYPE it = U_Address[cbeg].first; it != U_Address[cbeg].last; ++it) for (INMOST_DATA_ENUM_TYPE it = U_Address[cbeg].first; it != U_Address[cbeg].last; ++it)
{ {
EstU1[LU_Entries[it].first] = LU_Entries[it].second; EstU1[LU_Entries[it].first] = LU_Entries[it].second;
#if defined(ESTIMATOR_REFINE) #if defined(ESTIMATOR_REFINE)
EstU2[LU_Entries[it].first] = LU_Entries[it].second; EstU2[LU_Entries[it].first] = LU_Entries[it].second;
#endif #endif
} }
for (INMOST_DATA_ENUM_TYPE it = L_Address[cbeg].first; it != L_Address[cbeg].last; ++it) for (INMOST_DATA_ENUM_TYPE it = L_Address[cbeg].first; it != L_Address[cbeg].last; ++it)
{ {
EstL1[LU_Entries[it].first] = LU_Entries[it].second; EstL1[LU_Entries[it].first] = LU_Entries[it].second;
#if defined(ESTIMATOR_REFINE) #if defined(ESTIMATOR_REFINE)
EstL2[LU_Entries[it].first] = LU_Entries[it].second; EstL2[LU_Entries[it].first] = LU_Entries[it].second;
#endif #endif
} }
} }
} }
#endif #endif
nzLU += L_Address[cbeg].Size() + U_Address[cbeg].Size() + 1; nzLU += L_Address[cbeg].Size() + U_Address[cbeg].Size() + 1;
max_diag = min_diag = fabs(LU_Diag[cbeg]); max_diag = min_diag = fabs(LU_Diag[cbeg]);
...@@ -1785,7 +1939,7 @@ public: ...@@ -1785,7 +1939,7 @@ public:
// diagonal pivoting algorithm // // diagonal pivoting algorithm //
/////////////////////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////////////////////
#if defined(DIAGONAL_PIVOT_COND) #if defined(DIAGONAL_PIVOT_COND)
int no_swap_algorithm = 1; int no_swap_algorithm = 1;
#endif #endif
#if defined(DIAGONAL_PIVOT_TAU) #if defined(DIAGONAL_PIVOT_TAU)
if ( k < cend-1 && fabs(LU_Diag[k])*(cend-k) < DIAGONAL_PIVOT_TAU*mean_diag) if ( k < cend-1 && fabs(LU_Diag[k])*(cend-k) < DIAGONAL_PIVOT_TAU*mean_diag)
...@@ -1794,19 +1948,19 @@ public: ...@@ -1794,19 +1948,19 @@ public:
#endif #endif
{ {
#if defined(DIAGONAL_PIVOT_COND) #if defined(DIAGONAL_PIVOT_COND)
if( false ) if( false )
{ {
swap_algorithm: swap_algorithm:
no_swap_algorithm--; no_swap_algorithm--;
} }
#endif #endif
j = k+1; j = k+1;
for (i = k + 1; i < cend; i++) if (fabs(LU_Diag[j]) < fabs(LU_Diag[i])) j = i; for (i = k + 1; i < cend; i++) if (fabs(LU_Diag[j]) < fabs(LU_Diag[i])) j = i;
if( j > cend-1 ) if( j > cend-1 )
{ {
std::cout << "Oops! asking to exchange beyond matrix! k = " << k << " j = " << j << " cend = " << cend << std::endl; std::cout << "Oops! asking to exchange beyond matrix! k = " << k << " j = " << j << " cend = " << cend << std::endl;
} }
/* /*
if( fabs(LU_Diag[j]) < fabs(LU_Diag[k]) ) if( fabs(LU_Diag[j]) < fabs(LU_Diag[k]) )
{ {
...@@ -1816,7 +1970,7 @@ swap_algorithm: ...@@ -1816,7 +1970,7 @@ swap_algorithm:
else*/ else*/
{ {
//TODO!!! //TODO!!!
++swaps; ++swaps;
#if defined(REPORT_ILU) #if defined(REPORT_ILU)
//std::cout << "Detected, that there is a much better pivot, i'm " << k << " " << LU_Diag[k] << " other " << j << " " << LU_Diag[j] << std::endl; //std::cout << "Detected, that there is a much better pivot, i'm " << k << " " << LU_Diag[k] << " other " << j << " " << LU_Diag[j] << std::endl;
//std::cout << "Condition numbers: L " << NuL << " D " << NuD << " U " << NuU << std::endl; //std::cout << "Condition numbers: L " << NuL << " D " << NuD << " U " << NuU << std::endl;
...@@ -3032,10 +3186,14 @@ swap_algorithm: ...@@ -3032,10 +3186,14 @@ swap_algorithm:
printf("reorder %f (%6.2f%%)\n", treorder, 100.0*treorder/ttotal); printf("reorder %f (%6.2f%%)\n", treorder, 100.0*treorder/ttotal);
#if defined(REORDER_METIS_ND) #if defined(REORDER_METIS_ND)
printf("metis graph %f nd %f\n", tmetisgraph, tmetisnd); printf("metis graph %f nd %f\n", tmetisgraph, tmetisnd);
#endif
#if defined(REORDER_RCM)
printf("rcm graph %f reorder %f\n", trcmgraph, trcmorder);
#endif #endif
printf("reassamble %f (%6.2f%%)\n", treassamble, 100.0*treassamble / ttotal); printf("reassamble %f (%6.2f%%)\n", treassamble, 100.0*treassamble / ttotal);
printf("rescale %f (%6.2f%%)\n", trescale, 100.0*trescale / ttotal); printf("rescale %f (%6.2f%%)\n", trescale, 100.0*trescale / ttotal);
printf("factor %f (%6.2f%%)\n", tfactor, 100.0*tfactor / ttotal); printf("factor %f (%6.2f%%)\n", tfactor, 100.0*tfactor / ttotal);
printf("nnz A %d LU %d LU2 %d\n",nzA,nzLU,nzLU2);
#endif #endif
init = true; init = true;
/*