Actual source code: tcqmr.c

  1: #define PETSCKSP_DLL

  3: /*
  4:     This file contains an implementation of Tony Chan's transpose-free QMR.

  6:     Note: The vector dot products in the code have not been checked for the
  7:     complex numbers version, so most probably some are incorrect.
  8: */

 10:  #include include/private/kspimpl.h
 11:  #include src/ksp/ksp/impls/tcqmr/tcqmrp.h

 15: static PetscErrorCode KSPSolve_TCQMR(KSP ksp)
 16: {
 17:   PetscReal      rnorm0,rnorm,dp1,Gamma;
 18:   PetscScalar    theta,ep,cl1,sl1,cl,sl,sprod,tau_n1,f;
 19:   PetscScalar    deltmp,rho,beta,eptmp,ta,s,c,tau_n,delta;
 20:   PetscScalar    dp11,dp2,rhom1,alpha,tmp;

 24:   ksp->its = 0;

 26:   KSPInitialResidual(ksp,x,u,v,r,b);
 27:   VecNorm(r,NORM_2,&rnorm0);         /*  rnorm0 = ||r|| */

 29:   (*ksp->converged)(ksp,0,rnorm0,&ksp->reason,ksp->cnvP);
 30:   if (ksp->reason) return(0);

 32:   VecSet(um1,0.0);
 33:   VecCopy(r,u);
 34:   rnorm = rnorm0;
 35:   tmp   = 1.0/rnorm; VecScale(u,tmp);
 36:   VecSet(vm1,0.0);
 37:   VecCopy(u,v);
 38:   VecCopy(u,v0);
 39:   VecSet(pvec1,0.0);
 40:   VecSet(pvec2,0.0);
 41:   VecSet(p,0.0);
 42:   theta = 0.0;
 43:   ep    = 0.0;
 44:   cl1   = 0.0;
 45:   sl1   = 0.0;
 46:   cl    = 0.0;
 47:   sl    = 0.0;
 48:   sprod = 1.0;
 49:   tau_n1= rnorm0;
 50:   f     = 1.0;
 51:   Gamma = 1.0;
 52:   rhom1 = 1.0;

 54:   /*
 55:    CALCULATE SQUARED LANCZOS  vectors
 56:    */
 57:   (*ksp->converged)(ksp,ksp->its,rnorm,&ksp->reason,ksp->cnvP);
 58:   while (!ksp->reason){
 59:     KSPMonitor(ksp,ksp->its,rnorm);
 60:     ksp->its++;

 62:     KSP_PCApplyBAorAB(ksp,u,y,vtmp); /* y = A*u */
 63:     VecDot(v0,y,&dp11);
 64:     VecDot(v0,u,&dp2);
 65:     alpha  = dp11 / dp2;                          /* alpha = v0'*y/v0'*u */
 66:     deltmp = alpha;
 67:     VecCopy(y,z);
 68:     VecAXPY(z,-alpha,u); /* z = y - alpha u */
 69:     VecDot(v0,u,&rho);
 70:     beta   = rho / (f*rhom1);
 71:     rhom1  = rho;
 72:     VecCopy(z,utmp);    /* up1 = (A-alpha*I)*
 73:                                                  (z-2*beta*p) + f*beta*
 74:                                                  beta*um1 */
 75:     VecAXPY(utmp,-2.0*beta,p);
 76:     KSP_PCApplyBAorAB(ksp,utmp,up1,vtmp);
 77:     VecAXPY(up1,-alpha,utmp);
 78:     VecAXPY(up1,f*beta*beta,um1);
 79:     VecNorm(up1,NORM_2,&dp1);
 80:     f      = 1.0 / dp1;
 81:     VecScale(up1,f);
 82:     VecAYPX(p,-beta,z);   /* p = f*(z-beta*p) */
 83:     VecScale(p,f);
 84:     VecCopy(u,um1);
 85:     VecCopy(up1,u);
 86:     beta   = beta/Gamma;
 87:     eptmp  = beta;
 88:     KSP_PCApplyBAorAB(ksp,v,vp1,vtmp);
 89:     VecAXPY(vp1,-alpha,v);
 90:     VecAXPY(vp1,-beta,vm1);
 91:     VecNorm(vp1,NORM_2,&Gamma);
 92:     VecScale(vp1,1.0/Gamma);
 93:     VecCopy(v,vm1);
 94:     VecCopy(vp1,v);

 96:   /*
 97:      SOLVE  Ax = b
 98:    */
 99:   /* Apply last two Given's (Gl-1 and Gl) rotations to (beta,alpha,Gamma) */
100:     if (ksp->its > 2) {
101:       theta =  sl1*beta;
102:       eptmp = -cl1*beta;
103:     }
104:     if (ksp->its > 1) {
105:       ep     = -cl*eptmp + sl*alpha;
106:       deltmp = -sl*eptmp - cl*alpha;
107:     }
108:     if (PetscAbsReal(Gamma) > PetscAbsScalar(deltmp)) {
109:       ta = -deltmp / Gamma;
110:       s  = 1.0 / PetscSqrtScalar(1.0 + ta*ta);
111:       c  = s*ta;
112:     } else {
113:       ta = -Gamma/deltmp;
114:       c  = 1.0 / PetscSqrtScalar(1.0 + ta*ta);
115:       s  = c*ta;
116:     }

118:     delta  = -c*deltmp + s*Gamma;
119:     tau_n  = -c*tau_n1; tau_n1 = -s*tau_n1;
120:     VecCopy(vm1,pvec);
121:     VecAXPY(pvec,-theta,pvec2);
122:     VecAXPY(pvec,-ep,pvec1);
123:     VecScale(pvec,1.0/delta);
124:     VecAXPY(x,tau_n,pvec);
125:     cl1    = cl; sl1 = sl; cl = c; sl = s;

127:     VecCopy(pvec1,pvec2);
128:     VecCopy(pvec,pvec1);

130:     /* Compute the upper bound on the residual norm r (See QMR paper p. 13) */
131:     sprod = sprod*PetscAbsScalar(s);
132: #if defined(PETSC_USE_COMPLEX)
133:     rnorm = rnorm0 * sqrt((double)ksp->its+2.0) * PetscRealPart(sprod);
134: #else
135:     rnorm = rnorm0 * sqrt((double)ksp->its+2.0) * sprod;
136: #endif
137:     if (ksp->its >= ksp->max_it) {ksp->reason = KSP_DIVERGED_ITS; break;}
138:     (*ksp->converged)(ksp,ksp->its,rnorm,&ksp->reason,ksp->cnvP);
139:   }
140:   KSPMonitor(ksp,ksp->its,rnorm);
141:   KSPUnwindPreconditioner(ksp,x,vtmp);

143:   return(0);
144: }

148: static PetscErrorCode KSPSetUp_TCQMR(KSP ksp)
149: {

153:   if (ksp->pc_side == PC_SYMMETRIC){
154:     SETERRQ(PETSC_ERR_SUP,"no symmetric preconditioning for KSPTCQMR");
155:   }
156:   KSPDefaultGetWork(ksp,TCQMR_VECS);
157:   return(0);
158: }

160: /*MC
161:      KSPRTCQMR - A variant of QMR (quasi minimal residual) developed by Tony Chan

163:    Options Database Keys:
164: .   see KSPSolve()

166:    Level: beginner

168: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPTFQMR

170: M*/

175: PetscErrorCode  KSPCreate_TCQMR(KSP ksp)
176: {
178:   ksp->data                = (void*)0;
179:   ksp->pc_side             = PC_LEFT;
180:   ksp->ops->buildsolution  = KSPDefaultBuildSolution;
181:   ksp->ops->buildresidual  = KSPDefaultBuildResidual;
182:   ksp->ops->setup          = KSPSetUp_TCQMR;
183:   ksp->ops->solve          = KSPSolve_TCQMR;
184:   ksp->ops->destroy        = KSPDefaultDestroy;
185:   ksp->ops->setfromoptions = 0;
186:   ksp->ops->view           = 0;
187:   return(0);
188: }