Actual source code: minres.c

  1: /*$Id: minres.c,v 1.18 2001/08/10 18:14:38 bsmith Exp $*/
 2:  #include src/ksp/ksp/kspimpl.h

  4: typedef struct {
  5:   PetscReal haptol;
  6: } KSP_MINRES;

 10: int KSPSetUp_MINRES(KSP ksp)
 11: {


 16:   if (ksp->pc_side == PC_RIGHT) {
 17:     SETERRQ(2,"No right preconditioning for KSPMINRES");
 18:   } else if (ksp->pc_side == PC_SYMMETRIC) {
 19:     SETERRQ(2,"No symmetric preconditioning for KSPMINRES");
 20:   }

 22:   KSPDefaultGetWork(ksp,9);

 24:   return(0);
 25: }


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

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

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

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

 60:   ksp->its = 0;

 62:   VecSet(&zero,UOLD);         /*     u_old  <-   0   */
 63:   VecCopy(UOLD,VOLD);         /*     v_old  <-   0   */
 64:   VecCopy(UOLD,W);            /*     w      <-   0   */
 65:   VecCopy(UOLD,WOLD);         /*     w_old  <-   0   */

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

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

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

 86: #if !defined(PETSC_USE_COMPLEX)
 87:   if (dp < 0.0) SETERRQ(PETSC_ERR_KSP_BRKDWN,"Indefinite preconditioner");
 88: #endif
 89:   dp   = PetscSqrtScalar(dp);
 90:   beta = dp;                                        /*  beta <- sqrt(r'*z  */
 91:   eta  = beta;

 93:   VecCopy(R,V);
 94:   VecCopy(Z,U);
 95:   ibeta = 1.0 / beta;
 96:   VecScale(&ibeta,V);         /*    v <- r / beta     */
 97:   VecScale(&ibeta,U);         /*    u <- z / beta     */

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

101:   KSPLogResidualHistory(ksp,np);
102:   KSPMonitor(ksp,0,np);            /* call any registered monitor routines */
103:   ksp->rnorm = np;
104:   (*ksp->converged)(ksp,0,np,&ksp->reason,ksp->cnvP);  /* test for convergence */
105:   if (ksp->reason) return(0);

107:   i = 0;
108:   do {
109:      ksp->its = i+1;

111: /*   Lanczos  */

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

117:      malpha = - alpha;
118:      VecAXPY(&malpha,V,R);     /*  r <- r - alpha v     */
119:      VecAXPY(&malpha,U,Z);     /*  z <- z - alpha u     */
120:      mbeta = - beta;
121:      VecAXPY(&mbeta,VOLD,R);   /*  r <- r - beta v_old  */
122:      VecAXPY(&mbeta,UOLD,Z);   /*  z <- z - beta u_old  */

124:      betaold = beta;

126:      VecDot(R,Z,&dp);
127:      if (PetscAbsScalar(dp) < minres->haptol) {
128:        PetscLogInfo(ksp,"KSPSolve_MINRES:Detected happy breakdown %g tolerance %g\n",PetscAbsScalar(dp),minres->haptol);
129:        dp = PetscAbsScalar(dp); /* tiny number, can we use 0.0? */
130:      }

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

137: /*    QR factorisation    */

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

141:      rho0 = cold * alpha - coold * sold * betaold;
142:      rho1 = PetscSqrtScalar(rho0*rho0 + beta*beta);
143:      rho2 = sold * alpha + coold * cold * betaold;
144:      rho3 = soold * betaold;

146: /*     Givens rotation    */

148:      c = rho0 / rho1;
149:      s = beta / rho1;

151: /*    Update    */

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

164:      ceta = c * eta;
165:      VecAXPY(&ceta,W,X);      /*  x <- x + c eta w     */
166:      eta = - s * eta;

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

178:      ksp->rnorm = np;
179:      KSPLogResidualHistory(ksp,np);
180:      KSPMonitor(ksp,i+1,np);
181:      (*ksp->converged)(ksp,i+1,np,&ksp->reason,ksp->cnvP); /* test for convergence */
182:      if (ksp->reason) break;
183:      i++;
184:   } while (i<ksp->max_it);
185:   if (i == ksp->max_it) {
186:     ksp->reason = KSP_DIVERGED_ITS;
187:   }
188:   return(0);
189: }

191: /*MC
192:      KSPMINRES -     This code implements the MINRES (Minimum Residual) method. 
193:                  Reference: Paige & Saunders, 1975.

195:    Options Database Keys:
196: .   see KSPSolve()

198:    Level: beginner

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

202:    Notes: The operator and the preconditioner must be positive definite for this method

204: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPCG

206: M*/

208: EXTERN_C_BEGIN
211: int KSPCreate_MINRES(KSP ksp)
212: {
213:   KSP_MINRES *minres;


218:   ksp->pc_side   = PC_LEFT;
219:   PetscNew(KSP_MINRES,&minres);
220:   minres->haptol = 1.e-18;
221:   ksp->data      = (void*)minres;

223:   /*
224:        Sets the functions that are associated with this data structure 
225:        (in C++ this is the same as defining virtual functions)
226:   */
227:   ksp->ops->setup                = KSPSetUp_MINRES;
228:   ksp->ops->solve                = KSPSolve_MINRES;
229:   ksp->ops->destroy              = KSPDefaultDestroy;
230:   ksp->ops->setfromoptions       = 0;
231:   ksp->ops->buildsolution        = KSPDefaultBuildSolution;
232:   ksp->ops->buildresidual        = KSPDefaultBuildResidual;

234:   return(0);
235: }
236: EXTERN_C_END