00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043 #include "Ifpack_Euclid.h"
00044 #if defined(HAVE_EUCLID) && defined(HAVE_MPI)
00045
00046 #include "Ifpack_Utils.h"
00047 #include <algorithm>
00048 #include "Epetra_MpiComm.h"
00049 #include "Epetra_IntVector.h"
00050 #include "Epetra_Import.h"
00051 #include "Teuchos_ParameterList.hpp"
00052 #include "Teuchos_RCP.hpp"
00053 #include "getRow_dh.h"
00054
00055 using Teuchos::RCP;
00056 using Teuchos::rcp;
00057
00058 Ifpack_Euclid::Ifpack_Euclid(Epetra_CrsMatrix* A):
00059 A_(rcp(A,false)),
00060 UseTranspose_(false),
00061 Condest_(-1),
00062 IsInitialized_(false),
00063 IsComputed_(false),
00064 Label_(),
00065 NumInitialize_(0),
00066 NumCompute_(0),
00067 NumApplyInverse_(0),
00068 InitializeTime_(0.0),
00069 ComputeTime_(0.0),
00070 ApplyInverseTime_(0.0),
00071 ComputeFlops_(0.0),
00072 ApplyInverseFlops_(0.0),
00073 Time_(A_->Comm()),
00074 SetLevel_(1),
00075 SetBJ_(0),
00076 SetStats_(0),
00077 SetMem_(0),
00078 SetSparse_(0.0),
00079 SetRowScale_(0),
00080 SetILUT_(0.0)
00081 {
00082
00083
00084 for(int i = 0; i < A_->NumMyRows(); i++){
00085 int *indices;
00086 int len;
00087 A_->Graph().ExtractMyRowView(i, len, indices);
00088 for(int j = 0; j < len; j++){
00089 indices[j] = A_->GCID(indices[j]);
00090 }
00091 }
00092 }
00093
00094
00095 void Ifpack_Euclid::Destroy(){
00096
00097 if(IsComputed()){
00098 Euclid_dhDestroy(eu);
00099 }
00100
00101 if(IsInitialized()){
00102 Parser_dhDestroy(parser_dh);
00103 parser_dh = NULL;
00104 TimeLog_dhDestroy(tlog_dh);
00105 tlog_dh = NULL;
00106 Mem_dhDestroy(mem_dh);
00107 mem_dh = NULL;
00108 }
00109
00110 for(int i = 0; i < A_->NumMyRows(); i++){
00111 int *indices;
00112 int len;
00113 A_->Graph().ExtractMyRowView(i, len, indices);
00114 for(int j = 0; j < len; j++){
00115 indices[j] = A_->LCID(indices[j]);
00116 }
00117 }
00118 }
00119
00120
00121 int Ifpack_Euclid::Initialize(){
00122
00123 comm_dh = GetMpiComm();
00124 MPI_Comm_size(comm_dh, &np_dh);
00125 MPI_Comm_rank(comm_dh, &myid_dh);
00126 Time_.ResetStartTime();
00127 if(mem_dh == NULL){
00128 Mem_dhCreate(&mem_dh);
00129 }
00130 if (tlog_dh == NULL) {
00131 TimeLog_dhCreate(&tlog_dh);
00132 }
00133
00134 if (parser_dh == NULL) {
00135 Parser_dhCreate(&parser_dh);
00136 }
00137 Parser_dhInit(parser_dh, 0, NULL);
00138
00139 Euclid_dhCreate(&eu);
00140 IsInitialized_=true;
00141 NumInitialize_ = NumInitialize_ + 1;
00142 InitializeTime_ = InitializeTime_ + Time_.ElapsedTime();
00143 return 0;
00144 }
00145
00146
00147 int Ifpack_Euclid::SetParameters(Teuchos::ParameterList& list){
00148 List_ = list;
00149 SetLevel_ = list.get("SetLevel", (int)1);
00150 SetBJ_ = list.get("SetBJ", (int)0);
00151 SetStats_ = list.get("SetStats", (int)0);
00152 SetMem_ = list.get("SetMem", (int)0);
00153 SetSparse_ = list.get("SetSparse", (double)0.0);
00154 SetRowScale_ = list.get("SetRowScale", (int)0);
00155 SetILUT_ = list.get("SetILUT", (double)0.0);
00156 return 0;
00157 }
00158
00159
00160 int Ifpack_Euclid::SetParameter(string name, int value){
00161
00162 std::locale loc;
00163 for(size_t i = 0; i < name.length(); i++){
00164 name[i] = (char) tolower(name[i], loc);
00165 }
00166 if(name.compare("setlevel") == 0){
00167 SetLevel_ = value;
00168 } else if(name.compare("setbj") == 0){
00169 SetBJ_ = value;
00170 } else if(name.compare("setstats") == 0){
00171 SetStats_ = value;
00172 } else if(name.compare("setmem") == 0){
00173 SetMem_ = value;
00174 } else if(name.compare("setrowscale") == 0){
00175 SetRowScale_ = value;
00176 } else {
00177 cout << "\nThe string " << name << " is not an available option.\n";
00178 IFPACK_CHK_ERR(-1);
00179 }
00180 return 0;
00181 }
00182
00183
00184 int Ifpack_Euclid::SetParameter(string name, double value){
00185
00186 std::locale loc;
00187 for(size_t i; i < name.length(); i++){
00188 name[i] = (char) tolower(name[i], loc);
00189 }
00190 if(name.compare("setsparse") == 0){
00191 SetSparse_ = value;
00192 } else if(name.compare("setilut") == 0){
00193 SetILUT_ = value;
00194 } else {
00195 cout << "\nThe string " << name << " is not an available option.\n";
00196 IFPACK_CHK_ERR(-1);
00197 }
00198 return 0;
00199 }
00200
00201
00202 int Ifpack_Euclid::Compute(){
00203 if(IsInitialized() == false){
00204 IFPACK_CHK_ERR(Initialize());
00205 }
00206 Time_.ResetStartTime();
00207 sprintf(Label_, "IFPACK_Euclid (level=%d, bj=%d, stats=%d, mem=%d, sparse=%f, rowscale=%d, ilut=%f)",
00208 SetLevel_, SetBJ_, SetStats_, SetMem_, SetSparse_, SetRowScale_, SetILUT_);
00209
00210 eu->level = SetLevel_;
00211 if(SetBJ_ != 0){
00212 strcpy("bj", eu->algo_par);
00213 }
00214 if(SetSparse_ != 0.0){
00215 eu->sparseTolA = SetSparse_;
00216 }
00217 if(SetRowScale_ != 0){
00218 eu->isScaled = true;
00219 }
00220 if(SetILUT_ != 0.0){
00221 eu->droptol = SetILUT_;
00222 }
00223 if(SetStats_ != 0 || SetMem_ != 0){
00224 eu->logging = true;
00225 Parser_dhInsert(parser_dh, "-eu_stats", "1");
00226 }
00227
00228 eu->A = (void*) A_.get();
00229 eu->m = A_->NumMyRows();
00230 eu->n = A_->NumGlobalRows();
00231 Euclid_dhSetup(eu);
00232 IsComputed_ = true;
00233 NumCompute_ = NumCompute_ + 1;
00234 ComputeTime_ = ComputeTime_ + Time_.ElapsedTime();
00235 return 0;
00236 }
00237
00238
00239 int Ifpack_Euclid::ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const{
00240 if(IsComputed() == false){
00241 IFPACK_CHK_ERR(-1);
00242 }
00243 int NumVectors = X.NumVectors();
00244 if(NumVectors != Y.NumVectors()){
00245 IFPACK_CHK_ERR(-2);
00246 }
00247 Time_.ResetStartTime();
00248
00249 for(int vecNum = 0; vecNum < NumVectors; vecNum++){
00250 CallEuclid(X[vecNum], Y[vecNum]);
00251 }
00252 if(SetStats_ != 0){
00253 Euclid_dhPrintTestData(eu, stdout);
00254 }
00255 NumApplyInverse_ = NumApplyInverse_ + 1;
00256 ApplyInverseTime_ = ApplyInverseTime_ + Time_.ElapsedTime();
00257 return 0;
00258 }
00259
00260
00261 int Ifpack_Euclid::CallEuclid(double *x, double *y) const{
00262 Euclid_dhApply(eu, x, y);
00263 return 0;
00264 }
00265
00266
00267 ostream& operator << (ostream& os, const Ifpack_Euclid& A){
00268 if (!A.Comm().MyPID()) {
00269 os << endl;
00270 os << "================================================================================" << endl;
00271 os << "Ifpack_Euclid: " << A.Label () << endl << endl;
00272 os << "Using " << A.Comm().NumProc() << " processors." << endl;
00273 os << "Global number of rows = " << A.Matrix().NumGlobalRows() << endl;
00274 os << "Global number of nonzeros = " << A.Matrix().NumGlobalNonzeros() << endl;
00275 os << "Condition number estimate = " << A.Condest() << endl;
00276 os << endl;
00277 os << "Phase # calls Total Time (s) Total MFlops MFlops/s" << endl;
00278 os << "----- ------- -------------- ------------ --------" << endl;
00279 os << "Initialize() " << std::setw(5) << A.NumInitialize()
00280 << " " << std::setw(15) << A.InitializeTime()
00281 << " 0.0 0.0" << endl;
00282 os << "Compute() " << std::setw(5) << A.NumCompute()
00283 << " " << std::setw(15) << A.ComputeTime()
00284 << " " << std::setw(15) << 1.0e-6 * A.ComputeFlops();
00285 if (A.ComputeTime() != 0.0)
00286 os << " " << std::setw(15) << 1.0e-6 * A.ComputeFlops() / A.ComputeTime() << endl;
00287 else
00288 os << " " << std::setw(15) << 0.0 << endl;
00289 os << "ApplyInverse() " << std::setw(5) << A.NumApplyInverse()
00290 << " " << std::setw(15) << A.ApplyInverseTime()
00291 << " " << std::setw(15) << 1.0e-6 * A.ApplyInverseFlops();
00292 if (A.ApplyInverseTime() != 0.0)
00293 os << " " << std::setw(15) << 1.0e-6 * A.ApplyInverseFlops() / A.ApplyInverseTime() << endl;
00294 else
00295 os << " " << std::setw(15) << 0.0 << endl;
00296 os << "================================================================================" << endl;
00297 os << endl;
00298 }
00299 return os;
00300 }
00301
00302
00303 double Ifpack_Euclid::Condest(const Ifpack_CondestType CT,
00304 const int MaxIters,
00305 const double Tol,
00306 Epetra_RowMatrix* Matrix_in){
00307 if (!IsComputed())
00308 return(-1.0);
00309 return(Condest_);
00310 }
00311
00312 #endif // HAVE_EUCLID && HAVE_MPI