Actual source code: minres.c

  1: /*$Id: minres.c,v 1.18 2001/08/10 18:14:38 bsmith Exp $*/
  2: /*                       
  3:     This code implements the MINRES (Minimum Residual) method. 
  4:     Reference: Paige & Saunders, 1975.

  6:     Contributed by: Robert Scheichl: maprs@maths.bath.ac.uk

  8: */
 9:  #include src/sles/ksp/kspimpl.h

 11: typedef struct {
 12:   PetscReal haptol;
 13: } KSP_MINRES;

 15: int KSPSetUp_MINRES(KSP ksp)
 16: {


 21:   if (ksp->pc_side == PC_RIGHT) {
 22:     SETERRQ(2,"No right preconditioning for KSPMINRES");
 23:   } else if (ksp->pc_side == PC_SYMMETRIC) {
 24:     SETERRQ(2,"No symmetric preconditioning for KSPMINRES");
 25:   }

 27:   KSPDefaultGetWork(ksp,9);

 29:   return(0);
 30: }


 33: int  KSPSolve_MINRES(KSP ksp,int *its)
 34: {
 35:   int          ierr,i,maxit;
 36:   PetscScalar  alpha,malpha,beta,mbeta,ibeta,betaold,eta,c=1.0,ceta,cold=1.0,coold,s=0.0,sold=0.0,soold;
 37:   PetscScalar  rho0,rho1,irho1,rho2,mrho2,rho3,mrho3,mone = -1.0,zero = 0.0,dp = 0.0;
 38:   PetscReal    np;
 39:   Vec          X,B,R,Z,U,V,W,UOLD,VOLD,WOLD,WOOLD;
 40:   Mat          Amat,Pmat;
 41:   MatStructure pflag;
 42:   KSP_MINRES   *minres = (KSP_MINRES*)ksp->data;
 43:   PetscTruth   diagonalscale;

 46:   ierr    = PCDiagonalScale(ksp->B,&diagonalscale);
 47:   if (diagonalscale) SETERRQ1(1,"Krylov method %s does not support diagonal scaling",ksp->type_name);

 49:   maxit   = ksp->max_it;
 50:   X       = ksp->vec_sol;
 51:   B       = ksp->vec_rhs;
 52:   R       = ksp->work[0];
 53:   Z       = ksp->work[1];
 54:   U       = ksp->work[2];
 55:   V       = ksp->work[3];
 56:   W       = ksp->work[4];
 57:   UOLD    = ksp->work[5];
 58:   VOLD    = ksp->work[6];
 59:   WOLD    = ksp->work[7];
 60:   WOOLD   = ksp->work[8];

 62:   PCGetOperators(ksp->B,&Amat,&Pmat,&pflag);

 64:   ksp->its = 0;

 66:   VecSet(&zero,UOLD);         /*     u_old  <-   0   */
 67:   VecCopy(UOLD,VOLD);         /*     v_old  <-   0   */
 68:   VecCopy(UOLD,W);            /*     w      <-   0   */
 69:   VecCopy(UOLD,WOLD);         /*     w_old  <-   0   */

 71:   if (!ksp->guess_zero) {
 72:     KSP_MatMult(ksp,Amat,X,R); /*     r <- b - A*x    */
 73:     VecAYPX(&mone,B,R);
 74:   } else {
 75:     VecCopy(B,R);              /*     r <- b (x is 0) */
 76:   }

 78:   KSP_PCApply(ksp,ksp->B,R,Z); /*     z  <- B*r       */

 80:   VecDot(R,Z,&dp);
 81:   if (PetscAbsScalar(dp) < minres->haptol) {
 82:     PetscLogInfo(ksp,"KSPSolve_MINRES:Detected happy breakdown %g tolerance %gn",PetscAbsScalar(dp),minres->haptol);
 83:     dp = PetscAbsScalar(dp); /* tiny number, can't use 0.0, cause divided by below */
 84:     if (dp == 0.0) {
 85:       ksp->reason = KSP_CONVERGED_ATOL;
 86:       *its        = 0;
 87:       return(0);
 88:     }
 89:   }

 91: #if !defined(PETSC_USE_COMPLEX)
 92:   if (dp < 0.0) SETERRQ(PETSC_ERR_KSP_BRKDWN,"Indefinite preconditioner");
 93: #endif
 94:   dp   = PetscSqrtScalar(dp);
 95:   beta = dp;                                        /*  beta <- sqrt(r'*z  */
 96:   eta  = beta;

 98:   VecCopy(R,V);
 99:   VecCopy(Z,U);
100:   ibeta = 1.0 / beta;
101:   VecScale(&ibeta,V);         /*    v <- r / beta     */
102:   VecScale(&ibeta,U);         /*    u <- z / beta     */

104:   VecNorm(Z,NORM_2,&np);      /*   np <- ||z||        */

106:   (*ksp->converged)(ksp,0,np,&ksp->reason,ksp->cnvP);  /* test for convergence */
107:   if (ksp->reason) {*its =  0; return(0);}
108:   KSPLogResidualHistory(ksp,np);
109:   KSPMonitor(ksp,0,np);            /* call any registered monitor routines */
110:   ksp->rnorm = np;

112:   for (i=0; i<maxit; i++) {
113:      ksp->its = i+1;

115: /*   Lanczos  */

117:      KSP_MatMult(ksp,Amat,U,R);   /*      r <- A*u   */
118:      VecDot(U,R,&alpha);          /*  alpha <- r'*u  */
119:      KSP_PCApply(ksp,ksp->B,R,Z); /*      z <- B*r   */

121:      malpha = - alpha;
122:      VecAXPY(&malpha,V,R);     /*  r <- r - alpha v     */
123:      VecAXPY(&malpha,U,Z);     /*  z <- z - alpha u     */
124:      mbeta = - beta;
125:      VecAXPY(&mbeta,VOLD,R);   /*  r <- r - beta v_old  */
126:      VecAXPY(&mbeta,UOLD,Z);   /*  z <- z - beta u_old  */

128:      betaold = beta;

130:      VecDot(R,Z,&dp);
131:      if (PetscAbsScalar(dp) < minres->haptol) {
132:        PetscLogInfo(ksp,"KSPSolve_MINRES:Detected happy breakdown %g tolerance %gn",PetscAbsScalar(dp),minres->haptol);
133:        dp = PetscAbsScalar(dp); /* tiny number, can we use 0.0? */
134:      }

136: #if !defined(PETSC_USE_COMPLEX)
137:      if (dp < 0.0) SETERRQ1(PETSC_ERR_KSP_BRKDWN,"Indefinite preconditioner R'Z = %g",dp);
138: #endif
139:      beta = PetscSqrtScalar(dp);                               /*  beta <- sqrt(r'*z)   */

141: /*    QR factorisation    */

143:      coold = cold; cold = c; soold = sold; sold = s;

145:      rho0 = cold * alpha - coold * sold * betaold;
146:      rho1 = PetscSqrtScalar(rho0*rho0 + beta*beta);
147:      rho2 = sold * alpha + coold * cold * betaold;
148:      rho3 = soold * betaold;

150: /*     Givens rotation    */

152:      c = rho0 / rho1;
153:      s = beta / rho1;

155: /*    Update    */

157:      VecCopy(WOLD,WOOLD);     /*  w_oold <- w_old      */
158:      VecCopy(W,WOLD);         /*  w_old  <- w          */
159: 
160:      VecCopy(U,W);            /*  w      <- u          */
161:      mrho2 = - rho2;
162:      VecAXPY(&mrho2,WOLD,W);  /*  w <- w - rho2 w_old  */
163:      mrho3 = - rho3;
164:      VecAXPY(&mrho3,WOOLD,W); /*  w <- w - rho3 w_oold */
165:      irho1 = 1.0 / rho1;
166:      VecScale(&irho1,W);      /*  w <- w / rho1        */

168:      ceta = c * eta;
169:      VecAXPY(&ceta,W,X);      /*  x <- x + c eta w     */
170:      eta = - s * eta;

172:      VecCopy(V,VOLD);
173:      VecCopy(U,UOLD);
174:      VecCopy(R,V);
175:      VecCopy(Z,U);
176:      ibeta = 1.0 / beta;
177:      VecScale(&ibeta,V);      /*  v <- r / beta       */
178:      VecScale(&ibeta,U);      /*  u <- z / beta       */
179: 
180:      np = ksp->rnorm * PetscAbsScalar(s);

182:      ksp->rnorm = np;
183:      KSPLogResidualHistory(ksp,np);
184:      KSPMonitor(ksp,i+1,np);
185:      (*ksp->converged)(ksp,i+1,np,&ksp->reason,ksp->cnvP); /* test for convergence */
186:      if (ksp->reason) break;
187:   }
188:   if (i == maxit) {
189:     ksp->reason = KSP_DIVERGED_ITS;
190:   }
191:   *its = ksp->its;
192:   return(0);
193: }

195: EXTERN_C_BEGIN
196: int KSPCreate_MINRES(KSP ksp)
197: {
198:   KSP_MINRES *minres;


203:   ksp->pc_side   = PC_LEFT;
204:   ierr           = PetscNew(KSP_MINRES,&minres);
205:   minres->haptol = 1.e-18;
206:   ksp->data      = (void*)minres;

208:   /*
209:        Sets the functions that are associated with this data structure 
210:        (in C++ this is the same as defining virtual functions)
211:   */
212:   ksp->ops->setup                = KSPSetUp_MINRES;
213:   ksp->ops->solve                = KSPSolve_MINRES;
214:   ksp->ops->destroy              = KSPDefaultDestroy;
215:   ksp->ops->setfromoptions       = 0;
216:   ksp->ops->buildsolution        = KSPDefaultBuildSolution;
217:   ksp->ops->buildresidual        = KSPDefaultBuildResidual;

219:   return(0);
220: }
221: EXTERN_C_END