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 #ifndef IFPACK_SPARSECONTAINER_H
00044 #define IFPACK_SPARSECONTAINER_H
00045
00046 #include "Ifpack_Container.h"
00047 #include "Epetra_IntSerialDenseVector.h"
00048 #include "Epetra_MultiVector.h"
00049 #include "Epetra_Vector.h"
00050 #include "Epetra_Map.h"
00051 #include "Epetra_RowMatrix.h"
00052 #include "Epetra_CrsMatrix.h"
00053 #include "Epetra_LinearProblem.h"
00054 #include "Epetra_IntSerialDenseVector.h"
00055 #include "Teuchos_ParameterList.hpp"
00056 #include "Teuchos_RefCountPtr.hpp"
00057 #ifdef HAVE_MPI
00058 #include "Epetra_MpiComm.h"
00059 #else
00060 #include "Epetra_SerialComm.h"
00061 #endif
00062
00092 template<typename T>
00093 class Ifpack_SparseContainer : public Ifpack_Container {
00094
00095 public:
00096
00098
00099 Ifpack_SparseContainer(const int NumRows, const int NumVectors = 1);
00100
00102 Ifpack_SparseContainer(const Ifpack_SparseContainer<T>& rhs);
00103
00105 virtual ~Ifpack_SparseContainer();
00107
00109
00111 Ifpack_SparseContainer& operator=(const Ifpack_SparseContainer<T>& rhs);
00113
00115
00116 virtual int NumRows() const;
00117
00119 virtual int NumVectors() const
00120 {
00121 return(NumVectors_);
00122 }
00123
00125 virtual int SetNumVectors(const int NumVectors_in)
00126 {
00127 if (NumVectors_ != NumVectors_in)
00128 {
00129 NumVectors_=NumVectors_in;
00130 LHS_=Teuchos::rcp(new Epetra_MultiVector(*Map_,NumVectors_));
00131 RHS_=Teuchos::rcp(new Epetra_MultiVector(*Map_,NumVectors_));
00132 }
00133 return(0);
00134 }
00135
00137 virtual double& LHS(const int i, const int Vector = 0);
00138
00140 virtual double& RHS(const int i, const int Vector = 0);
00141
00143
00152 virtual int& ID(const int i);
00153
00155 virtual int SetMatrixElement(const int row, const int col,
00156 const double value);
00157
00158
00160 virtual bool IsInitialized() const
00161 {
00162 return(IsInitialized_);
00163 }
00164
00166 virtual bool IsComputed() const
00167 {
00168 return(IsComputed_);
00169 }
00170
00172 virtual int SetParameters(Teuchos::ParameterList& List);
00173
00175 virtual const char* Label() const
00176 {
00177 return(Label_.c_str());
00178 }
00179
00181 Teuchos::RCP<const Epetra_Map> Map() const
00182 {
00183 return(Map_);
00184 }
00185
00187 Teuchos::RCP<const Epetra_MultiVector> LHS() const
00188 {
00189 return(LHS_);
00190 }
00191
00193 Teuchos::RCP<const Epetra_MultiVector> RHS() const
00194 {
00195 return(RHS_);
00196 }
00197
00199 Teuchos::RCP<const Epetra_CrsMatrix> Matrix() const
00200 {
00201 return(Matrix_);
00202 }
00203
00205 const Epetra_IntSerialDenseVector& ID() const
00206 {
00207 return(GID_);
00208 }
00209
00211 Teuchos::RCP<const T> Inverse() const
00212 {
00213 return(Inverse_);
00214 }
00216
00218
00225 virtual int Initialize();
00227 virtual int Compute(const Epetra_RowMatrix& Matrix_in);
00229 virtual int Apply();
00230
00232 virtual int ApplyInverse();
00233
00235
00237
00238 virtual int Destroy();
00240
00242 virtual double InitializeFlops() const
00243 {
00244 if (Inverse_ == Teuchos::null)
00245 return (0.0);
00246 else
00247 return(Inverse_->InitializeFlops());
00248 }
00249
00251 virtual double ComputeFlops() const
00252 {
00253 if (Inverse_ == Teuchos::null)
00254 return (0.0);
00255 else
00256 return(Inverse_->ComputeFlops());
00257 }
00258
00260 virtual double ApplyFlops() const
00261 {
00262 return(ApplyFlops_);
00263 }
00264
00266 virtual double ApplyInverseFlops() const
00267 {
00268 if (Inverse_ == Teuchos::null)
00269 return (0.0);
00270 else
00271 return(Inverse_->ApplyInverseFlops());
00272 }
00273
00275 virtual ostream& Print(std::ostream& os) const;
00276
00277 private:
00278
00280 virtual int Extract(const Epetra_RowMatrix& Matrix_in);
00281
00283 int NumRows_;
00285 int NumVectors_;
00287 Teuchos::RefCountPtr<Epetra_Map> Map_;
00289 Teuchos::RefCountPtr<Epetra_CrsMatrix> Matrix_;
00291 Teuchos::RefCountPtr<Epetra_MultiVector> LHS_;
00293 Teuchos::RefCountPtr<Epetra_MultiVector> RHS_;
00295 Epetra_IntSerialDenseVector GID_;
00297 bool IsInitialized_;
00299 bool IsComputed_;
00301 Teuchos::RefCountPtr<Epetra_Comm> SerialComm_;
00303 Teuchos::RefCountPtr<T> Inverse_;
00305 string Label_;
00306 Teuchos::ParameterList List_;
00307 double ApplyFlops_;
00308
00309 };
00310
00311
00312 template<typename T>
00313 Ifpack_SparseContainer<T>::
00314 Ifpack_SparseContainer(const int NumRows_in, const int NumVectors_in) :
00315 NumRows_(NumRows_in),
00316 NumVectors_(NumVectors_in),
00317 IsInitialized_(false),
00318 IsComputed_(false),
00319 ApplyFlops_(0.0)
00320 {
00321
00322 #ifdef HAVE_MPI
00323 SerialComm_ = Teuchos::rcp( new Epetra_MpiComm(MPI_COMM_SELF) );
00324 #else
00325 SerialComm_ = Teuchos::rcp( new Epetra_SerialComm );
00326 #endif
00327
00328 }
00329
00330
00331 template<typename T>
00332 Ifpack_SparseContainer<T>::
00333 Ifpack_SparseContainer(const Ifpack_SparseContainer<T>& rhs) :
00334 NumRows_(rhs.NumRows()),
00335 NumVectors_(rhs.NumVectors()),
00336 IsInitialized_(rhs.IsInitialized()),
00337 IsComputed_(rhs.IsComputed())
00338 {
00339
00340 #ifdef HAVE_MPI
00341 SerialComm_ = Teuchos::rcp( new Epetra_MpiComm(MPI_COMM_SELF) );
00342 #else
00343 SerialComm_ = Teuchos::rcp( new Epetra_SerialComm );
00344 #endif
00345
00346 if (!rhs.Map().is_null())
00347 Map_ = Teuchos::rcp( new Epetra_Map(*rhs.Map()) );
00348
00349 if (!rhs.Matrix().is_null())
00350 Matrix_ = Teuchos::rcp( new Epetra_CrsMatrix(*rhs.Matrix()) );
00351
00352 if (!rhs.LHS().is_null())
00353 LHS_ = Teuchos::rcp( new Epetra_MultiVector(*rhs.LHS()) );
00354
00355 if (!rhs.RHS().is_null())
00356 RHS_ = Teuchos::rcp( new Epetra_MultiVector(*rhs.RHS()) );
00357
00358 }
00359
00360 template<typename T>
00361 Ifpack_SparseContainer<T>::~Ifpack_SparseContainer()
00362 {
00363 Destroy();
00364 }
00365
00366
00367 template<typename T>
00368 int Ifpack_SparseContainer<T>::NumRows() const
00369 {
00370 if (IsInitialized() == false)
00371 return(0);
00372 else
00373 return(NumRows_);
00374 }
00375
00376
00377 template<typename T>
00378 int Ifpack_SparseContainer<T>::Initialize()
00379 {
00380
00381 if (IsInitialized_ == true)
00382 Destroy();
00383
00384 IsInitialized_ = false;
00385
00386 #if !defined(EPETRA_NO_32BIT_GLOBAL_INDICES) || !defined(EPETRA_NO_64BIT_GLOBAL_INDICES)
00387 Map_ = Teuchos::rcp( new Epetra_Map(NumRows_,0,*SerialComm_) );
00388 #endif
00389
00390 LHS_ = Teuchos::rcp( new Epetra_MultiVector(*Map_,NumVectors_) );
00391 RHS_ = Teuchos::rcp( new Epetra_MultiVector(*Map_,NumVectors_) );
00392 GID_.Reshape(NumRows_,1);
00393
00394 Matrix_ = Teuchos::rcp( new Epetra_CrsMatrix(Copy,*Map_,0) );
00395
00396
00397 Inverse_ = Teuchos::rcp( new T(Matrix_.get()) );
00398
00399 if (Inverse_ == Teuchos::null)
00400 IFPACK_CHK_ERR(-5);
00401
00402 IFPACK_CHK_ERR(Inverse_->SetParameters(List_));
00403
00404
00405
00406
00407
00408 Label_ = "Ifpack_SparseContainer";
00409
00410 IsInitialized_ = true;
00411 return(0);
00412
00413 }
00414
00415
00416 template<typename T>
00417 double& Ifpack_SparseContainer<T>::LHS(const int i, const int Vector)
00418 {
00419 return(((*LHS_)(Vector))->Values()[i]);
00420 }
00421
00422
00423 template<typename T>
00424 double& Ifpack_SparseContainer<T>::RHS(const int i, const int Vector)
00425 {
00426 return(((*RHS_)(Vector))->Values()[i]);
00427 }
00428
00429
00430 template<typename T>
00431 int Ifpack_SparseContainer<T>::
00432 SetMatrixElement(const int row, const int col, const double value)
00433 {
00434 if (!IsInitialized())
00435 IFPACK_CHK_ERR(-3);
00436
00437 if ((row < 0) || (row >= NumRows())) {
00438 IFPACK_CHK_ERR(-2);
00439 }
00440
00441 if ((col < 0) || (col >= NumRows())) {
00442 IFPACK_CHK_ERR(-2);
00443 }
00444
00445 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
00446 if(Matrix_->RowMatrixRowMap().GlobalIndicesInt()) {
00447 int ierr = Matrix_->InsertGlobalValues((int)row,1,(double*)&value,(int*)&col);
00448 if (ierr < 0) {
00449 ierr = Matrix_->SumIntoGlobalValues((int)row,1,(double*)&value,(int*)&col);
00450 if (ierr < 0)
00451 IFPACK_CHK_ERR(-1);
00452 }
00453 }
00454 else
00455 #endif
00456 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
00457 if(Matrix_->RowMatrixRowMap().GlobalIndicesLongLong()) {
00458 long long col_LL = col;
00459 int ierr = Matrix_->InsertGlobalValues(row,1,(double*)&value,&col_LL);
00460 if (ierr < 0) {
00461 ierr = Matrix_->SumIntoGlobalValues(row,1,(double*)&value,&col_LL);
00462 if (ierr < 0)
00463 IFPACK_CHK_ERR(-1);
00464 }
00465 }
00466 else
00467 #endif
00468 throw "Ifpack_SparseContainer<T>::SetMatrixElement: GlobalIndices type unknown";
00469
00470 return(0);
00471
00472 }
00473
00474
00475 template<typename T>
00476 int Ifpack_SparseContainer<T>::Compute(const Epetra_RowMatrix& Matrix_in)
00477 {
00478
00479 IsComputed_ = false;
00480 if (!IsInitialized()) {
00481 IFPACK_CHK_ERR(Initialize());
00482 }
00483
00484
00485 IFPACK_CHK_ERR(Extract(Matrix_in));
00486
00487
00488 IFPACK_CHK_ERR(Inverse_->Initialize());
00489
00490
00491 IFPACK_CHK_ERR(Inverse_->Compute());
00492
00493 Label_ = "Ifpack_SparseContainer";
00494
00495 IsComputed_ = true;
00496
00497 return(0);
00498 }
00499
00500
00501 template<typename T>
00502 int Ifpack_SparseContainer<T>::Apply()
00503 {
00504 if (IsComputed() == false) {
00505 IFPACK_CHK_ERR(-3);
00506 }
00507
00508 IFPACK_CHK_ERR(Matrix_->Apply(*RHS_, *LHS_));
00509
00510 ApplyFlops_ += 2 * Matrix_->NumGlobalNonzeros64();
00511 return(0);
00512 }
00513
00514
00515 template<typename T>
00516 int Ifpack_SparseContainer<T>::ApplyInverse()
00517 {
00518 if (!IsComputed())
00519 IFPACK_CHK_ERR(-1);
00520
00521 IFPACK_CHK_ERR(Inverse_->ApplyInverse(*RHS_, *LHS_));
00522
00523 return(0);
00524 }
00525
00526
00527
00528 template<typename T>
00529 int Ifpack_SparseContainer<T>::Destroy()
00530 {
00531 IsInitialized_ = false;
00532 IsComputed_ = false;
00533 return(0);
00534 }
00535
00536
00537 template<typename T>
00538 int& Ifpack_SparseContainer<T>::ID(const int i)
00539 {
00540 return(GID_[i]);
00541 }
00542
00543
00544 template<typename T>
00545 int Ifpack_SparseContainer<T>::
00546 SetParameters(Teuchos::ParameterList& List)
00547 {
00548 List_ = List;
00549 return(0);
00550 }
00551
00552
00553
00554 template<typename T>
00555 int Ifpack_SparseContainer<T>::Extract(const Epetra_RowMatrix& Matrix_in)
00556 {
00557
00558 for (int j = 0 ; j < NumRows_ ; ++j) {
00559
00560 if (ID(j) == -1)
00561 IFPACK_CHK_ERR(-1);
00562
00563 if (ID(j) > Matrix_in.NumMyRows())
00564 IFPACK_CHK_ERR(-1);
00565 }
00566
00567 int Length = Matrix_in.MaxNumEntries();
00568 std::vector<double> Values;
00569 Values.resize(Length);
00570 std::vector<int> Indices;
00571 Indices.resize(Length);
00572
00573 for (int j = 0 ; j < NumRows_ ; ++j) {
00574
00575 int LRID = ID(j);
00576
00577 int NumEntries;
00578
00579 int ierr =
00580 Matrix_in.ExtractMyRowCopy(LRID, Length, NumEntries,
00581 &Values[0], &Indices[0]);
00582 IFPACK_CHK_ERR(ierr);
00583
00584 for (int k = 0 ; k < NumEntries ; ++k) {
00585
00586 int LCID = Indices[k];
00587
00588
00589 if (LCID >= Matrix_in.NumMyRows())
00590 continue;
00591
00592
00593
00594
00595 int jj = -1;
00596 for (int kk = 0 ; kk < NumRows_ ; ++kk)
00597 if (ID(kk) == LCID)
00598 jj = kk;
00599
00600 if (jj != -1)
00601 SetMatrixElement(j,jj,Values[k]);
00602
00603 }
00604 }
00605
00606 IFPACK_CHK_ERR(Matrix_->FillComplete());
00607
00608 return(0);
00609 }
00610
00611
00612 template<typename T>
00613 ostream& Ifpack_SparseContainer<T>::Print(ostream & os) const
00614 {
00615 os << "================================================================================" << endl;
00616 os << "Ifpack_SparseContainer" << endl;
00617 os << "Number of rows = " << NumRows() << endl;
00618 os << "Number of vectors = " << NumVectors() << endl;
00619 os << "IsInitialized() = " << IsInitialized() << endl;
00620 os << "IsComputed() = " << IsComputed() << endl;
00621 os << "Flops in Initialize() = " << InitializeFlops() << endl;
00622 os << "Flops in Compute() = " << ComputeFlops() << endl;
00623 os << "Flops in ApplyInverse() = " << ApplyInverseFlops() << endl;
00624 os << "================================================================================" << endl;
00625 os << endl;
00626
00627 return(os);
00628 }
00629 #endif // IFPACK_SPARSECONTAINER_H