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_SORa.h"
00045 #include "Ifpack.h"
00046 #include "Ifpack_Utils.h"
00047 #include "Teuchos_ParameterList.hpp"
00048 #include "Teuchos_RefCountPtr.hpp"
00049 #include "Epetra_Import.h"
00050 #include "Epetra_Export.h"
00051 #include "Epetra_IntVector.h"
00052
00053 using Teuchos::RefCountPtr;
00054 using Teuchos::rcp;
00055
00056
00057
00058 #ifdef HAVE_IFPACK_EPETRAEXT
00059 #include "EpetraExt_MatrixMatrix.h"
00060 #include "EpetraExt_RowMatrixOut.h"
00061 #include "EpetraExt_MultiVectorOut.h"
00062 #include "EpetraExt_Transpose_RowMatrix.h"
00063
00064
00065 #define ABS(x) ((x)>=0?(x):-(x))
00066 #define MIN(x,y) ((x)<(y)?(x):(y))
00067 #define MAX(x,y) ((x)>(y)?(x):(y))
00068
00069
00070 int* FindLocalDiricheltRowsFromOnesAndZeros(const Epetra_CrsMatrix & Matrix, int &numBCRows);
00071 void Apply_BCsToMatrixRowsAndColumns(const int *dirichletRows, int numBCRows,const Epetra_IntVector &dirichletColumns,const Epetra_CrsMatrix & Matrix);
00072 Epetra_IntVector * FindLocalDirichletColumnsFromRows(const int *dirichletRows, int numBCRows,const Epetra_CrsMatrix & Matrix);
00073
00074 Ifpack_SORa::Ifpack_SORa(Epetra_RowMatrix* A):
00075 IsInitialized_(false),
00076 IsComputed_(false),
00077 Label_(),
00078 Alpha_(1.5),
00079 Gamma_(1.0),
00080 NumSweeps_(1),
00081 IsParallel_(false),
00082 HaveOAZBoundaries_(false),
00083 UseInterprocDamping_(false),
00084 UseGlobalDamping_(false),
00085 LambdaMax_(-1.0),
00086 Time_(A->Comm())
00087 {
00088 Epetra_CrsMatrix *Acrs=dynamic_cast<Epetra_CrsMatrix*>(A);
00089 if(Acrs) Acrs_=rcp(Acrs,false);
00090 A_=rcp(A,false);
00091 }
00092
00093 void Ifpack_SORa::Destroy(){
00094 }
00095
00096
00097
00098 int Ifpack_SORa::Initialize(){
00099 Alpha_ = List_.get("sora: alpha", Alpha_);
00100 Gamma_ = List_.get("sora: gamma",Gamma_);
00101 NumSweeps_ = List_.get("sora: sweeps",NumSweeps_);
00102 HaveOAZBoundaries_= List_.get("sora: oaz boundaries", HaveOAZBoundaries_);
00103 UseInterprocDamping_ = List_.get("sora: use interproc damping",UseInterprocDamping_);
00104 UseGlobalDamping_ = List_.get("sora: use global damping",UseGlobalDamping_);
00105
00106 if (A_->Comm().NumProc() != 1) IsParallel_ = true;
00107 else {
00108 IsParallel_ = false;
00109
00110 UseInterprocDamping_=false;
00111 }
00112
00113
00114 IsInitialized_=true;
00115 NumInitialize_++;
00116 return 0;
00117 }
00118
00119 int Ifpack_SORa::SetParameters(Teuchos::ParameterList& parameterlist){
00120 List_=parameterlist;
00121 return 0;
00122 }
00123
00124
00125 template<typename int_type>
00126 int Ifpack_SORa::TCompute(){
00127 Epetra_Map *RowMap=const_cast<Epetra_Map*>(&A_->RowMatrixRowMap());
00128 Epetra_Vector Adiag(*RowMap);
00129 Epetra_CrsMatrix *Askew2=0, *Aherm2=0,*W=0;
00130 int *rowptr_s,*colind_s,*rowptr_h,*colind_h;
00131 double *vals_s,*vals_h;
00132 bool RowMatrixMode=(Acrs_==Teuchos::null);
00133
00134
00135 sprintf(Label_, "IFPACK SORa (alpha=%5.2f gamma=%5.2f)",Alpha_,Gamma_);
00136
00137 if(RowMatrixMode){
00138 if(!A_->Comm().MyPID()) cout<<"SORa: RowMatrix mode enabled"<<endl;
00139
00140 Epetra_RowMatrix *Arow=&*A_;
00141 Epetra_Map *ColMap=const_cast<Epetra_Map*>(&A_->RowMatrixColMap());
00142
00143 int Nmax=A_->MaxNumEntries();
00144 int length;
00145 std::vector<double> values(Nmax);
00146 Epetra_CrsMatrix *Acrs=new Epetra_CrsMatrix(Copy,*RowMap,Nmax);
00147
00148 std::vector<int> indices_local(Nmax);
00149 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
00150 if(RowMap->GlobalIndicesInt()) {
00151 for(int i=0;i<Arow->NumMyRows();i++) {
00152 Arow->ExtractMyRowCopy(i,Nmax,length,&values[0],&indices_local[0]);
00153 for(int j=0;j<length;j++)
00154 indices_local[j]= ColMap->GID(indices_local[j]);
00155 Acrs->InsertGlobalValues(RowMap->GID(i),length,&values[0],&indices_local[0]);
00156 }
00157 }
00158 else
00159 #endif
00160 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
00161 if(RowMap->GlobalIndicesLongLong()) {
00162 std::vector<int_type> indices_global(Nmax);
00163 for(int i=0;i<Arow->NumMyRows();i++) {
00164 Arow->ExtractMyRowCopy(i,Nmax,length,&values[0],&indices_local[0]);
00165 for(int j=0;j<length;j++)
00166 indices_global[j]= (int_type) ColMap->GID64(indices_local[j]);
00167 Acrs->InsertGlobalValues((int_type) RowMap->GID64(i),length,&values[0],&indices_global[0]);
00168 }
00169 }
00170 else
00171 #endif
00172 throw "Ifpack_SORa::TCompute: GlobalIndices type unknown";
00173
00174 Acrs->FillComplete(A_->OperatorDomainMap(),A_->OperatorRangeMap());
00175 Acrs_=rcp(Acrs,true);
00176 }
00177
00178
00179
00180
00181 IFPACK_CHK_ERR(EpetraExt::MatrixMatrix::Add(*Acrs_,false,1,*Acrs_,true,-1,Askew2));
00182 Askew2->FillComplete();
00183
00184
00185 IFPACK_CHK_ERR(EpetraExt::MatrixMatrix::Add(*Acrs_,false,1,*Acrs_,true,1,Aherm2));
00186 Aherm2->FillComplete();
00187
00188 int nnz2=Askew2->NumMyNonzeros();
00189 int N=Askew2->NumMyRows();
00190
00191
00192
00193 IFPACK_CHK_ERR(Askew2->ExtractCrsDataPointers(rowptr_s,colind_s,vals_s));
00194 IFPACK_CHK_ERR(Aherm2->ExtractCrsDataPointers(rowptr_h,colind_h,vals_h));
00195
00196
00197 #define SANITY_CHECK
00198 #ifdef SANITY_CHECK
00199 for(int i=0;i<N;i++)
00200 if(rowptr_s[i]!=rowptr_h[i]) IFPACK_CHK_ERR(-2);
00201 for(int i=0;i<nnz2;i++)
00202 if(colind_s[i]!=colind_h[i]) IFPACK_CHK_ERR(-3);
00203 #endif
00204
00205
00206
00207 if(HaveOAZBoundaries_){
00208 int numBCRows;
00209 int* dirRows=FindLocalDiricheltRowsFromOnesAndZeros(*Acrs_,numBCRows);
00210 Epetra_IntVector* dirCols=FindLocalDirichletColumnsFromRows(dirRows,numBCRows,*Aherm2);
00211 Apply_BCsToMatrixRowsAndColumns(dirRows,numBCRows,*dirCols,*Aherm2);
00212 Apply_BCsToMatrixRowsAndColumns(dirRows,numBCRows,*dirCols,*Askew2);
00213 delete [] dirRows;
00214 delete dirCols;
00215 }
00216
00217
00218 A_->ExtractDiagonalCopy(Adiag);
00219
00220
00221 Epetra_Vector *Wdiag = new Epetra_Vector(*RowMap);
00222
00223
00224
00225 int maxentries=Askew2->MaxNumEntries();
00226 int_type* gids=new int_type [maxentries];
00227 double* newvals=new double[maxentries];
00228 W=new Epetra_CrsMatrix(Copy,*RowMap,0);
00229 for(int i=0;i<N;i++){
00230
00231 int_type rowgid = (int_type) Acrs_->GRID64(i);
00232 double c_data=0.0;
00233 double ipdamp=0.0;
00234 int idx=0;
00235
00236 for(int j=rowptr_s[i];j<rowptr_s[i+1];j++){
00237 int_type colgid = (int_type) Askew2->GCID64(colind_s[j]);
00238 c_data+=fabs(vals_s[j]);
00239 if(rowgid>colgid){
00240
00241 if(colind_s[j] < N) {
00242 gids[idx]=colgid;
00243 newvals[idx]=vals_h[j]/2 + Alpha_ * vals_s[j]/2;
00244 idx++;
00245 }
00246 else{
00247 ipdamp+=fabs(vals_h[j]/2 + Alpha_ * vals_s[j]/2);
00248 }
00249 }
00250 }
00251 if(idx>0)
00252 IFPACK_CHK_ERR(W->InsertGlobalValues(rowgid,idx,newvals,gids));
00253
00254
00255 double w_val= c_data*Alpha_*Gamma_/4 + Adiag[Acrs_->LRID(rowgid)];
00256 if(UseInterprocDamping_) w_val+=ipdamp;
00257
00258 W->InsertGlobalValues(rowgid,1,&w_val,&rowgid);
00259 IFPACK_CHK_ERR(Wdiag->ReplaceGlobalValues(1,&w_val,&rowgid));
00260 }
00261 W->FillComplete(A_->OperatorDomainMap(),A_->OperatorRangeMap());
00262 W_=rcp(W);
00263 Wdiag_=rcp(Wdiag);
00264
00265
00266 IsComputed_=true;
00267
00268
00269 if(UseGlobalDamping_) {
00270 PowerMethod(10,LambdaMax_);
00271 if(!A_->Comm().MyPID()) printf("SORa: Global damping parameter = %6.4e (lmax=%6.4e)\n",GetOmega(),LambdaMax_);
00272 }
00273
00274
00275 delete [] gids;
00276 delete [] newvals;
00277 delete Aherm2;
00278 delete Askew2;
00279 if(RowMatrixMode) {
00280 Acrs_=Teuchos::null;
00281 }
00282
00283
00284 NumCompute_++;
00285 ComputeTime_ += Time_.ElapsedTime();
00286 return 0;
00287 }
00288
00289 int Ifpack_SORa::Compute(){
00290 if(!IsInitialized_) Initialize();
00291 int ret_val = 0;
00292 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
00293 if(A_->RowMatrixRowMap().GlobalIndicesInt()) {
00294 ret_val = TCompute<int>();
00295 }
00296 else
00297 #endif
00298 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
00299 if(A_->RowMatrixRowMap().GlobalIndicesLongLong()) {
00300 ret_val = TCompute<long long>();
00301 }
00302 else
00303 #endif
00304 throw "Ifpack_SORa::Compute: GlobalIndices type unknown";
00305
00306 return ret_val;
00307 }
00308
00309
00310 int Ifpack_SORa::ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const{
00311 if(!IsComputed_) return -1;
00312 Time_.ResetStartTime();
00313 bool initial_guess_is_zero=false;
00314 int NumMyRows=W_->NumMyRows();
00315 int NumVectors = X.NumVectors();
00316 Epetra_MultiVector Temp(A_->RowMatrixRowMap(),NumVectors);
00317
00318 double omega=GetOmega();
00319
00320
00321 Teuchos::RefCountPtr<const Epetra_MultiVector> Xcopy;
00322 if (X.Pointers()[0] == Y.Pointers()[0]){
00323 Xcopy = Teuchos::rcp( new Epetra_MultiVector(X) );
00324
00325 Y.Scale(0.0);
00326 initial_guess_is_zero=true;
00327 }
00328 else
00329 Xcopy = Teuchos::rcp( &X, false );
00330
00331 Teuchos::RefCountPtr< Epetra_MultiVector > T2;
00332
00333
00334
00335
00336
00337 if (IsParallel_ && W_->Importer()) T2 = Teuchos::rcp( new Epetra_MultiVector(W_->Importer()->TargetMap(),NumVectors,true));
00338 else T2 = Teuchos::rcp( new Epetra_MultiVector(A_->RowMatrixRowMap(),NumVectors,true));
00339
00340
00341 int* rowptr,*colind;
00342 double *values;
00343 double **t_ptr,** y_ptr, ** t2_ptr, **x_ptr,*d_ptr;
00344 T2->ExtractView(&t2_ptr);
00345 Y.ExtractView(&y_ptr);
00346 Temp.ExtractView(&t_ptr);
00347 Xcopy->ExtractView(&x_ptr);
00348 Wdiag_->ExtractView(&d_ptr);
00349 IFPACK_CHK_ERR(W_->ExtractCrsDataPointers(rowptr,colind,values));
00350
00351
00352 for(int i=0; i<NumSweeps_; i++){
00353
00354 if(!initial_guess_is_zero || i > 0) {
00355 A_->Apply(Y,Temp);
00356 Temp.Update(1.0,*Xcopy,-1.0);
00357 }
00358 else
00359 Temp.Update(1.0,*Xcopy,0.0);
00360
00361
00362
00363
00364
00365
00366 for(int j=0; j<NumMyRows; j++){
00367 double diag=d_ptr[j];
00368 for (int m=0 ; m<NumVectors; m++) {
00369 double dtmp=0.0;
00370
00371 t2_ptr[m][j]=0.0;
00372 for(int k=rowptr[j];k<rowptr[j+1];k++){
00373 dtmp+= values[k]*t2_ptr[m][colind[k]];
00374 }
00375
00376 t2_ptr[m][j] = (t_ptr[m][j]- dtmp)/diag;
00377 y_ptr[m][j] += omega*t2_ptr[m][j];
00378 }
00379 }
00380 }
00381
00382
00383 NumApplyInverse_++;
00384 ApplyInverseTime_ += Time_.ElapsedTime();
00385 return 0;
00386 }
00387
00388
00389 ostream& Ifpack_SORa::Print(ostream& os) const{
00390 os<<"Ifpack_SORa"<<endl;
00391 os<<" alpha = "<<Alpha_<<endl;
00392 os<<" gamma = "<<Gamma_<<endl;
00393 os<<endl;
00394 return os;
00395 }
00396
00397
00398 double Ifpack_SORa::Condest(const Ifpack_CondestType CT,
00399 const int MaxIters,
00400 const double Tol,
00401 Epetra_RowMatrix* Matrix_in){
00402 return -1.0;
00403 }
00404
00405
00406
00407
00408
00409
00410 inline int* FindLocalDiricheltRowsFromOnesAndZeros(const Epetra_CrsMatrix & Matrix, int &numBCRows){
00411 int *dirichletRows = new int[Matrix.NumMyRows()];
00412 numBCRows = 0;
00413 for (int i=0; i<Matrix.NumMyRows(); i++) {
00414 int numEntries, *cols;
00415 double *vals;
00416 int ierr = Matrix.ExtractMyRowView(i,numEntries,vals,cols);
00417 if (ierr == 0) {
00418 int nz=0;
00419 for (int j=0; j<numEntries; j++) if (vals[j]!=0) nz++;
00420 if (nz==1) dirichletRows[numBCRows++] = i;
00421
00422 if(nz==2) dirichletRows[numBCRows++] = i;
00423 }
00424 }
00425 return dirichletRows;
00426 }
00427
00428
00429
00431 inline Epetra_IntVector * FindLocalDirichletColumnsFromRows(const int *dirichletRows, int numBCRows,const Epetra_CrsMatrix & Matrix){
00432
00433
00434
00435 Epetra_IntVector ZeroRows(Matrix.RowMap());
00436 Epetra_IntVector *ZeroCols=new Epetra_IntVector(Matrix.ColMap());
00437 ZeroRows.PutValue(0);
00438 ZeroCols->PutValue(0);
00439
00440
00441 if(Matrix.RowMap().SameAs(Matrix.ColMap())){
00442 for (int i=0; i < numBCRows; i++)
00443 (*ZeroCols)[dirichletRows[i]]=1;
00444 return ZeroCols;
00445 }
00446
00447
00448 for (int i=0; i < numBCRows; i++)
00449 ZeroRows[dirichletRows[i]]=1;
00450
00451
00452 if(Matrix.RowMap().SameAs(Matrix.DomainMap())){
00453
00454 ZeroCols->Import(ZeroRows,*Matrix.Importer(),Insert);
00455 }
00456 else{
00457
00458 Epetra_Import Importer(Matrix.ColMap(),Matrix.RowMap());
00459 ZeroCols->Import(ZeroRows,Importer,Insert);
00460 }
00461 return ZeroCols;
00462 }
00463
00464
00465
00466 inline void Apply_BCsToMatrixRowsAndColumns(const int *dirichletRows, int numBCRows,const Epetra_IntVector &dirichletColumns,const Epetra_CrsMatrix & Matrix){
00467
00468
00469
00470
00471 for(int i=0;i<numBCRows;i++){
00472 int numEntries, *cols;
00473 double *vals;
00474 Matrix.ExtractMyRowView(dirichletRows[i],numEntries,vals,cols);
00475 for (int j=0; j<numEntries; j++) vals[j]=0.0;
00476 }
00477
00478
00479 for (int i=0; i < Matrix.NumMyRows(); i++) {
00480 int numEntries;
00481 double *vals;
00482 int *cols;
00483 Matrix.ExtractMyRowView(i,numEntries,vals,cols);
00484 for (int j=0; j < numEntries; j++) {
00485 if (dirichletColumns[ cols[j] ] > 0) vals[j] = 0.0;
00486 }
00487 }
00488 }
00489
00490
00491
00492
00493
00494 int Ifpack_SORa::
00495 PowerMethod(const int MaximumIterations, double& lambda_max)
00496 {
00497
00498 lambda_max = 0.0;
00499 double RQ_top, RQ_bottom, norm;
00500 Epetra_Vector x(A_->OperatorDomainMap());
00501 Epetra_Vector y(A_->OperatorRangeMap());
00502 Epetra_Vector z(A_->OperatorRangeMap());
00503 x.Random();
00504 x.Norm2(&norm);
00505 if (norm == 0.0) IFPACK_CHK_ERR(-1);
00506
00507 x.Scale(1.0 / norm);
00508
00509
00510 int NumSweepsBackup=NumSweeps_;
00511 bool UseGlobalDampingBackup=UseGlobalDamping_;
00512 NumSweeps_=1;UseGlobalDamping_=false;
00513
00514 for (int iter = 0; iter < MaximumIterations; ++iter)
00515 {
00516 y.PutScalar(0.0);
00517 A_->Apply(x, z);
00518 ApplyInverse(z,y);
00519 y.Update(1.0,x,-1.0);
00520
00521
00522 y.Dot(x, &RQ_top);
00523 x.Dot(x, &RQ_bottom);
00524 lambda_max = RQ_top / RQ_bottom;
00525 y.Norm2(&norm);
00526 if (norm == 0.0) IFPACK_CHK_ERR(-1);
00527 x.Update(1.0 / norm, y, 0.0);
00528
00529 }
00530
00531
00532
00533
00534
00535 NumSweeps_=NumSweepsBackup;
00536 UseGlobalDamping_=UseGlobalDampingBackup;
00537
00538 return(0);
00539 }
00540
00541 #endif