Actual source code: minres.c

petsc-3.4.5 2014-06-29
  2: #include <petsc-private/kspimpl.h>

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

 10: PetscErrorCode KSPSetUp_MINRES(KSP ksp)
 11: {

 15:   if (ksp->pc_side == PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"No right preconditioning for KSPMINRES");
 16:   else if (ksp->pc_side == PC_SYMMETRIC) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"No symmetric preconditioning for KSPMINRES");
 17:   KSPSetWorkVecs(ksp,9);
 18:   return(0);
 19: }


 24: PetscErrorCode  KSPSolve_MINRES(KSP ksp)
 25: {
 27:   PetscInt       i;
 28:   PetscScalar    alpha,beta,ibeta,betaold,eta,c=1.0,ceta,cold=1.0,coold,s=0.0,sold=0.0,soold;
 29:   PetscScalar    rho0,rho1,irho1,rho2,mrho2,rho3,mrho3,dp = 0.0;
 30:   PetscReal      np;
 31:   Vec            X,B,R,Z,U,V,W,UOLD,VOLD,WOLD,WOOLD;
 32:   Mat            Amat,Pmat;
 33:   MatStructure   pflag;
 34:   KSP_MINRES     *minres = (KSP_MINRES*)ksp->data;
 35:   PetscBool      diagonalscale;

 38:   PCGetDiagonalScale(ksp->pc,&diagonalscale);
 39:   if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);

 41:   X     = ksp->vec_sol;
 42:   B     = ksp->vec_rhs;
 43:   R     = ksp->work[0];
 44:   Z     = ksp->work[1];
 45:   U     = ksp->work[2];
 46:   V     = ksp->work[3];
 47:   W     = ksp->work[4];
 48:   UOLD  = ksp->work[5];
 49:   VOLD  = ksp->work[6];
 50:   WOLD  = ksp->work[7];
 51:   WOOLD = ksp->work[8];

 53:   PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);

 55:   ksp->its = 0;

 57:   VecSet(UOLD,0.0);          /*     u_old  <-   0   */
 58:   VecCopy(UOLD,VOLD);         /*     v_old  <-   0   */
 59:   VecCopy(UOLD,W);            /*     w      <-   0   */
 60:   VecCopy(UOLD,WOLD);         /*     w_old  <-   0   */

 62:   if (!ksp->guess_zero) {
 63:     KSP_MatMult(ksp,Amat,X,R); /*     r <- b - A*x    */
 64:     VecAYPX(R,-1.0,B);
 65:   } else {
 66:     VecCopy(B,R);              /*     r <- b (x is 0) */
 67:   }

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

 71:   VecDot(R,Z,&dp);
 72:   if (PetscRealPart(dp) < minres->haptol) {
 73:     PetscInfo2(ksp,"Detected indefinite operator %G tolerance %G\n",PetscRealPart(dp),minres->haptol);
 74:     ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
 75:     return(0);
 76:   }

 78:   dp   = PetscAbsScalar(dp);
 79:   dp   = PetscSqrtScalar(dp);
 80:   beta = dp;                                        /*  beta <- sqrt(r'*z  */
 81:   eta  = beta;

 83:   VecCopy(R,V);
 84:   VecCopy(Z,U);
 85:   ibeta = 1.0 / beta;
 86:   VecScale(V,ibeta);        /*    v <- r / beta     */
 87:   VecScale(U,ibeta);        /*    u <- z / beta     */

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

 91:   KSPLogResidualHistory(ksp,np);
 92:   KSPMonitor(ksp,0,np);
 93:   ksp->rnorm = np;
 94:   (*ksp->converged)(ksp,0,np,&ksp->reason,ksp->cnvP); /* test for convergence */
 95:   if (ksp->reason) return(0);

 97:   i = 0;
 98:   do {
 99:     ksp->its = i+1;

101: /*   Lanczos  */

103:     KSP_MatMult(ksp,Amat,U,R);   /*      r <- A*u   */
104:     VecDot(U,R,&alpha);          /*  alpha <- r'*u  */
105:     KSP_PCApply(ksp,R,Z); /*      z <- B*r   */

107:     VecAXPY(R,-alpha,V);     /*  r <- r - alpha v     */
108:     VecAXPY(Z,-alpha,U);     /*  z <- z - alpha u     */
109:     VecAXPY(R,-beta,VOLD);   /*  r <- r - beta v_old  */
110:     VecAXPY(Z,-beta,UOLD);   /*  z <- z - beta u_old  */

112:     betaold = beta;

114:     VecDot(R,Z,&dp);
115:     if ( PetscRealPart(dp) < minres->haptol) {
116:       PetscInfo2(ksp,"Detected indefinite operator %G tolerance %G\n",PetscRealPart(dp),minres->haptol);
117:       ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
118:       break;
119:     }

121:     dp   = PetscAbsScalar(dp);
122:     beta = PetscSqrtScalar(dp);                               /*  beta <- sqrt(r'*z)   */

124: /*    QR factorisation    */

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

128:     rho0 = cold * alpha - coold * sold * betaold;
129:     rho1 = PetscSqrtScalar(rho0*rho0 + beta*beta);
130:     rho2 = sold * alpha + coold * cold * betaold;
131:     rho3 = soold * betaold;

133: /*     Givens rotation    */

135:     c = rho0 / rho1;
136:     s = beta / rho1;

138: /*    Update    */

140:     VecCopy(WOLD,WOOLD);     /*  w_oold <- w_old      */
141:     VecCopy(W,WOLD);         /*  w_old  <- w          */

143:     VecCopy(U,W);           /*  w      <- u          */
144:     mrho2 = -rho2;
145:     VecAXPY(W,mrho2,WOLD); /*  w <- w - rho2 w_old  */
146:     mrho3 = -rho3;
147:     VecAXPY(W,mrho3,WOOLD); /*  w <- w - rho3 w_oold */
148:     irho1 = 1.0 / rho1;
149:     VecScale(W,irho1);     /*  w <- w / rho1        */

151:     ceta = c * eta;
152:     VecAXPY(X,ceta,W);      /*  x <- x + c eta w     */
153:     eta  = -s * eta;

155:     VecCopy(V,VOLD);
156:     VecCopy(U,UOLD);
157:     VecCopy(R,V);
158:     VecCopy(Z,U);
159:     ibeta = 1.0 / beta;
160:     VecScale(V,ibeta);     /*  v <- r / beta       */
161:     VecScale(U,ibeta);     /*  u <- z / beta       */

163:     np = ksp->rnorm * PetscAbsScalar(s);

165:     ksp->rnorm = np;
166:     KSPLogResidualHistory(ksp,np);
167:     KSPMonitor(ksp,i+1,np);
168:     (*ksp->converged)(ksp,i+1,np,&ksp->reason,ksp->cnvP); /* test for convergence */
169:     if (ksp->reason) break;
170:     i++;
171:   } while (i<ksp->max_it);
172:   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
173:   return(0);
174: }

176: /*MC
177:      KSPMINRES - This code implements the MINRES (Minimum Residual) method.

179:    Options Database Keys:
180: .   see KSPSolve()

182:    Level: beginner

184:    Notes: The operator and the preconditioner must be symmetric and the preconditioner must
185:           be positive definite for this method.
186:           Supports only left preconditioning.

188:    Reference: Paige & Saunders, 1975.

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

192: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPCG, KSPCR
193: M*/
196: PETSC_EXTERN PetscErrorCode KSPCreate_MINRES(KSP ksp)
197: {
198:   KSP_MINRES     *minres;

202:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
203:   PetscNewLog(ksp,KSP_MINRES,&minres);
204:   minres->haptol = 1.e-18;
205:   ksp->data      = (void*)minres;

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