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
00044 #ifndef IFPACK_BLOCKPRECONDITIONER_H
00045 #define IFPACK_BLOCKPRECONDITIONER_H
00046
00047 #include "Ifpack_ConfigDefs.h"
00048 #include "Ifpack_Preconditioner.h"
00049 #include "Ifpack_Partitioner.h"
00050 #include "Ifpack_LinePartitioner.h"
00051 #include "Ifpack_LinearPartitioner.h"
00052 #include "Ifpack_GreedyPartitioner.h"
00053 #include "Ifpack_METISPartitioner.h"
00054 #include "Ifpack_EquationPartitioner.h"
00055 #include "Ifpack_UserPartitioner.h"
00056 #include "Ifpack_Graph_Epetra_RowMatrix.h"
00057 #include "Ifpack_DenseContainer.h"
00058 #include "Ifpack_Utils.h"
00059 #include "Teuchos_ParameterList.hpp"
00060 #include "Teuchos_RefCountPtr.hpp"
00061 #include "Epetra_RowMatrix.h"
00062 #include "Epetra_MultiVector.h"
00063 #include "Epetra_Vector.h"
00064 #include "Epetra_Time.h"
00065 #include "Epetra_Import.h"
00066
00067 static const int IFPACK_JACOBI = 0;
00068 static const int IFPACK_GS = 1;
00069 static const int IFPACK_SGS = 2;
00070
00071
00073
00136 template<typename T>
00137 class Ifpack_BlockRelaxation : public Ifpack_Preconditioner {
00138
00139 public:
00140
00142
00143
00148 Ifpack_BlockRelaxation(const Epetra_RowMatrix* Matrix);
00149
00150 virtual ~Ifpack_BlockRelaxation();
00151
00153
00155
00157
00165 virtual int Apply(const Epetra_MultiVector& X,
00166 Epetra_MultiVector& Y) const;
00167
00169
00178 virtual int ApplyInverse(const Epetra_MultiVector& X,
00179 Epetra_MultiVector& Y) const;
00180
00182 virtual double NormInf() const
00183 {
00184 return(-1.0);
00185 }
00187
00189
00190 virtual int SetUseTranspose(bool UseTranspose_in)
00191 {
00192 if (UseTranspose_in)
00193 IFPACK_CHK_ERR(-98);
00194 return(0);
00195 }
00196
00197 virtual const char* Label() const;
00198
00200 virtual bool UseTranspose() const
00201 {
00202 return(false);
00203 }
00204
00206 virtual bool HasNormInf() const
00207 {
00208 return(false);
00209 }
00210
00212 virtual const Epetra_Comm & Comm() const;
00213
00215 virtual const Epetra_Map & OperatorDomainMap() const;
00216
00218 virtual const Epetra_Map & OperatorRangeMap() const;
00220
00222 int NumLocalBlocks() const
00223 {
00224 return(NumLocalBlocks_);
00225 }
00226
00228 virtual bool IsInitialized() const
00229 {
00230 return(IsInitialized_);
00231 }
00232
00234 virtual bool IsComputed() const
00235 {
00236 return(IsComputed_);
00237 }
00238
00240 virtual int SetParameters(Teuchos::ParameterList& List);
00241
00243 virtual int Initialize();
00244
00246 virtual int Compute();
00247
00248 virtual const Epetra_RowMatrix& Matrix() const
00249 {
00250 return(*Matrix_);
00251 }
00252
00253 virtual double Condest(const Ifpack_CondestType CT = Ifpack_Cheap,
00254 const int MaxIters = 1550,
00255 const double Tol = 1e-9,
00256 Epetra_RowMatrix* Matrix_in = 0)
00257 {
00258 return(-1.0);
00259 }
00260
00261 virtual double Condest() const
00262 {
00263 return(-1.0);
00264 }
00265
00266 std::ostream& Print(std::ostream& os) const;
00267
00269 virtual int NumInitialize() const
00270 {
00271 return(NumInitialize_);
00272 }
00273
00275 virtual int NumCompute() const
00276 {
00277 return(NumCompute_);
00278 }
00279
00281 virtual int NumApplyInverse() const
00282 {
00283 return(NumApplyInverse_);
00284 }
00285
00287 virtual double InitializeTime() const
00288 {
00289 return(InitializeTime_);
00290 }
00291
00293 virtual double ComputeTime() const
00294 {
00295 return(ComputeTime_);
00296 }
00297
00299 virtual double ApplyInverseTime() const
00300 {
00301 return(ApplyInverseTime_);
00302 }
00303
00305 virtual double InitializeFlops() const
00306 {
00307 #ifdef IFPACK_FLOPCOUNTERS
00308 if (Containers_.size() == 0)
00309 return(0.0);
00310
00311
00312
00313
00314 double total = InitializeFlops_;
00315 for (unsigned int i = 0 ; i < Containers_.size() ; ++i)
00316 if(Containers_[i]) total += Containers_[i]->InitializeFlops();
00317 return(total);
00318 #else
00319 return(0.0);
00320 #endif
00321 }
00322
00323 virtual double ComputeFlops() const
00324 {
00325 #ifdef IFPACK_FLOPCOUNTERS
00326 if (Containers_.size() == 0)
00327 return(0.0);
00328
00329 double total = ComputeFlops_;
00330 for (unsigned int i = 0 ; i < Containers_.size() ; ++i)
00331 if(Containers_[i]) total += Containers_[i]->ComputeFlops();
00332 return(total);
00333 #else
00334 return(0.0);
00335 #endif
00336 }
00337
00338 virtual double ApplyInverseFlops() const
00339 {
00340 #ifdef IFPACK_FLOPCOUNTERS
00341 if (Containers_.size() == 0)
00342 return(0.0);
00343
00344 double total = ApplyInverseFlops_;
00345 for (unsigned int i = 0 ; i < Containers_.size() ; ++i) {
00346 if(Containers_[i]) total += Containers_[i]->ApplyInverseFlops();
00347 }
00348 return(total);
00349 #else
00350 return(0.0);
00351 #endif
00352 }
00353
00354 private:
00355
00357 Ifpack_BlockRelaxation(const Ifpack_BlockRelaxation& rhs);
00358
00360 Ifpack_BlockRelaxation & operator=(const Ifpack_BlockRelaxation& rhs)
00361 {
00362 return(*this);
00363 }
00364
00365 virtual int ApplyInverseJacobi(const Epetra_MultiVector& X,
00366 Epetra_MultiVector& Y) const;
00367
00368 virtual int DoJacobi(const Epetra_MultiVector& X,
00369 Epetra_MultiVector& Y) const;
00370
00371 virtual int ApplyInverseGS(const Epetra_MultiVector& X,
00372 Epetra_MultiVector& Y) const;
00373
00374 virtual int DoGaussSeidel(Epetra_MultiVector& X,
00375 Epetra_MultiVector& Y) const;
00376
00377 virtual int ApplyInverseSGS(const Epetra_MultiVector& X,
00378 Epetra_MultiVector& Y) const;
00379
00380 virtual int DoSGS(const Epetra_MultiVector& X,
00381 Epetra_MultiVector& Xtmp,
00382 Epetra_MultiVector& Y) const;
00383
00384 int ExtractSubmatrices();
00385
00386
00387
00389 bool IsInitialized_;
00391 bool IsComputed_;
00393 int NumInitialize_;
00395 int NumCompute_;
00397 mutable int NumApplyInverse_;
00399 double InitializeTime_;
00401 double ComputeTime_;
00403 mutable double ApplyInverseTime_;
00405 double InitializeFlops_;
00407 double ComputeFlops_;
00409 mutable double ApplyInverseFlops_;
00410
00411
00412
00414 int NumSweeps_;
00416 double DampingFactor_;
00418 int NumLocalBlocks_;
00420 Teuchos::ParameterList List_;
00421
00422
00423
00426 Teuchos::RefCountPtr< const Epetra_RowMatrix > Matrix_;
00427 mutable std::vector<Teuchos::RefCountPtr<T> > Containers_;
00428 Epetra_Vector Diagonal_ ;
00429
00431 Teuchos::RefCountPtr<Ifpack_Partitioner> Partitioner_;
00432 string PartitionerType_;
00433 int PrecType_;
00435 string Label_;
00437 bool ZeroStartingSolution_;
00438 Teuchos::RefCountPtr<Ifpack_Graph> Graph_;
00440 Teuchos::RefCountPtr<Epetra_Vector> W_;
00441
00442 int OverlapLevel_;
00443 mutable Epetra_Time Time_;
00444 bool IsParallel_;
00445 Teuchos::RefCountPtr<Epetra_Import> Importer_;
00446
00447
00448 };
00449
00450
00451 template<typename T>
00452 Ifpack_BlockRelaxation<T>::
00453 Ifpack_BlockRelaxation(const Epetra_RowMatrix* Matrix_in) :
00454 IsInitialized_(false),
00455 IsComputed_(false),
00456 NumInitialize_(0),
00457 NumCompute_(0),
00458 NumApplyInverse_(0),
00459 InitializeTime_(0.0),
00460 ComputeTime_(0.0),
00461 ApplyInverseTime_(0.0),
00462 InitializeFlops_(0.0),
00463 ComputeFlops_(0.0),
00464 ApplyInverseFlops_(0.0),
00465 NumSweeps_(1),
00466 DampingFactor_(1.0),
00467 NumLocalBlocks_(1),
00468 Matrix_(Teuchos::rcp(Matrix_in,false)),
00469 Diagonal_( Matrix_in->Map()),
00470 PartitionerType_("greedy"),
00471 PrecType_(IFPACK_JACOBI),
00472 ZeroStartingSolution_(true),
00473 OverlapLevel_(0),
00474 Time_(Comm()),
00475 IsParallel_(false)
00476 {
00477 if (Matrix_in->Comm().NumProc() != 1)
00478 IsParallel_ = true;
00479 }
00480
00481
00482 template<typename T>
00483 Ifpack_BlockRelaxation<T>::~Ifpack_BlockRelaxation()
00484 {
00485 }
00486
00487
00488 template<typename T>
00489 const char* Ifpack_BlockRelaxation<T>::Label() const
00490 {
00491 return(Label_.c_str());
00492 }
00493
00494
00495 template<typename T>
00496 int Ifpack_BlockRelaxation<T>::
00497 Apply(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const
00498 {
00499 int ierr = Matrix().Apply(X,Y);
00500 IFPACK_RETURN(ierr);
00501 }
00502
00503
00504 template<typename T>
00505 const Epetra_Comm& Ifpack_BlockRelaxation<T>::
00506 Comm() const
00507 {
00508 return(Matrix().Comm());
00509 }
00510
00511
00512 template<typename T>
00513 const Epetra_Map& Ifpack_BlockRelaxation<T>::
00514 OperatorDomainMap() const
00515 {
00516 return(Matrix().OperatorDomainMap());
00517 }
00518
00519
00520 template<typename T>
00521 const Epetra_Map& Ifpack_BlockRelaxation<T>::
00522 OperatorRangeMap() const
00523 {
00524 return(Matrix().OperatorRangeMap());
00525 }
00526
00527
00528 template<typename T>
00529 int Ifpack_BlockRelaxation<T>::ExtractSubmatrices()
00530 {
00531
00532 if (Partitioner_ == Teuchos::null)
00533 IFPACK_CHK_ERR(-3);
00534
00535 NumLocalBlocks_ = Partitioner_->NumLocalParts();
00536
00537 Containers_.resize(NumLocalBlocks());
00538
00539 Diagonal_ = Epetra_Vector(Matrix_->Map());
00540 Matrix_->ExtractDiagonalCopy(Diagonal_);
00541
00542 for (int i = 0 ; i < NumLocalBlocks() ; ++i) {
00543
00544 int rows = Partitioner_->NumRowsInPart(i);
00545
00546
00547
00548 if( rows != 1 ) {
00549 Containers_[i] = Teuchos::rcp( new T(rows) );
00550
00551 IFPACK_CHK_ERR(Containers_[i]->SetParameters(List_));
00552 IFPACK_CHK_ERR(Containers_[i]->Initialize());
00553
00554
00555
00556 for (int j = 0 ; j < rows ; ++j) {
00557 int LRID = (*Partitioner_)(i,j);
00558 Containers_[i]->ID(j) = LRID;
00559 }
00560
00561 IFPACK_CHK_ERR(Containers_[i]->Compute(*Matrix_));
00562 }
00563
00564 }
00565
00566 int issing = 0;
00567
00568 for (int i = 0 ; i < NumLocalBlocks() ; ++i)
00569 issing += (int) ( Partitioner_->NumRowsInPart(i) == 1);
00570 printf( " %d of %d containers are singleton \n",issing,NumLocalBlocks());
00571
00572 return(0);
00573 }
00574
00575
00576 template<typename T>
00577 int Ifpack_BlockRelaxation<T>::Compute()
00578 {
00579
00580 if (!IsInitialized())
00581 IFPACK_CHK_ERR(Initialize());
00582
00583 Time_.ResetStartTime();
00584
00585 IsComputed_ = false;
00586
00587 if (Matrix().NumGlobalRows64() != Matrix().NumGlobalCols64())
00588 IFPACK_CHK_ERR(-2);
00589
00590 IFPACK_CHK_ERR(ExtractSubmatrices());
00591
00592 if (IsParallel_ && PrecType_ != IFPACK_JACOBI) {
00593
00594 Importer_ = Teuchos::rcp( new Epetra_Import(Matrix().RowMatrixColMap(),
00595 Matrix().RowMatrixRowMap()) );
00596
00597 if (Importer_ == Teuchos::null) IFPACK_CHK_ERR(-5);
00598 }
00599 IsComputed_ = true;
00600 ComputeTime_ += Time_.ElapsedTime();
00601 ++NumCompute_;
00602
00603 return(0);
00604
00605 }
00606
00607
00608 template<typename T>
00609 int Ifpack_BlockRelaxation<T>::
00610 ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const
00611 {
00612 if (!IsComputed())
00613 IFPACK_CHK_ERR(-3);
00614
00615 if (X.NumVectors() != Y.NumVectors())
00616 IFPACK_CHK_ERR(-2);
00617
00618 Time_.ResetStartTime();
00619
00620
00621
00622 Teuchos::RefCountPtr< const Epetra_MultiVector > Xcopy;
00623 if (X.Pointers()[0] == Y.Pointers()[0])
00624 Xcopy = Teuchos::rcp( new Epetra_MultiVector(X) );
00625 else
00626 Xcopy = Teuchos::rcp( &X, false );
00627
00628 switch (PrecType_) {
00629 case IFPACK_JACOBI:
00630 IFPACK_CHK_ERR(ApplyInverseJacobi(*Xcopy,Y));
00631 break;
00632 case IFPACK_GS:
00633 IFPACK_CHK_ERR(ApplyInverseGS(*Xcopy,Y));
00634 break;
00635 case IFPACK_SGS:
00636 IFPACK_CHK_ERR(ApplyInverseSGS(*Xcopy,Y));
00637 break;
00638 }
00639
00640 ApplyInverseTime_ += Time_.ElapsedTime();
00641 ++NumApplyInverse_;
00642
00643 return(0);
00644 }
00645
00646
00647
00648
00649
00650 template<typename T>
00651 int Ifpack_BlockRelaxation<T>::
00652 ApplyInverseJacobi(const Epetra_MultiVector& X,
00653 Epetra_MultiVector& Y) const
00654 {
00655
00656 if (ZeroStartingSolution_)
00657 Y.PutScalar(0.0);
00658
00659
00660 if (NumSweeps_ == 1 && ZeroStartingSolution_) {
00661 int ierr = DoJacobi(X,Y);
00662 IFPACK_RETURN(ierr);
00663 }
00664
00665 Epetra_MultiVector AX(Y);
00666
00667 for (int j = 0; j < NumSweeps_ ; j++) {
00668 IFPACK_CHK_ERR(Apply(Y,AX));
00669 ApplyInverseFlops_ += X.NumVectors() * 2 * Matrix_->NumGlobalNonzeros64();
00670 IFPACK_CHK_ERR(AX.Update(1.0,X,-1.0));
00671 ApplyInverseFlops_ += X.NumVectors() * 2 * Matrix_->NumGlobalRows64();
00672 IFPACK_CHK_ERR(DoJacobi(AX,Y));
00673
00674 }
00675
00676
00677 return(0);
00678 }
00679
00680
00681 template<typename T>
00682 int Ifpack_BlockRelaxation<T>::
00683 DoJacobi(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const
00684 {
00685 int NumVectors = X.NumVectors();
00686
00687 if (OverlapLevel_ == 0) {
00688
00689 for (int i = 0 ; i < NumLocalBlocks() ; ++i) {
00690
00691 int rows = Partitioner_->NumRowsInPart(i);
00692
00693 if (rows == 0)
00694 continue;
00695
00696 if(rows != 1) {
00697 int LID;
00698
00699
00700 for (int j = 0 ; j < Partitioner_->NumRowsInPart(i) ; ++j) {
00701 LID = Containers_[i]->ID(j);
00702 for (int k = 0 ; k < NumVectors ; ++k) {
00703 Containers_[i]->RHS(j,k) = X[k][LID];
00704 }
00705 }
00706
00707
00708
00709
00710
00711 IFPACK_CHK_ERR(Containers_[i]->ApplyInverse());
00712
00713
00714 for (int j = 0 ; j < Partitioner_->NumRowsInPart(i) ; ++j) {
00715 LID = Containers_[i]->ID(j);
00716 for (int k = 0 ; k < NumVectors ; ++k) {
00717 Y[k][LID] += DampingFactor_ * Containers_[i]->LHS(j,k);
00718 }
00719 }
00720 }
00721 else {
00722
00723 int LRID = (*Partitioner_)(i,0);
00724 double b = X[0][LRID];
00725 double a = Diagonal_[LRID];
00726 Y[0][LRID] += DampingFactor_* b/a;
00727 }
00728
00729
00730 #ifdef IFPACK_FLOPCOUNTERS
00731 ApplyInverseFlops_ += NumVectors * 2 * Matrix_->NumGlobalRows();
00732 #endif
00733
00734 }
00735 }
00736 else {
00737
00738 for (int i = 0 ; i < NumLocalBlocks() ; ++i) {
00739
00740 int rows = Partitioner_->NumRowsInPart(i);
00741
00742
00743 if (rows == 0)
00744 continue;
00745 if(rows != 1) {
00746 int LID;
00747
00748
00749 for (int j = 0 ; j < Partitioner_->NumRowsInPart(i) ; ++j) {
00750 LID = Containers_[i]->ID(j);
00751 for (int k = 0 ; k < NumVectors ; ++k) {
00752 Containers_[i]->RHS(j,k) = (*W_)[LID] * X[k][LID];
00753 }
00754 }
00755
00756
00757
00758 IFPACK_CHK_ERR(Containers_[i]->ApplyInverse());
00759
00760
00761 for (int j = 0 ; j < Partitioner_->NumRowsInPart(i) ; ++j) {
00762 LID = Containers_[i]->ID(j);
00763 for (int k = 0 ; k < NumVectors ; ++k) {
00764 Y[k][LID] += DampingFactor_ * (*W_)[LID] * Containers_[i]->LHS(j,k);
00765 }
00766 }
00767 }
00768 else {
00769 int LRID = (*Partitioner_)(i,0);
00770 double w = (*W_)[LRID];
00771 double b = w * X[0][LRID];
00772 double a = Diagonal_[LRID];
00773
00774 Y[0][LRID] += DampingFactor_ * w * b / a;
00775 }
00776 }
00777
00778
00779
00780 #ifdef IFPACK_FLOPCOUNTERS
00781 ApplyInverseFlops_ += NumVectors * 4 * Matrix_->NumGlobalRows();
00782 #endif
00783 }
00784
00785 return(0);
00786 }
00787
00788
00789 template<typename T>
00790 int Ifpack_BlockRelaxation<T>::
00791 ApplyInverseGS(const Epetra_MultiVector& X,
00792 Epetra_MultiVector& Y) const
00793 {
00794
00795 if (ZeroStartingSolution_)
00796 Y.PutScalar(0.0);
00797
00798 Epetra_MultiVector Xcopy(X);
00799 for (int j = 0; j < NumSweeps_ ; j++) {
00800 IFPACK_CHK_ERR(DoGaussSeidel(Xcopy,Y));
00801 if (j != NumSweeps_ - 1)
00802 Xcopy = X;
00803 }
00804
00805 return(0);
00806
00807 }
00808
00809
00810 template<typename T>
00811 int Ifpack_BlockRelaxation<T>::
00812 DoGaussSeidel(Epetra_MultiVector& X, Epetra_MultiVector& Y) const
00813 {
00814
00815
00816
00817 int Length = Matrix().MaxNumEntries();
00818 std::vector<int> Indices(Length);
00819 std::vector<double> Values(Length);
00820
00821 int NumMyRows = Matrix().NumMyRows();
00822 int NumVectors = X.NumVectors();
00823
00824
00825
00826
00827 Teuchos::RefCountPtr< Epetra_MultiVector > Y2;
00828 if (IsParallel_)
00829 Y2 = Teuchos::rcp( new Epetra_MultiVector(Importer_->TargetMap(), NumVectors) );
00830 else
00831 Y2 = Teuchos::rcp( &Y, false );
00832
00833 double** y_ptr;
00834 double** y2_ptr;
00835 Y.ExtractView(&y_ptr);
00836 Y2->ExtractView(&y2_ptr);
00837
00838
00839 if (IsParallel_)
00840 IFPACK_CHK_ERR(Y2->Import(Y,*Importer_,Insert));
00841
00842 for (int i = 0 ; i < NumLocalBlocks() ; ++i) {
00843 int rows = Partitioner_->NumRowsInPart(i);
00844
00845
00846 if (rows!=1 && Containers_[i]->NumRows() == 0)
00847 continue;
00848
00849 int LID;
00850
00851
00852
00853 for (int j = 0 ; j < Partitioner_->NumRowsInPart(i) ; ++j) {
00854 LID = (*Partitioner_)(i,j);
00855
00856 int NumEntries;
00857 IFPACK_CHK_ERR(Matrix().ExtractMyRowCopy(LID, Length,NumEntries,
00858 &Values[0], &Indices[0]));
00859
00860 for (int k = 0 ; k < NumEntries ; ++k) {
00861 int col = Indices[k];
00862
00863 for (int kk = 0 ; kk < NumVectors ; ++kk) {
00864 X[kk][LID] -= Values[k] * y2_ptr[kk][col];
00865 }
00866 }
00867 }
00868
00869 if(rows != 1) {
00870
00871
00872 for (int j = 0 ; j < Partitioner_->NumRowsInPart(i) ; ++j) {
00873 LID = Containers_[i]->ID(j);
00874 for (int k = 0 ; k < NumVectors ; ++k) {
00875 Containers_[i]->RHS(j,k) = X[k][LID];
00876 }
00877 }
00878
00879 IFPACK_CHK_ERR(Containers_[i]->ApplyInverse());
00880 #ifdef IFPACK_FLOPCOUNTERS
00881 ApplyInverseFlops_ += Containers_[i]->ApplyInverseFlops();
00882 #endif
00883
00884 for (int j = 0 ; j < Partitioner_->NumRowsInPart(i) ; ++j) {
00885 LID = Containers_[i]->ID(j);
00886 for (int k = 0 ; k < NumVectors ; ++k) {
00887 double temp = DampingFactor_ * Containers_[i]->LHS(j,k);
00888 y2_ptr[k][LID] += temp;
00889 }
00890 }
00891 }
00892 else {
00893 int LRID = (*Partitioner_)(i,0);
00894 double b = X[0][LRID];
00895 double a = Diagonal_[LRID];
00896 y2_ptr[0][LRID]+= DampingFactor_* b/a;
00897 }
00898 }
00899
00900
00901
00902 #ifdef IFPACK_FLOPCOUNTERS
00903 ApplyInverseFlops_ += NumVectors * 2 * Matrix_->NumGlobalNonzeros();
00904 ApplyInverseFlops_ += NumVectors * 2 * Matrix_->NumGlobalRows();
00905 #endif
00906
00907
00908
00909 if (IsParallel_)
00910 for (int m = 0 ; m < NumVectors ; ++m)
00911 for (int i = 0 ; i < NumMyRows ; ++i)
00912 y_ptr[m][i] = y2_ptr[m][i];
00913
00914 return(0);
00915 }
00916
00917
00918 template<typename T>
00919 int Ifpack_BlockRelaxation<T>::
00920 ApplyInverseSGS(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const
00921 {
00922
00923 if (ZeroStartingSolution_)
00924 Y.PutScalar(0.0);
00925
00926 Epetra_MultiVector Xcopy(X);
00927 for (int j = 0; j < NumSweeps_ ; j++) {
00928 IFPACK_CHK_ERR(DoSGS(X,Xcopy,Y));
00929 if (j != NumSweeps_ - 1)
00930 Xcopy = X;
00931 }
00932 return(0);
00933 }
00934
00935
00936 template<typename T>
00937 int Ifpack_BlockRelaxation<T>::
00938 DoSGS(const Epetra_MultiVector& X, Epetra_MultiVector& Xcopy,
00939 Epetra_MultiVector& Y) const
00940 {
00941
00942 int NumMyRows = Matrix().NumMyRows();
00943 int NumVectors = X.NumVectors();
00944
00945 int Length = Matrix().MaxNumEntries();
00946 std::vector<int> Indices;
00947 std::vector<double> Values;
00948 Indices.resize(Length);
00949 Values.resize(Length);
00950
00951
00952
00953
00954 Teuchos::RefCountPtr< Epetra_MultiVector > Y2;
00955 if (IsParallel_)
00956 Y2 = Teuchos::rcp( new Epetra_MultiVector(Importer_->TargetMap(), NumVectors) );
00957 else
00958 Y2 = Teuchos::rcp( &Y, false );
00959
00960 double** y_ptr;
00961 double** y2_ptr;
00962 Y.ExtractView(&y_ptr);
00963 Y2->ExtractView(&y2_ptr);
00964
00965
00966 if (IsParallel_)
00967 IFPACK_CHK_ERR(Y2->Import(Y,*Importer_,Insert));
00968
00969 for (int i = 0 ; i < NumLocalBlocks() ; ++i) {
00970 int rows = Partitioner_->NumRowsInPart(i);
00971
00972 if (rows !=1 && Containers_[i]->NumRows() == 0)
00973 continue;
00974
00975 int LID;
00976
00977
00978
00979 for (int j = 0 ; j < Partitioner_->NumRowsInPart(i) ; ++j) {
00980 LID = (*Partitioner_)(i,j);
00981 int NumEntries;
00982 IFPACK_CHK_ERR(Matrix().ExtractMyRowCopy(LID, Length,NumEntries,
00983 &Values[0], &Indices[0]));
00984
00985 for (int k = 0 ; k < NumEntries ; ++k) {
00986 int col = Indices[k];
00987
00988 for (int kk = 0 ; kk < NumVectors ; ++kk) {
00989 Xcopy[kk][LID] -= Values[k] * y2_ptr[kk][col];
00990 }
00991 }
00992 }
00993
00994
00995
00996 if(rows != 1) {
00997 for (int j = 0 ; j < Partitioner_->NumRowsInPart(i) ; ++j) {
00998 LID = Containers_[i]->ID(j);
00999 for (int k = 0 ; k < NumVectors ; ++k) {
01000 Containers_[i]->RHS(j,k) = Xcopy[k][LID];
01001 }
01002 }
01003
01004 IFPACK_CHK_ERR(Containers_[i]->ApplyInverse());
01005 #ifdef IFPACK_FLOPCOUNTERS
01006 ApplyInverseFlops_ += Containers_[i]->ApplyInverseFlops();
01007 #endif
01008
01009 for (int j = 0 ; j < Partitioner_->NumRowsInPart(i) ; ++j) {
01010 LID = Containers_[i]->ID(j);
01011 for (int k = 0 ; k < NumVectors ; ++k) {
01012 y2_ptr[k][LID] += DampingFactor_ * Containers_[i]->LHS(j,k);
01013 }
01014 }
01015 }
01016 else {
01017 int LRID = (*Partitioner_)(i,0);
01018 double b = Xcopy[0][LRID];
01019 double a = Diagonal_[LRID];
01020 y2_ptr[0][LRID]+= DampingFactor_* b/a;
01021 }
01022 }
01023
01024 #ifdef IFPACK_FLOPCOUNTERS
01025
01026 ApplyInverseFlops_ += NumVectors * 2 * Matrix_->NumGlobalNonzeros();
01027 ApplyInverseFlops_ += NumVectors * 2 * Matrix_->NumGlobalRows();
01028 #endif
01029
01030 Xcopy = X;
01031
01032 for (int i = NumLocalBlocks() - 1; i >=0 ; --i) {
01033 int rows = Partitioner_->NumRowsInPart(i);
01034 if (rows != 1 &&Containers_[i]->NumRows() == 0)
01035 continue;
01036
01037 int LID;
01038
01039
01040
01041 for (int j = 0 ; j < Partitioner_->NumRowsInPart(i) ; ++j) {
01042 LID = (*Partitioner_)(i,j);
01043
01044 int NumEntries;
01045 IFPACK_CHK_ERR(Matrix().ExtractMyRowCopy(LID, Length,NumEntries,
01046 &Values[0], &Indices[0]));
01047
01048 for (int k = 0 ; k < NumEntries ; ++k) {
01049 int col = Indices[k];
01050
01051 for (int kk = 0 ; kk < NumVectors ; ++kk) {
01052 Xcopy[kk][LID] -= Values[k] * y2_ptr[kk][col];
01053 }
01054 }
01055 }
01056
01057
01058 if(rows != 1) {
01059 for (int j = 0 ; j < Partitioner_->NumRowsInPart(i) ; ++j) {
01060 LID = Containers_[i]->ID(j);
01061 for (int k = 0 ; k < NumVectors ; ++k) {
01062 Containers_[i]->RHS(j,k) = Xcopy[k][LID];
01063 }
01064 }
01065
01066 IFPACK_CHK_ERR(Containers_[i]->ApplyInverse());
01067 #ifdef IFPACK_FLOPCOUNTERS
01068 ApplyInverseFlops_ += Containers_[i]->ApplyInverseFlops();
01069 #endif
01070
01071 for (int j = 0 ; j < Partitioner_->NumRowsInPart(i) ; ++j) {
01072 LID = Containers_[i]->ID(j);
01073 for (int k = 0 ; k < NumVectors ; ++k) {
01074 y2_ptr[k][LID] += DampingFactor_ * Containers_[i]->LHS(j,k);
01075 }
01076 }
01077 }
01078 else {
01079 int LRID = (*Partitioner_)(i,0);
01080 double b = Xcopy[0][LRID];
01081 double a = Diagonal_[LRID];
01082 y2_ptr[0][LRID]+= DampingFactor_* b/a;
01083 }
01084 }
01085
01086 #ifdef IFPACK_FLOPCOUNTERS
01087
01088 ApplyInverseFlops_ += NumVectors * 2 * Matrix_->NumGlobalNonzeros();
01089 ApplyInverseFlops_ += NumVectors * 2 * Matrix_->NumGlobalRows();
01090 #endif
01091
01092
01093
01094 if (IsParallel_)
01095 for (int m = 0 ; m < NumVectors ; ++m)
01096 for (int i = 0 ; i < NumMyRows ; ++i)
01097 y_ptr[m][i] = y2_ptr[m][i];
01098
01099 return(0);
01100 }
01101
01102
01103 template<typename T>
01104 ostream& Ifpack_BlockRelaxation<T>::Print(ostream & os) const
01105 {
01106
01107 string PT;
01108 if (PrecType_ == IFPACK_JACOBI)
01109 PT = "Jacobi";
01110 else if (PrecType_ == IFPACK_GS)
01111 PT = "Gauss-Seidel";
01112 else if (PrecType_ == IFPACK_SGS)
01113 PT = "symmetric Gauss-Seidel";
01114
01115 if (!Comm().MyPID()) {
01116 os << endl;
01117 os << "================================================================================" << endl;
01118 os << "Ifpack_BlockRelaxation, " << PT << endl;
01119 os << "Sweeps = " << NumSweeps_ << endl;
01120 os << "Damping factor = " << DampingFactor_;
01121 if (ZeroStartingSolution_)
01122 os << ", using zero starting solution" << endl;
01123 else
01124 os << ", using input starting solution" << endl;
01125 os << "Number of local blocks = " << Partitioner_->NumLocalParts() << endl;
01126
01127 os << "Global number of rows = " << Matrix_->NumGlobalRows64() << endl;
01128 os << endl;
01129 os << "Phase # calls Total Time (s) Total MFlops MFlops/s" << endl;
01130 os << "----- ------- -------------- ------------ --------" << endl;
01131 os << "Initialize() " << std::setw(5) << NumInitialize()
01132 << " " << std::setw(15) << InitializeTime()
01133 << " " << std::setw(15) << 1.0e-6 * InitializeFlops();
01134 if (InitializeTime() != 0.0)
01135 os << " " << std::setw(15) << 1.0e-6 * InitializeFlops() / InitializeTime() << endl;
01136 else
01137 os << " " << std::setw(15) << 0.0 << endl;
01138 os << "Compute() " << std::setw(5) << NumCompute()
01139 << " " << std::setw(15) << ComputeTime()
01140 << " " << std::setw(15) << 1.0e-6 * ComputeFlops();
01141 if (ComputeTime() != 0.0)
01142 os << " " << std::setw(15) << 1.0e-6 * ComputeFlops() / ComputeTime() << endl;
01143 else
01144 os << " " << std::setw(15) << 0.0 << endl;
01145 os << "ApplyInverse() " << std::setw(5) << NumApplyInverse()
01146 << " " << std::setw(15) << ApplyInverseTime()
01147 << " " << std::setw(15) << 1.0e-6 * ApplyInverseFlops();
01148 if (ApplyInverseTime() != 0.0)
01149 os << " " << std::setw(15) << 1.0e-6 * ApplyInverseFlops() / ApplyInverseTime() << endl;
01150 else
01151 os << " " << std::setw(15) << 0.0 << endl;
01152 os << "================================================================================" << endl;
01153 os << endl;
01154 }
01155
01156 return(os);
01157 }
01158
01159
01160 template<typename T>
01161 int Ifpack_BlockRelaxation<T>::SetParameters(Teuchos::ParameterList& List)
01162 {
01163
01164 string PT;
01165 if (PrecType_ == IFPACK_JACOBI)
01166 PT = "Jacobi";
01167 else if (PrecType_ == IFPACK_GS)
01168 PT = "Gauss-Seidel";
01169 else if (PrecType_ == IFPACK_SGS)
01170 PT = "symmetric Gauss-Seidel";
01171
01172 PT = List.get("relaxation: type", PT);
01173
01174 if (PT == "Jacobi") {
01175 PrecType_ = IFPACK_JACOBI;
01176 }
01177 else if (PT == "Gauss-Seidel") {
01178 PrecType_ = IFPACK_GS;
01179 }
01180 else if (PT == "symmetric Gauss-Seidel") {
01181 PrecType_ = IFPACK_SGS;
01182 } else {
01183 cerr << "Option `relaxation: type' has an incorrect value ("
01184 << PT << ")'" << endl;
01185 cerr << "(file " << __FILE__ << ", line " << __LINE__ << ")" << endl;
01186 exit(EXIT_FAILURE);
01187 }
01188
01189 NumSweeps_ = List.get("relaxation: sweeps", NumSweeps_);
01190 DampingFactor_ = List.get("relaxation: damping factor",
01191 DampingFactor_);
01192 ZeroStartingSolution_ = List.get("relaxation: zero starting solution",
01193 ZeroStartingSolution_);
01194 PartitionerType_ = List.get("partitioner: type",
01195 PartitionerType_);
01196 NumLocalBlocks_ = List.get("partitioner: local parts",
01197 NumLocalBlocks_);
01198
01199 OverlapLevel_ = List.get("partitioner: overlap",
01200 OverlapLevel_);
01201
01202
01203 if (PrecType_ != IFPACK_JACOBI)
01204 OverlapLevel_ = 0;
01205 if (NumLocalBlocks_ < 0)
01206 NumLocalBlocks_ = Matrix().NumMyRows() / (-NumLocalBlocks_);
01207
01208
01209
01210
01211 List_ = List;
01212
01213
01214 string PT2;
01215 if (PrecType_ == IFPACK_JACOBI)
01216 PT2 = "BJ";
01217 else if (PrecType_ == IFPACK_GS)
01218 PT2 = "BGS";
01219 else if (PrecType_ == IFPACK_SGS)
01220 PT2 = "BSGS";
01221 Label_ = "IFPACK (" + PT2 + ", sweeps="
01222 + Ifpack_toString(NumSweeps_) + ", damping="
01223 + Ifpack_toString(DampingFactor_) + ", blocks="
01224 + Ifpack_toString(NumLocalBlocks()) + ")";
01225
01226 return(0);
01227 }
01228
01229
01230 template<typename T>
01231 int Ifpack_BlockRelaxation<T>::Initialize()
01232 {
01233 IsInitialized_ = false;
01234 Time_.ResetStartTime();
01235
01236 Graph_ = Teuchos::rcp( new Ifpack_Graph_Epetra_RowMatrix(Teuchos::rcp(&Matrix(),false)) );
01237 if (Graph_ == Teuchos::null) IFPACK_CHK_ERR(-5);
01238
01239 if (PartitionerType_ == "linear")
01240 Partitioner_ = Teuchos::rcp( new Ifpack_LinearPartitioner(&*Graph_) );
01241 else if (PartitionerType_ == "greedy")
01242 Partitioner_ = Teuchos::rcp( new Ifpack_GreedyPartitioner(&*Graph_) );
01243 else if (PartitionerType_ == "metis")
01244 Partitioner_ = Teuchos::rcp( new Ifpack_METISPartitioner(&*Graph_) );
01245 else if (PartitionerType_ == "equation")
01246 Partitioner_ = Teuchos::rcp( new Ifpack_EquationPartitioner(&*Graph_) );
01247 else if (PartitionerType_ == "user")
01248 Partitioner_ = Teuchos::rcp( new Ifpack_UserPartitioner(&*Graph_) );
01249 else if (PartitionerType_ == "line")
01250 Partitioner_ = Teuchos::rcp( new Ifpack_LinePartitioner(&*Graph_) );
01251 else
01252 IFPACK_CHK_ERR(-2);
01253
01254 if (Partitioner_ == Teuchos::null) IFPACK_CHK_ERR(-5);
01255
01256
01257 IFPACK_CHK_ERR(Partitioner_->SetParameters(List_));
01258 IFPACK_CHK_ERR(Partitioner_->Compute());
01259
01260
01261 NumLocalBlocks_ = Partitioner_->NumLocalParts();
01262
01263
01264 W_ = Teuchos::rcp( new Epetra_Vector(Matrix().RowMatrixRowMap()) );
01265 W_->PutScalar(0.0);
01266
01267 for (int i = 0 ; i < NumLocalBlocks() ; ++i) {
01268
01269 for (int j = 0 ; j < Partitioner_->NumRowsInPart(i) ; ++j) {
01270 int LID = (*Partitioner_)(i,j);
01271 (*W_)[LID]++;
01272 }
01273 }
01274 W_->Reciprocal(*W_);
01275
01276
01277 if (PartitionerType_ == "line") {
01278
01279 string PT2;
01280 if (PrecType_ == IFPACK_JACOBI)
01281 PT2 = "BJ";
01282 else if (PrecType_ == IFPACK_GS)
01283 PT2 = "BGS";
01284 else if (PrecType_ == IFPACK_SGS)
01285 PT2 = "BSGS";
01286 Label_ = "IFPACK (" + PT2 + ", auto-line, sweeps="
01287 + Ifpack_toString(NumSweeps_) + ", damping="
01288 + Ifpack_toString(DampingFactor_) + ", blocks="
01289 + Ifpack_toString(NumLocalBlocks()) + ")";
01290 }
01291
01292 InitializeTime_ += Time_.ElapsedTime();
01293 IsInitialized_ = true;
01294 ++NumInitialize_;
01295
01296 return(0);
01297 }
01298
01299
01300 #endif // IFPACK_BLOCKPRECONDITIONER_H