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_ConfigDefs.h"
00044 #ifdef HAVE_IFPACK_SPARSKIT
00045 #include "Ifpack_Preconditioner.h"
00046 #include "Ifpack_SPARSKIT.h"
00047 #include "Ifpack_Condest.h"
00048 #include "Ifpack_Utils.h"
00049 #include "Epetra_Comm.h"
00050 #include "Epetra_Map.h"
00051 #include "Epetra_RowMatrix.h"
00052 #include "Epetra_Vector.h"
00053 #include "Epetra_MultiVector.h"
00054 #include "Epetra_Util.h"
00055 #include "Teuchos_ParameterList.hpp"
00056
00057 #define F77_ILUT F77_FUNC(ilut, ILUT)
00058 #define F77_ILUD F77_FUNC(ilud, ILUD)
00059 #define F77_ILUTP F77_FUNC(ilutp, ILUTP)
00060 #define F77_ILUDP F77_FUNC(iludp, ILUDP)
00061 #define F77_ILUK F77_FUNC(iluk, ILUK)
00062 #define F77_ILU0 F77_FUNC(ilu0, ILU0)
00063 #define F77_LUSOL F77_FUNC(lusol, LUSOL)
00064
00065 extern "C" {
00066 void F77_ILUT(int *, double*, int*, int*, int*, double*,
00067 double*, int*, int*, int*, double*, int*, int*);
00068 void F77_ILUD(int *, double*, int*, int*, double*, double*,
00069 double*, int*, int*, int*, double*, int*, int*);
00070 void F77_ILUTP(int *, double*, int*, int*, int*, double*, double*, int*,
00071 double*, int*, int*, int*, double*, int*, int*, int*);
00072 void F77_ILUDP(int *, double*, int*, int*, double*, double*, double*, int*,
00073 double*, int*, int*, int*, double*, int*, int*, int*);
00074 void F77_ILUK(int *, double*, int*, int*, int*,
00075 double*, int*, int*, int*, int*, double*, int*, int*);
00076 void F77_ILU0(int*, double*, int*, int*, double*, int*, int*, int*, int*);
00077 void F77_LUSOL(int *, double*, double*, double*, int*, int*);
00078 }
00079
00080
00081
00082 Ifpack_SPARSKIT::Ifpack_SPARSKIT(Epetra_RowMatrix* A) :
00083 A_(*A),
00084 Comm_(A->Comm()),
00085 UseTranspose_(false),
00086 lfil_(0),
00087 droptol_(0.0),
00088 tol_(0.0),
00089 permtol_(0.1),
00090 alph_(0.0),
00091 mbloc_(-1),
00092 Type_("ILUT"),
00093 Condest_(-1.0),
00094 IsInitialized_(false),
00095 IsComputed_(false),
00096 NumInitialize_(0),
00097 NumCompute_(0),
00098 NumApplyInverse_(0),
00099 InitializeTime_(0.0),
00100 ComputeTime_(0),
00101 ApplyInverseTime_(0),
00102 ComputeFlops_(0.0),
00103 ApplyInverseFlops_(0.0)
00104 {
00105 Teuchos::ParameterList List;
00106 SetParameters(List);
00107 }
00108
00109
00110 Ifpack_SPARSKIT::~Ifpack_SPARSKIT()
00111 {
00112 }
00113
00114
00115 int Ifpack_SPARSKIT::SetParameters(Teuchos::ParameterList& List)
00116 {
00117 lfil_ = List.get("fact: sparskit: lfil", lfil_);
00118 tol_ = List.get("fact: sparskit: tol", tol_);
00119 droptol_ = List.get("fact: sparskit: droptol", droptol_);
00120 permtol_ = List.get("fact: sparskit: permtol", permtol_);
00121 alph_ = List.get("fact: sparskit: alph", alph_);
00122 mbloc_ = List.get("fact: sparskit: mbloc", mbloc_);
00123 Type_ = List.get("fact: sparskit: type", Type_);
00124
00125
00126 Label_ = "IFPACK SPARSKIT (Type=" + Type_ + ", fill=" +
00127 Ifpack_toString(lfil_) + ")";
00128
00129 return(0);
00130 }
00131
00132
00133 int Ifpack_SPARSKIT::Initialize()
00134 {
00135 IsInitialized_ = true;
00136 IsComputed_ = false;
00137 return(0);
00138 }
00139
00140
00141 int Ifpack_SPARSKIT::Compute()
00142 {
00143 if (!IsInitialized())
00144 IFPACK_CHK_ERR(Initialize());
00145
00146 IsComputed_ = false;
00147
00148
00149
00150
00151
00152
00153 int n = Matrix().NumMyRows();
00154 int nnz = Matrix().NumMyNonzeros();
00155
00156 std::vector<double> a(nnz);
00157 std::vector<int> ja(nnz);
00158 std::vector<int> ia(n + 1);
00159
00160 const int MaxNumEntries = Matrix().MaxNumEntries();
00161
00162 std::vector<double> Values(MaxNumEntries);
00163 std::vector<int> Indices(MaxNumEntries);
00164
00165 int count = 0;
00166
00167 ia[0] = 1;
00168 for (int i = 0 ; i < n ; ++i)
00169 {
00170 int NumEntries;
00171 int NumMyEntries = 0;
00172 Matrix().ExtractMyRowCopy(i, MaxNumEntries, NumEntries, &Values[0],
00173 &Indices[0]);
00174
00175
00176
00177
00178 for (int j = 0 ; j < NumEntries ; ++j)
00179 {
00180 if (Indices[j] < n)
00181 {
00182 a[count] = Values[j];
00183 ja[count] = Indices[j] + 1;
00184 ++count;
00185 ++NumMyEntries;
00186 }
00187 }
00188 ia[i + 1] = ia[i] + NumMyEntries;
00189 }
00190
00191 if (mbloc_ == -1) mbloc_ = n;
00192
00193 int iwk;
00194
00195 if (Type_ == "ILUT" || Type_ == "ILUTP" || Type_ == "ILUD" ||
00196 Type_ == "ILUDP")
00197 iwk = nnz + 2 * lfil_ * n;
00198 else if (Type_ == "ILUK")
00199 iwk = (2 * lfil_ + 1) * nnz + 1;
00200 else if (Type_ == "ILU0")
00201 iwk = nnz+2;
00202
00203 int ierr = 0;
00204
00205 alu_.resize(iwk);
00206 jlu_.resize(iwk);
00207 ju_.resize(n + 1);
00208
00209 std::vector<int> jw(n + 1);
00210 std::vector<double> w(n + 1);
00211
00212 if (Type_ == "ILUT")
00213 {
00214 jw.resize(2 * n);
00215 F77_ILUT(&n, &a[0], &ja[0], &ia[0], &lfil_, &droptol_,
00216 &alu_[0], &jlu_[0], &ju_[0], &iwk, &w[0], &jw[0], &ierr);
00217 }
00218 else if (Type_ == "ILUD")
00219 {
00220 jw.resize(2 * n);
00221 F77_ILUD(&n, &a[0], &ja[0], &ia[0], &alph_, &tol_,
00222 &alu_[0], &jlu_[0], &ju_[0], &iwk, &w[0], &jw[0], &ierr);
00223 }
00224 else if (Type_ == "ILUTP")
00225 {
00226 jw.resize(2 * n);
00227 iperm_.resize(2 * n);
00228 F77_ILUTP(&n, &a[0], &ja[0], &ia[0], &lfil_, &droptol_, &permtol_,
00229 &mbloc_, &alu_[0], &jlu_[0], &ju_[0], &iwk, &w[0], &jw[0],
00230 &iperm_[0], &ierr);
00231 for (int i = 0 ; i < n ; ++i)
00232 iperm_[i]--;
00233 }
00234 else if (Type_ == "ILUDP")
00235 {
00236 jw.resize(2 * n);
00237 iperm_.resize(2 * n);
00238 F77_ILUDP(&n, &a[0], &ja[0], &ia[0], &alph_, &droptol_, &permtol_,
00239 &mbloc_, &alu_[0], &jlu_[0], &ju_[0], &n, &w[0], &jw[0],
00240 &iperm_[0], &ierr);
00241 for (int i = 0 ; i < n ; ++i)
00242 iperm_[i]--;
00243 }
00244 else if (Type_ == "ILUK")
00245 {
00246 std::vector<int> levs(iwk);
00247 jw.resize(3 * n);
00248 F77_ILUK(&n, &a[0], &ja[0], &ia[0], &lfil_,
00249 &alu_[0], &jlu_[0], &ju_[0], &levs[0], &iwk, &w[0], &jw[0], &ierr);
00250 }
00251 else if (Type_ == "ILU0")
00252 {
00253
00254 jw.resize(2 * n);
00255 F77_ILU0(&n, &a[0], &ja[0], &ia[0],
00256 &alu_[0], &jlu_[0], &ju_[0], &jw[0], &ierr);
00257 }
00258 IFPACK_CHK_ERR(ierr);
00259
00260 IsComputed_ = true;
00261 return(0);
00262 }
00263
00264
00265
00266 int Ifpack_SPARSKIT::ApplyInverse(const Epetra_MultiVector& X,
00267 Epetra_MultiVector& Y) const
00268 {
00269 if (!IsComputed())
00270 IFPACK_CHK_ERR(-3);
00271
00272 if (X.NumVectors() != Y.NumVectors())
00273 IFPACK_CHK_ERR(-2);
00274
00275 int n = Matrix().NumMyRows();
00276
00277 for (int i = 0 ; i < X.NumVectors() ; ++i)
00278 F77_LUSOL(&n, (double*)X(i)->Values(), Y(i)->Values(), (double*)&alu_[0],
00279 (int*)&jlu_[0], (int*)&ju_[0]);
00280
00281
00282 if (Type_ == "ILUTP" || Type_ == "ILUDP")
00283 {
00284 std::vector<double> tmp(n);
00285 for (int j = 0 ; j < n ; ++j)
00286 tmp[iperm_[j]] = Y[0][j];
00287 for (int j = 0 ; j < n ; ++j)
00288 Y[0][j] = tmp[j];
00289 }
00290
00291 ++NumApplyInverse_;
00292 return(0);
00293
00294 }
00295
00296
00297 double Ifpack_SPARSKIT::Condest(const Ifpack_CondestType CT,
00298 const int MaxIters, const double Tol,
00299 Epetra_RowMatrix* Matrix)
00300 {
00301 if (!IsComputed())
00302 return(-1.0);
00303
00304 if (Condest_ == -1.0)
00305 Condest_ = Ifpack_Condest(*this, CT, MaxIters, Tol, Matrix);
00306
00307 return(Condest_);
00308 }
00309
00310
00311 std::ostream&
00312 Ifpack_SPARSKIT::Print(std::ostream& os) const
00313 {
00314 if (!Comm().MyPID()) {
00315 os << endl;
00316 os << "================================================================================" << endl;
00317 os << "Ifpack_SPARSKIT: " << Label() << endl << endl;
00318 os << "Condition number estimate = " << Condest() << endl;
00319 os << "Global number of rows = " << A_.NumGlobalRows() << endl;
00320 os << endl;
00321 os << "Phase # calls Total Time (s) Total MFlops MFlops/s" << endl;
00322 os << "----- ------- -------------- ------------ --------" << endl;
00323 os << "Initialize() " << std::setw(5) << NumInitialize()
00324 << " " << std::setw(15) << InitializeTime()
00325 << " 0.0 0.0" << endl;
00326 os << "Compute() " << std::setw(5) << NumCompute()
00327 << " " << std::setw(15) << ComputeTime()
00328 << " " << std::setw(15) << 1.0e-6 * ComputeFlops();
00329 if (ComputeTime() != 0.0)
00330 os << " " << std::setw(15) << 1.0e-6 * ComputeFlops() / ComputeTime() << endl;
00331 else
00332 os << " " << std::setw(15) << 0.0 << endl;
00333 os << "ApplyInverse() " << std::setw(5) << NumApplyInverse()
00334 << " " << std::setw(15) << ApplyInverseTime()
00335 << " " << std::setw(15) << 1.0e-6 * ApplyInverseFlops();
00336 if (ApplyInverseTime() != 0.0)
00337 os << " " << std::setw(15) << 1.0e-6 * ApplyInverseFlops() / ApplyInverseTime() << endl;
00338 else
00339 os << " " << std::setw(15) << 0.0 << endl;
00340 os << "================================================================================" << endl;
00341 os << endl;
00342 }
00343
00344
00345 return(os);
00346 }
00347
00348 #endif // HAVE_IFPACK_SPARSKIT