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 #include "Ifpack_IHSS.h"
00045 #include "Ifpack.h"
00046 #include "Ifpack_Utils.h"
00047 #include "Teuchos_ParameterList.hpp"
00048 #include "Teuchos_RefCountPtr.hpp"
00049
00050
00051
00052 using Teuchos::RefCountPtr;
00053 using Teuchos::rcp;
00054
00055
00056 #ifdef HAVE_IFPACK_EPETRAEXT
00057 #include "EpetraExt_MatrixMatrix.h"
00058
00059
00060 Ifpack_IHSS::Ifpack_IHSS(Epetra_RowMatrix* A):
00061 IsInitialized_(false),
00062 IsComputed_(false),
00063 Label_(),
00064 EigMaxIters_(10),
00065 EigRatio_(30.0),
00066 LambdaMax_(-1.0),
00067 Alpha_(-1.0),
00068 NumSweeps_(1),
00069
00070 Time_(A->Comm())
00071 {
00072 Epetra_CrsMatrix *Acrs=dynamic_cast<Epetra_CrsMatrix*>(A);
00073 TEUCHOS_TEST_FOR_EXCEPT(!Acrs)
00074 A_=rcp(Acrs,false);
00075 }
00076
00077 void Ifpack_IHSS::Destroy(){
00078 }
00079
00080
00081
00082 int Ifpack_IHSS::Initialize(){
00083 EigMaxIters_ = List_.get("ihss: eigenvalue max iterations",EigMaxIters_);
00084 EigRatio_ = List_.get("ihss: ratio eigenvalue", EigRatio_);
00085 NumSweeps_ = List_.get("ihss: sweeps",NumSweeps_);
00086
00087
00088 IsInitialized_=true;
00089 NumInitialize_++;
00090 return 0;
00091 }
00092
00093 int Ifpack_IHSS::SetParameters(Teuchos::ParameterList& parameterlist){
00094 List_=parameterlist;
00095 return 0;
00096 }
00097
00098
00099 int Ifpack_IHSS::Compute(){
00100 if(!IsInitialized_) Initialize();
00101
00102 int rv;
00103 Ifpack Factory;
00104 Epetra_CrsMatrix *Askew=0,*Aherm=0;
00105 Ifpack_Preconditioner *Pskew=0, *Pherm=0;
00106 Time_.ResetStartTime();
00107
00108
00109 rv=EpetraExt::MatrixMatrix::Add(*A_,false,.5,*A_,true,.5,Aherm);
00110 Aherm->FillComplete();
00111 if(rv) IFPACK_CHK_ERR(-1);
00112
00113
00114 Epetra_Vector avec(Aherm->RowMap());
00115 IFPACK_CHK_ERR(Aherm->ExtractDiagonalCopy(avec));
00116
00117
00118
00119
00120
00121
00122
00123 avec.MaxValue(&Alpha_);
00124
00125
00126 for(int i=0;i<Aherm->NumMyRows();i++) avec[i]+=Alpha_;
00127 IFPACK_CHK_ERR(Aherm->ReplaceDiagonalValues(avec));
00128 Aherm_=rcp(Aherm);
00129
00130
00131 Askew=new Epetra_CrsMatrix(Copy,A_->RowMap(),0);
00132 rv=EpetraExt::MatrixMatrix::Add(*A_,false,.5,*A_,true,-.5,Askew);
00133 if(rv) IFPACK_CHK_ERR(-2);
00134
00135 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
00136 if(Askew->RowMap().GlobalIndicesInt()) {
00137 for(int i=0;i<Askew->NumMyRows();i++) {
00138 int gid=Askew->GRID(i);
00139 Askew->InsertGlobalValues(gid,1,&Alpha_,&gid);
00140 }
00141 }
00142 else
00143 #endif
00144 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
00145 if(Askew->RowMap().GlobalIndicesLongLong()) {
00146 for(int i=0;i<Askew->NumMyRows();i++) {
00147 long long gid=Askew->GRID64(i);
00148 Askew->InsertGlobalValues(gid,1,&Alpha_,&gid);
00149 }
00150 }
00151 else
00152 #endif
00153 throw "Ifpack_IHSS::Compute: Unable to deduce GlobalIndices type";
00154
00155 Askew->FillComplete();
00156 Askew_=rcp(Askew);
00157
00158
00159 Teuchos::ParameterList PLh=List_.sublist("ihss: hermetian list");
00160 string htype=List_.get("ihss: hermetian type","ILU");
00161 Pherm= Factory.Create(htype, Aherm);
00162 Pherm->SetParameters(PLh);
00163 IFPACK_CHK_ERR(Pherm->Compute());
00164 Pherm_=rcp(Pherm);
00165
00166
00167 Teuchos::ParameterList PLs=List_.sublist("ihss: skew hermetian list");
00168 string stype=List_.get("ihss: skew hermetian type","ILU");
00169 Pskew= Factory.Create(stype, Askew);
00170 Pskew->SetParameters(PLs);
00171 IFPACK_CHK_ERR(Pskew->Compute());
00172 Pskew_=rcp(Pskew);
00173
00174
00175 sprintf(Label_, "IFPACK IHSS (H,S)=(%s/%s)",htype.c_str(),stype.c_str());
00176
00177
00178 IsComputed_=true;
00179 NumCompute_++;
00180 ComputeTime_ += Time_.ElapsedTime();
00181 return 0;
00182 }
00183
00184
00185 int Ifpack_IHSS::ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const{
00186 if(!IsComputed_) return -1;
00187 Time_.ResetStartTime();
00188 bool initial_guess_is_zero=false;
00189
00190
00191
00192 Teuchos::RefCountPtr<const Epetra_MultiVector> Xcopy;
00193 Epetra_MultiVector Temp(X);
00194 if (X.Pointers()[0] == Y.Pointers()[0]){
00195 Xcopy = Teuchos::rcp( new Epetra_MultiVector(X) );
00196
00197 Y.Scale(0.0);
00198 initial_guess_is_zero=true;
00199 }
00200 else
00201 Xcopy = Teuchos::rcp( &X, false );
00202
00203 Epetra_MultiVector T1(Y),T2(*Xcopy);
00204
00205
00206
00207
00208
00209
00210
00211
00212
00213
00214
00215 for(int i=0;i<NumSweeps_;i++){
00216
00217 if(!initial_guess_is_zero || i >0 ){
00218 Askew_->Apply(Y,T1);
00219 T2.Update(2*Alpha_,Y,-1,T1,1);
00220 }
00221 Pherm_->ApplyInverse(T2,Y);
00222
00223
00224 Aherm_->Apply(Y,T1);
00225 T2.Scale(1.0,*Xcopy);
00226 T2.Update(2*Alpha_,Y,-1,T1,1.0);
00227 Pskew_->ApplyInverse(T2,Y);
00228 }
00229
00230
00231 NumApplyInverse_++;
00232 ApplyInverseTime_ += Time_.ElapsedTime();
00233 return 0;
00234 }
00235
00236
00237 ostream& Ifpack_IHSS::Print(ostream& os) const{
00238 os<<"Ifpack_IHSS"<<endl;
00239 os<<"-Hermetian preconditioner"<<endl;
00240 os<<"-Skew Hermetian preconditioner"<<endl;
00241 os<<endl;
00242 return os;
00243 }
00244
00245
00246 double Ifpack_IHSS::Condest(const Ifpack_CondestType CT,
00247 const int MaxIters,
00248 const double Tol,
00249 Epetra_RowMatrix* Matrix_in){
00250 return -1.0;
00251 }
00252
00253
00254 int Ifpack_IHSS::PowerMethod(Epetra_Operator * Op,const int MaximumIterations, double& lambda_max)
00255 {
00256
00257
00258 lambda_max = 0.0;
00259 double RQ_top, RQ_bottom, norm;
00260 Epetra_Vector x(Op->OperatorDomainMap());
00261 Epetra_Vector y(Op->OperatorRangeMap());
00262 x.Random();
00263 x.Norm2(&norm);
00264 if (norm == 0.0) IFPACK_CHK_ERR(-1);
00265
00266 x.Scale(1.0 / norm);
00267
00268 for (int iter = 0; iter < MaximumIterations; ++iter)
00269 {
00270 Op->Apply(x, y);
00271 IFPACK_CHK_ERR(y.Dot(x, &RQ_top));
00272 IFPACK_CHK_ERR(x.Dot(x, &RQ_bottom));
00273 lambda_max = RQ_top / RQ_bottom;
00274 IFPACK_CHK_ERR(y.Norm2(&norm));
00275 if (norm == 0.0) IFPACK_CHK_ERR(-1);
00276 IFPACK_CHK_ERR(x.Update(1.0 / norm, y, 0.0));
00277 }
00278
00279 return(0);
00280 }
00281
00282 #endif