Actual source code: symmlq.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_SYMMLQ;

  8: PetscErrorCode KSPSetUp_SYMMLQ(KSP ksp)
  9: {

 13:   KSPSetWorkVecs(ksp,9);
 14:   return(0);
 15: }

 17: PetscErrorCode  KSPSolve_SYMMLQ(KSP ksp)
 18: {
 20:   PetscInt       i;
 21:   PetscScalar    alpha,beta,ibeta,betaold,beta1,ceta = 0,ceta_oold = 0.0, ceta_old = 0.0,ceta_bar;
 22:   PetscScalar    c  = 1.0,cold=1.0,s=0.0,sold=0.0,coold,soold,rho0,rho1,rho2,rho3;
 23:   PetscScalar    dp = 0.0;
 24:   PetscReal      np,s_prod;
 25:   Vec            X,B,R,Z,U,V,W,UOLD,VOLD,Wbar;
 26:   Mat            Amat,Pmat;
 27:   KSP_SYMMLQ     *symmlq = (KSP_SYMMLQ*)ksp->data;
 28:   PetscBool      diagonalscale;

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

 34:   X    = ksp->vec_sol;
 35:   B    = ksp->vec_rhs;
 36:   R    = ksp->work[0];
 37:   Z    = ksp->work[1];
 38:   U    = ksp->work[2];
 39:   V    = ksp->work[3];
 40:   W    = ksp->work[4];
 41:   UOLD = ksp->work[5];
 42:   VOLD = ksp->work[6];
 43:   Wbar = ksp->work[7];

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

 47:   ksp->its = 0;

 49:   VecSet(UOLD,0.0);           /* u_old <- zeros;  */
 50:   VecCopy(UOLD,VOLD);          /* v_old <- u_old;  */
 51:   VecCopy(UOLD,W);             /* w     <- u_old;  */
 52:   VecCopy(UOLD,Wbar);          /* w_bar <- u_old;  */
 53:   if (!ksp->guess_zero) {
 54:     KSP_MatMult(ksp,Amat,X,R); /*     r <- b - A*x */
 55:     VecAYPX(R,-1.0,B);
 56:   } else {
 57:     VecCopy(B,R);              /*     r <- b (x is 0) */
 58:   }

 60:   KSP_PCApply(ksp,R,Z); /* z  <- B*r       */
 61:   VecDot(R,Z,&dp);             /* dp = r'*z;      */
 62:   KSPCheckDot(ksp,dp);
 63:   if (PetscAbsScalar(dp) < symmlq->haptol) {
 64:     PetscInfo2(ksp,"Detected happy breakdown %g tolerance %g\n",(double)PetscAbsScalar(dp),(double)symmlq->haptol);
 65:     ksp->rnorm  = 0.0;  /* what should we really put here? */
 66:     ksp->reason = KSP_CONVERGED_HAPPY_BREAKDOWN;  /* bugfix proposed by Lourens (lourens.vanzanen@shell.com) */
 67:     return(0);
 68:   }

 70: #if !defined(PETSC_USE_COMPLEX)
 71:   if (dp < 0.0) {
 72:     ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
 73:     return(0);
 74:   }
 75: #endif
 76:   dp     = PetscSqrtScalar(dp);
 77:   beta   = dp;                         /*  beta <- sqrt(r'*z)  */
 78:   beta1  = beta;
 79:   s_prod = PetscAbsScalar(beta1);

 81:   VecCopy(R,V); /* v <- r; */
 82:   VecCopy(Z,U); /* u <- z; */
 83:   ibeta = 1.0 / beta;
 84:   VecScale(V,ibeta);    /* v <- ibeta*v; */
 85:   VecScale(U,ibeta);    /* u <- ibeta*u; */
 86:   VecCopy(U,Wbar);       /* w_bar <- u;   */
 87:   VecNorm(Z,NORM_2,&np);     /*   np <- ||z||        */
 88:   KSPCheckNorm(ksp,np);
 89:   KSPLogResidualHistory(ksp,np);
 90:   KSPMonitor(ksp,0,np);
 91:   ksp->rnorm = np;
 92:   (*ksp->converged)(ksp,0,np,&ksp->reason,ksp->cnvP); /* test for convergence */
 93:   if (ksp->reason) return(0);

 95:   i = 0; ceta = 0.;
 96:   do {
 97:     ksp->its = i+1;

 99:     /*    Update    */
100:     if (ksp->its > 1) {
101:       VecCopy(V,VOLD);  /* v_old <- v; */
102:       VecCopy(U,UOLD);  /* u_old <- u; */

104:       VecCopy(R,V);
105:       VecScale(V,1.0/beta); /* v <- ibeta*r; */
106:       VecCopy(Z,U);
107:       VecScale(U,1.0/beta); /* u <- ibeta*z; */

109:       VecCopy(Wbar,W);
110:       VecScale(W,c);
111:       VecAXPY(W,s,U);   /* w  <- c*w_bar + s*u;    (w_k) */
112:       VecScale(Wbar,-s);
113:       VecAXPY(Wbar,c,U); /* w_bar <- -s*w_bar + c*u; (w_bar_(k+1)) */
114:       VecAXPY(X,ceta,W); /* x <- x + ceta * w;       (xL_k)  */

116:       ceta_oold = ceta_old;
117:       ceta_old  = ceta;
118:     }

120:     /*   Lanczos  */
121:     KSP_MatMult(ksp,Amat,U,R);   /*  r     <- Amat*u; */
122:     VecDot(U,R,&alpha);          /*  alpha <- u'*r;   */
123:     KSP_PCApply(ksp,R,Z); /*      z <- B*r;    */

125:     VecAXPY(R,-alpha,V);   /*  r <- r - alpha* v;  */
126:     VecAXPY(Z,-alpha,U);   /*  z <- z - alpha* u;  */
127:     VecAXPY(R,-beta,VOLD); /*  r <- r - beta * v_old; */
128:     VecAXPY(Z,-beta,UOLD); /*  z <- z - beta * u_old; */
129:     betaold = beta;                                /* beta_k                  */
130:     VecDot(R,Z,&dp);       /* dp <- r'*z;             */
131:     KSPCheckDot(ksp,dp);
132:     if (PetscAbsScalar(dp) < symmlq->haptol) {
133:       PetscInfo2(ksp,"Detected happy breakdown %g tolerance %g\n",(double)PetscAbsScalar(dp),(double)symmlq->haptol);
134:       dp   = 0.0;
135:     }

137: #if !defined(PETSC_USE_COMPLEX)
138:     if (dp < 0.0) {
139:       ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
140:       break;
141:     }
142: #endif
143:     beta = PetscSqrtScalar(dp);                    /*  beta = sqrt(dp); */

145:     /*    QR factorization    */
146:     coold = cold; cold = c; soold = sold; sold = s;
147:     rho0  = cold * alpha - coold * sold * betaold;   /* gamma_bar */
148:     rho1  = PetscSqrtScalar(rho0*rho0 + beta*beta);  /* gamma     */
149:     rho2  = sold * alpha + coold * cold * betaold;   /* delta     */
150:     rho3  = soold * betaold;                         /* epsilon   */

152:     /* Givens rotation: [c -s; s c] (different from the Reference!) */
153:     c = rho0 / rho1; s = beta / rho1;

155:     if (ksp->its==1) ceta = beta1/rho1;
156:     else ceta = -(rho2*ceta_old + rho3*ceta_oold)/rho1;

158:     s_prod = s_prod*PetscAbsScalar(s);
159:     if (c == 0.0) np = s_prod*1.e16;
160:     else np = s_prod/PetscAbsScalar(c);       /* residual norm for xc_k (CGNORM) */

162:     ksp->rnorm = np;
163:     KSPLogResidualHistory(ksp,np);
164:     KSPMonitor(ksp,i+1,np);
165:     (*ksp->converged)(ksp,i+1,np,&ksp->reason,ksp->cnvP); /* test for convergence */
166:     if (ksp->reason) break;
167:     i++;
168:   } while (i<ksp->max_it);

170:   /* move to the CG point: xc_(k+1) */
171:   if (c == 0.0) ceta_bar = ceta*1.e15;
172:   else ceta_bar = ceta/c;

174:   VecAXPY(X,ceta_bar,Wbar); /* x <- x + ceta_bar*w_bar */

176:   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
177:   return(0);
178: }

180: /*MC
181:      KSPSYMMLQ -  This code implements the SYMMLQ method.

183:    Options Database Keys:
184: .   see KSPSolve()

186:    Level: beginner

188:    Notes:
189:     The operator and the preconditioner must be symmetric for this method. The
190:           preconditioner must be POSITIVE-DEFINITE.

192:           Supports only left preconditioning.

194:    Reference: Paige & Saunders, 1975.

196: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP
197: M*/
198: PETSC_EXTERN PetscErrorCode KSPCreate_SYMMLQ(KSP ksp)
199: {
200:   KSP_SYMMLQ     *symmlq;

204:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);

206:   PetscNewLog(ksp,&symmlq);
207:   symmlq->haptol = 1.e-18;
208:   ksp->data      = (void*)symmlq;

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