00001 #ifdef CH_LANG_CC
00002
00003
00004
00005
00006
00007
00008
00009 #endif
00010
00011
00012
00013 #ifndef _GMRESSOLVER_H_
00014 #define _GMRESSOLVER_H_
00015
00016 #include "LinearSolver.H"
00017 #include "parstream.H"
00018 #include "CH_Timer.H"
00019 #include "NamespaceHeader.H"
00020
00021
00022
00023
00024
00025
00026
00027
00028 template <class T>
00029 class GMRESSolver : public LinearSolver<T>
00030 {
00031 public:
00032
00033 GMRESSolver();
00034
00035 virtual ~GMRESSolver();
00036
00037 void clearData();
00038
00039 virtual void setHomogeneous(bool a_homogeneous)
00040 {
00041 m_homogeneous = a_homogeneous;
00042 }
00043
00044
00045 virtual void setConvergenceMetrics(Real a_metric,
00046 Real a_tolerance);
00047
00048
00049
00050
00051
00052
00053
00054 virtual void define(LinearOp<T>* a_op, bool a_homogeneous=true);
00055
00056
00057 virtual void solve(T& a_phi, const T& a_rhs);
00058
00059
00060
00061
00062
00063
00064 bool m_homogeneous;
00065
00066
00067
00068
00069
00070 LinearOp<T>* m_op;
00071
00072
00073
00074
00075
00076 private:
00077 int m_restrtLen;
00078 public:
00079 int restartLen()const
00080 {
00081 return m_restrtLen;
00082 }
00083 void setRestartLen(int mm);
00084
00085
00086
00087
00088 int m_imax;
00089
00090
00091
00092
00093
00094 int m_verbosity;
00095
00096
00097
00098
00099
00100 Real m_eps;
00101
00102
00103
00104
00105
00106 Real m_reps;
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116 int m_exitStatus;
00117
00118
00119
00120
00121
00122 Real m_small;
00123
00124 int m_normType;
00125
00126 private:
00127 void allocate();
00128 void CycleGMRES( T &xx,
00129 int &reason, int &itcount, Real &rnorm0,
00130 const bool avoidnorms = false);
00131
00132 void ResidualGMRES( T &a_vv, const T &a_xx,
00133 const T &a_bb, T &a_temp );
00134
00135 void BuildGMRESSoln( Real nrs[], T &a_xx, const int it,
00136 const T vv_0[] );
00137
00138 void UpdateGMRESHessenberg( const int it, bool hapend, Real &res );
00139
00140
00141
00142 void TwoUnmodifiedGramSchmidtOrthogonalization( const int it );
00143
00144 void ApplyAB( T &a_dest, const T &a_xx, T &a_temp ) const;
00145
00146
00147 Real *m_data;
00148 Real *m_hes, *m_hh, *m_d, *m_ee, *m_dd;
00149 T *m_work_arr;
00150 };
00151
00152
00153
00154
00155
00156 template <class T>
00157 void GMRESSolver<T>::allocate()
00158 {
00159 int max_k = m_restrtLen;
00160 int hh = (max_k + 2) * (max_k + 1);
00161 int hes = (max_k + 1) * (max_k + 1);
00162 int rs = (max_k + 2);
00163 int cc = (max_k + 1);
00164 int size = (hh + hes + rs + 2*cc + 1);
00165
00166 m_data = new Real[size];
00167 m_hh = m_data;
00168 m_hes = m_data + hh;
00169 m_d = m_hes + hes;
00170 m_ee = m_d + rs;
00171 m_dd = m_ee + cc;
00172
00173 m_work_arr = 0;
00174 }
00175
00176 template <class T>
00177 void GMRESSolver<T>::clearData()
00178 {
00179 delete [] m_data;
00180 }
00181
00182 template <class T>
00183 void GMRESSolver<T>::setRestartLen(int mm)
00184 {
00185 CH_assert(mm>0);
00186 clearData();
00187 m_restrtLen = mm;
00188 allocate();
00189 }
00190
00191 template <class T>
00192 GMRESSolver<T>::GMRESSolver()
00193 :m_homogeneous(false),
00194 m_op(NULL),
00195 m_restrtLen(30),
00196 m_imax(1000),
00197 m_verbosity(3),
00198 m_eps(1.0E-50),
00199 m_reps(1.0E-12),
00200 m_exitStatus(-1),
00201 m_small(1.0E-30),
00202 m_normType(2)
00203 {
00204 allocate();
00205 }
00206
00207 template <class T>
00208 GMRESSolver<T>::~GMRESSolver()
00209 {
00210 m_op = NULL;
00211 clearData();
00212 }
00213
00214 template <class T>
00215 void GMRESSolver<T>::define(LinearOp<T>* a_operator, bool a_homogeneous)
00216 {
00217 m_homogeneous = a_homogeneous;
00218 m_op = a_operator;
00219 }
00220
00221
00222
00223
00224 #define HH(a,b) (m_hh + (b)*(m_restrtLen+2) + (a))
00225 #define HES(a,b) (m_hes + (b)*(m_restrtLen+1) + (a))
00226 #define CC(a) (m_ee + (a))
00227 #define SS(a) (m_dd + (a))
00228 #define GRS(a) (m_d + (a))
00229
00230
00231 #define VEC_OFFSET 2
00232 #define VEC_TEMP_RHS m_work_arr[0]
00233 #define VEC_TEMP_LHS m_work_arr[1]
00234 #define VEC_VV(i) m_work_arr[VEC_OFFSET + i]
00235
00236
00237
00238
00239 template <class T>
00240 void GMRESSolver<T>::solve( T& a_xx, const T& a_bb )
00241 {
00242 CH_TIMERS("GMRESSolver::solve");
00243
00244 CH_TIMER("GMRESSolver::solve::Initialize",timeInitialize);
00245 CH_TIMER("GMRESSolver::solve::MainLoop",timeMainLoop);
00246 CH_TIMER("GMRESSolver::solve::Cleanup",timeCleanup);
00247
00248 CH_START(timeInitialize);
00249
00250 if (m_verbosity >= 3)
00251 {
00252 pout() << "GMRESSolver::solve" << endl;
00253 }
00254
00255 const int nwork = VEC_OFFSET + m_restrtLen + 1;
00256 m_work_arr = new T[nwork];
00257 m_op->create(VEC_TEMP_RHS, a_bb);
00258 m_op->create(VEC_TEMP_LHS, a_xx);
00259 for (int i=VEC_OFFSET;i<nwork;i++)
00260 {
00261 m_op->create(m_work_arr[i], a_bb);
00262 }
00263
00264 int it;
00265 Real rnorm0 = 0.0;
00266 T &vv_0 = VEC_VV(0);
00267
00268
00269 m_op->assign( vv_0, a_bb );
00270 m_op->setToZero( a_xx );
00271
00272 CH_STOP(timeInitialize);
00273
00274
00275
00276
00277
00278
00279
00280
00281
00282 CH_START(timeMainLoop);
00283
00284
00285 it = 0; m_exitStatus = -1;
00286 CycleGMRES( a_xx, m_exitStatus, it, rnorm0 );
00287
00288 while ( m_exitStatus==-1 && it < m_imax )
00289 {
00290 if (m_verbosity >= 3)
00291 {
00292 pout() << "*";
00293 }
00294 ResidualGMRES( vv_0, a_xx, a_bb, VEC_TEMP_RHS );
00295 CycleGMRES( a_xx, m_exitStatus, it, rnorm0 );
00296 }
00297 if (m_exitStatus==-1 && it >= m_imax)
00298 {
00299 m_exitStatus = 3;
00300 }
00301
00302 CH_STOP(timeMainLoop);
00303
00304 CH_START(timeCleanup);
00305
00306
00307 for (int i=0;i<nwork;i++)
00308 {
00309 m_op->clear(m_work_arr[i]);
00310 }
00311 delete [] m_work_arr; m_work_arr = 0;
00312
00313 if (m_verbosity >= 3)
00314 {
00315 pout() << "GMRESSolver::solve done, status = " << m_exitStatus << endl;
00316 }
00317
00318 CH_STOP(timeCleanup);
00319 }
00320
00321 #define CONVERGED(r0,r) (r<r0*m_reps || r<m_eps)
00322
00323 template <class T>
00324 void GMRESSolver<T>::CycleGMRES( T &a_xx,
00325 int &a_reason, int &a_itcount, Real &a_rnorm0,
00326 const bool a_avoidnorms )
00327 {
00328 CH_assert(m_op != 0);
00329 Real res_norm,res,hapbnd,tt;
00330 int it;
00331 bool hapend = false;
00332 T &vv_0 = VEC_VV(0);
00333
00334
00335
00336 res_norm = m_op->norm( vv_0, m_normType );
00337 res = res_norm;
00338 *GRS(0) = res_norm;
00339
00340
00341 if ( res == 0. )
00342 {
00343 if (m_verbosity >= 3)
00344 {
00345 pout() << "GMRESSolver::solve zero residual!!!" << endl;
00346 }
00347 a_reason = 1;
00348 return;
00349 }
00350
00351 tt = 1./res;
00352 m_op->scale( vv_0, tt);
00353
00354 if ( a_itcount == 0 ) a_rnorm0 = res;
00355 bool conv = CONVERGED( a_rnorm0, res );
00356 a_reason = conv ? 1 : -1;
00357 it = 0;
00358 while ( a_reason == -1 && it < m_restrtLen && a_itcount < m_imax )
00359 {
00360 if ( m_verbosity>=4 && (it!=0 || a_itcount==0) )
00361 {
00362 pout() << a_itcount << ") GMRES residual = " << res << endl;
00363 }
00364
00365 const T &vv_it = VEC_VV(it);
00366 T &vv_it1 = VEC_VV(it+1);
00367
00368
00369 {
00370 T &Mb = VEC_TEMP_LHS;
00371 ApplyAB( vv_it1, vv_it, Mb );
00372 }
00373
00374 TwoUnmodifiedGramSchmidtOrthogonalization(it);
00375
00376
00377
00378 tt = m_op->norm( vv_it1, m_normType );
00379
00380
00381
00382 hapbnd = 1.e-99;
00383 if (tt < hapbnd)
00384 {
00385 if ( m_verbosity>=3 )
00386 {
00387 pout() << "Detected happy breakdown, "<<it<<") current hapbnd = "<< tt << endl;
00388 }
00389 hapend = true;
00390 }
00391 else
00392 {
00393 m_op->scale( vv_it1, 1./tt);
00394 }
00395
00396
00397 *HH(it+1,it) = tt; *HES(it+1,it) = tt;
00398
00399 UpdateGMRESHessenberg( it, hapend, res );
00400
00401
00402 it++; (a_itcount)++;
00403 conv = CONVERGED( a_rnorm0, res );
00404 a_reason = conv ? 1 : -1;
00405
00406 if ( hapend )
00407 {
00408 if ( !conv && m_verbosity>=3)
00409 {
00410 pout() << "You reached the happy break down, but convergence was not indicated. Residual norm=" << res << endl;
00411 }
00412 break;
00413 }
00414 }
00415
00416 if ( (a_reason!=1 || a_itcount>=m_imax) && m_verbosity>=4 )
00417 {
00418 pout() << a_itcount << ") GMRES residual = " << res << endl;
00419 }
00420
00421
00422
00423
00424
00425
00426
00427 BuildGMRESSoln( GRS(0), a_xx, it-1, &VEC_VV(0) );
00428 }
00429
00430 template <class T>
00431 void GMRESSolver<T>::ResidualGMRES( T &a_vv, const T &a_xx,
00432 const T &a_bb, T &a_temp_rhs )
00433 {
00434 CH_assert(m_op != 0);
00435 m_op->applyOp( a_temp_rhs, a_xx, true);
00436 m_op->assign( a_vv, a_bb );
00437 m_op->incr( a_vv, a_temp_rhs, -1.0);
00438 }
00439 template <class T>
00440 void GMRESSolver<T>::BuildGMRESSoln( Real nrs[], T &a_xx, const int it,
00441 const T vv_0[] )
00442 {
00443 Real tt;
00444 int ii,k,j;
00445
00446
00447
00448
00449 if (it < 0)
00450 {
00451 return;
00452 }
00453 if (*HH(it,it) == 0.0)
00454 {
00455 pout() << "HH(it,it) is identically zero!!! GRS(it) = " << *GRS(it) << endl;
00456 }
00457 if (*HH(it,it) != 0.0)
00458 {
00459 nrs[it] = *GRS(it) / *HH(it,it);
00460 }
00461 else
00462 {
00463 nrs[it] = 0.0;
00464 }
00465
00466 for (ii=1; ii<=it; ii++)
00467 {
00468 k = it - ii;
00469 tt = *GRS(k);
00470 for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
00471 nrs[k] = tt / *HH(k,k);
00472 }
00473
00474
00475 T &temp = VEC_TEMP_RHS;
00476 m_op->setToZero(temp);
00477
00478 for (ii=0; ii<it+1; ii++)
00479 {
00480 m_op->incr(temp, vv_0[ii], nrs[ii]);
00481 }
00482
00483
00484
00485 T &temp_matop = VEC_TEMP_LHS;
00486
00487 m_op->preCond( temp_matop, temp );
00488 m_op->incr( a_xx, temp_matop, 1.0 );
00489 }
00490
00491
00492
00493
00494
00495
00496
00497
00498
00499
00500
00501
00502
00503 template <class T>
00504 void GMRESSolver<T>::UpdateGMRESHessenberg( const int it, bool hapend, Real &res )
00505 {
00506 Real *hh,*cc,*ss,tt;
00507 int j;
00508
00509 hh = HH(0,it);
00510 cc = CC(0);
00511 ss = SS(0);
00512
00513
00514
00515 for (j=1; j<=it; j++)
00516 {
00517 tt = *hh;
00518 *hh = *cc * tt + *ss * *(hh+1);
00519 hh++;
00520 *hh = *cc++ * *hh - (*ss++ * tt);
00521 }
00522
00523
00524
00525
00526
00527
00528
00529 if ( !hapend )
00530 {
00531 tt = sqrt( *hh * *hh + *(hh+1) * *(hh+1) );
00532 if (tt == 0.0)
00533 {
00534 pout() << "Your matrix or preconditioner is the null operator\n";
00535 }
00536 *cc = *hh / tt;
00537 *ss = *(hh+1) / tt;
00538 *GRS(it+1) = - (*ss * *GRS(it));
00539 *GRS(it) = *cc * *GRS(it);
00540 *hh = *cc * *hh + *ss * *(hh+1);
00541 res = Abs( *GRS(it+1) );
00542 }
00543 else
00544 {
00545
00546
00547
00548
00549
00550
00551
00552 res = 0.0;
00553 }
00554 }
00555
00556
00557
00558
00559
00560
00561
00562
00563
00564
00565
00566
00567
00568
00569
00570
00571
00572
00573
00574
00575
00576
00577
00578
00579
00580
00581
00582
00583
00584
00585
00586
00587
00588
00589
00590
00591
00592
00593
00594
00595
00596
00597
00598
00599
00600
00601
00602
00603
00604
00605
00606
00607
00608
00609
00610
00611
00612
00613
00614
00615
00616
00617
00618
00619
00620
00621
00622
00623
00624
00625 template <class T>
00626 void GMRESSolver<T>::TwoUnmodifiedGramSchmidtOrthogonalization( const int it )
00627 {
00628 Real *hh,*hes,*lhh = 0;
00629 T &vv_1 = VEC_VV(it+1);
00630 const T *vv_0 = &(VEC_VV(0));
00631
00632
00633 lhh = new Real[it+1];
00634
00635
00636 hh = HH(0,it);
00637 hes = HES(0,it);
00638
00639
00640 for (int j=0; j<=it; j++)
00641 {
00642 hh[j] = 0.0;
00643 hes[j] = 0.0;
00644 }
00645
00646 for ( int ncnt = 0 ; ncnt < 2 ; ncnt++ )
00647 {
00648
00649
00650
00651
00652
00653
00654
00655
00656
00657 m_op->mDotProduct(vv_1, it+1, vv_0, lhh);
00658
00659
00660
00661
00662
00663 for (int j=0; j<=it; j++) lhh[j] = - lhh[j];
00664
00665 for (int j=0; j<=it; j++)
00666 {
00667 m_op->incr(vv_1, vv_0[j], lhh[j]);
00668 }
00669 for (int j=0; j<=it; j++)
00670 {
00671 hh[j] -= lhh[j];
00672 hes[j] += lhh[j];
00673 }
00674 }
00675
00676 delete [] lhh;
00677 }
00678
00679
00680
00681
00682 template <class T>
00683 void GMRESSolver<T>::ApplyAB( T &a_dest, const T &a_xx, T &a_tmp_lhs ) const
00684 {
00685 m_op->preCond( a_tmp_lhs, a_xx );
00686 m_op->applyOp( a_dest, a_tmp_lhs, true);
00687 }
00688
00689 template <class T>
00690 void GMRESSolver<T>::setConvergenceMetrics(Real a_metric,
00691 Real a_tolerance)
00692 {
00693 m_eps = a_tolerance;
00694 }
00695
00696 #include "NamespaceFooter.H"
00697 #endif