Actual source code: minres.c

petsc-3.11.4 2019-09-28
Report Typos and Errors

  2:  #include <petsc/private/kspimpl.h>

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

  8: static PetscErrorCode KSPSetUp_MINRES(KSP ksp)
  9: {

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


 20: static PetscErrorCode  KSPSolve_MINRES(KSP ksp)
 21: {
 22:   PetscErrorCode    ierr;
 23:   PetscInt          i;
 24:   PetscScalar       alpha,beta,ibeta,betaold,eta,c=1.0,ceta,cold=1.0,coold,s=0.0,sold=0.0,soold;
 25:   PetscScalar       rho0,rho1,irho1,rho2,rho3,dp = 0.0;
 26:   const PetscScalar none = -1.0;
 27:   PetscReal         np;
 28:   Vec               X,B,R,Z,U,V,W,UOLD,VOLD,WOLD,WOOLD;
 29:   Mat               Amat,Pmat;
 30:   KSP_MINRES        *minres = (KSP_MINRES*)ksp->data;
 31:   PetscBool         diagonalscale;

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

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

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

 51:   ksp->its = 0;

 53:   VecSet(UOLD,0.0);          /*     u_old  <-   0   */
 54:   VecSet(VOLD,0.0);         /*     v_old  <-   0   */
 55:   VecSet(W,0.0);            /*     w      <-   0   */
 56:   VecSet(WOLD,0.0);         /*     w_old  <-   0   */

 58:   if (!ksp->guess_zero) {
 59:     KSP_MatMult(ksp,Amat,X,R); /*     r <- b - A*x    */
 60:     VecAYPX(R,-1.0,B);
 61:   } else {
 62:     VecCopy(B,R);              /*     r <- b (x is 0) */
 63:   }
 64:   KSP_PCApply(ksp,R,Z);       /*     z  <- B*r       */
 65:   VecNorm(Z,NORM_2,&np);      /*   np <- ||z||        */
 66:   KSPCheckNorm(ksp,np);
 67:   VecDot(R,Z,&dp);
 68:   KSPCheckDot(ksp,dp);

 70:   if (PetscRealPart(dp) < minres->haptol && np > minres->haptol) {
 71:     if (ksp->errorifnotconverged) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_CONV_FAILED,"Detected indefinite operator %g tolerance %g",(double)PetscRealPart(dp),(double)minres->haptol);
 72:     PetscInfo2(ksp,"Detected indefinite operator %g tolerance %g\n",(double)PetscRealPart(dp),(double)minres->haptol);
 73:     ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
 74:     return(0);
 75:   }

 77:   KSPLogResidualHistory(ksp,np);
 78:   KSPMonitor(ksp,0,np);
 79:   ksp->rnorm = np;
 80:   (*ksp->converged)(ksp,0,np,&ksp->reason,ksp->cnvP); /* test for convergence */
 81:   if (ksp->reason) return(0);

 83:   dp   = PetscAbsScalar(dp);
 84:   dp   = PetscSqrtScalar(dp);
 85:   beta = dp;                                        /*  beta <- sqrt(r'*z  */
 86:   eta  = beta;

 88:   VecCopy(R,V);
 89:   VecCopy(Z,U);
 90:   ibeta = 1.0 / beta;
 91:   VecScale(V,ibeta);        /*    v <- r / beta     */
 92:   VecScale(U,ibeta);        /*    u <- z / beta     */

 94:   i = 0;
 95:   do {
 96:     ksp->its = i+1;

 98:     /*   Lanczos  */

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

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

109:     betaold = beta;

111:     VecDot(R,Z,&dp);
112:     KSPCheckDot(ksp,dp);
113:     dp   = PetscAbsScalar(dp);
114:     beta = PetscSqrtScalar(dp);                               /*  beta <- sqrt(r'*z)   */

116:     /*    QR factorisation    */

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

120:     rho0 = cold * alpha - coold * sold * betaold;
121:     rho1 = PetscSqrtScalar(rho0*rho0 + beta*beta);
122:     rho2 = sold * alpha + coold * cold * betaold;
123:     rho3 = soold * betaold;

125:     /*     Givens rotation    */

127:     c = rho0 / rho1;
128:     s = beta / rho1;

130:     /*    Update    */

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

135:     VecCopy(U,W);           /*  w      <- u          */
136:     VecAXPY(W,-rho2,WOLD); /*  w <- w - rho2 w_old  */
137:     VecAXPY(W,-rho3,WOOLD); /*  w <- w - rho3 w_oold */
138:     irho1 = 1.0 / rho1;
139:     VecScale(W,irho1);     /*  w <- w / rho1        */

141:     ceta = c * eta;
142:     VecAXPY(X,ceta,W);      /*  x <- x + c eta w     */

144:     /*
145:         when dp is really small we have either convergence or an indefinite operator so compute true
146:         residual norm to check for convergence
147:     */
148:     if (PetscRealPart(dp) < minres->haptol) {
149:       PetscInfo2(ksp,"Possible indefinite operator %g tolerance %g\n",(double)PetscRealPart(dp),(double)minres->haptol);
150:       KSP_MatMult(ksp,Amat,X,VOLD);
151:       VecAXPY(VOLD,none,B);
152:       VecNorm(VOLD,NORM_2,&np);
153:       KSPCheckNorm(ksp,np);
154:     } else {
155:       /* otherwise compute new residual norm via recurrence relation */
156:       np = ksp->rnorm * PetscAbsScalar(s);
157:     }

159:     ksp->rnorm = np;
160:     KSPLogResidualHistory(ksp,np);
161:     KSPMonitor(ksp,i+1,np);
162:     (*ksp->converged)(ksp,i+1,np,&ksp->reason,ksp->cnvP); /* test for convergence */
163:     if (ksp->reason) break;

165:     if (PetscRealPart(dp) < minres->haptol) {
166:       if (ksp->errorifnotconverged) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_CONV_FAILED,"Detected indefinite operator %g tolerance %g",(double)PetscRealPart(dp),(double)minres->haptol);
167:       PetscInfo2(ksp,"Detected indefinite operator %g tolerance %g\n",(double)PetscRealPart(dp),(double)minres->haptol);
168:       ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
169:       break;
170:     }

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

181:     i++;
182:   } while (i<ksp->max_it);
183:   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
184:   return(0);
185: }

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

190:    Options Database Keys:
191: .   see KSPSolve()

193:    Level: beginner

195:    Notes:
196:     The operator and the preconditioner must be symmetric and the preconditioner must
197:           be positive definite for this method.
198:           Supports only left preconditioning.

200:    Reference: Paige & Saunders, 1975.

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

204: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPCG, KSPCR
205: M*/
206: PETSC_EXTERN PetscErrorCode KSPCreate_MINRES(KSP ksp)
207: {
208:   KSP_MINRES     *minres;

212:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
213:   PetscNewLog(ksp,&minres);

215:   /* this parameter is arbitrary; but e-50 didn't work for __float128 in one example */
216: #if defined(PETSC_USE_REAL___FLOAT128)
217:   minres->haptol = 1.e-100;
218: #elif defined(PETSC_USE_REAL_SINGLE)
219:   minres->haptol = 1.e-25;
220: #else
221:   minres->haptol = 1.e-50;
222: #endif
223:   ksp->data      = (void*)minres;

225:   /*
226:        Sets the functions that are associated with this data structure
227:        (in C++ this is the same as defining virtual functions)
228:   */
229:   ksp->ops->setup          = KSPSetUp_MINRES;
230:   ksp->ops->solve          = KSPSolve_MINRES;
231:   ksp->ops->destroy        = KSPDestroyDefault;
232:   ksp->ops->setfromoptions = 0;
233:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
234:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
235:   return(0);
236: }