Actual source code: qcg.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

  2:  #include <../src/ksp/ksp/impls/qcg/qcgimpl.h>

  4: static PetscErrorCode KSPQCGQuadraticRoots(Vec,Vec,PetscReal,PetscReal*,PetscReal*);

  6: /*@
  7:     KSPQCGSetTrustRegionRadius - Sets the radius of the trust region.

  9:     Logically Collective on ksp

 11:     Input Parameters:
 12: +   ksp   - the iterative context
 13: -   delta - the trust region radius (Infinity is the default)

 15:     Options Database Key:
 16: .   -ksp_qcg_trustregionradius <delta>

 18:     Level: advanced

 20: @*/
 21: PetscErrorCode  KSPQCGSetTrustRegionRadius(KSP ksp,PetscReal delta)
 22: {

 27:   if (delta < 0.0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Tolerance must be non-negative");
 28:   PetscTryMethod(ksp,"KSPQCGSetTrustRegionRadius_C",(KSP,PetscReal),(ksp,delta));
 29:   return(0);
 30: }

 32: /*@
 33:     KSPQCGGetTrialStepNorm - Gets the norm of a trial step vector.  The WCG step may be
 34:     constrained, so this is not necessarily the length of the ultimate step taken in QCG.

 36:     Not Collective

 38:     Input Parameter:
 39: .   ksp - the iterative context

 41:     Output Parameter:
 42: .   tsnorm - the norm

 44:     Level: advanced
 45: @*/
 46: PetscErrorCode  KSPQCGGetTrialStepNorm(KSP ksp,PetscReal *tsnorm)
 47: {

 52:   PetscUseMethod(ksp,"KSPQCGGetTrialStepNorm_C",(KSP,PetscReal*),(ksp,tsnorm));
 53:   return(0);
 54: }

 56: /*@
 57:     KSPQCGGetQuadratic - Gets the value of the quadratic function, evaluated at the new iterate:

 59:        q(s) = g^T * s + 0.5 * s^T * H * s

 61:     which satisfies the Euclidian Norm trust region constraint

 63:        || D * s || <= delta,

 65:     where

 67:      delta is the trust region radius,
 68:      g is the gradient vector, and
 69:      H is Hessian matrix,
 70:      D is a scaling matrix.

 72:     Collective on ksp

 74:     Input Parameter:
 75: .   ksp - the iterative context

 77:     Output Parameter:
 78: .   quadratic - the quadratic function evaluated at the new iterate

 80:     Level: advanced
 81: @*/
 82: PetscErrorCode  KSPQCGGetQuadratic(KSP ksp,PetscReal *quadratic)
 83: {

 88:   PetscUseMethod(ksp,"KSPQCGGetQuadratic_C",(KSP,PetscReal*),(ksp,quadratic));
 89:   return(0);
 90: }


 93: PetscErrorCode KSPSolve_QCG(KSP ksp)
 94: {
 95: /*
 96:    Correpondence with documentation above:
 97:       B = g = gradient,
 98:       X = s = step
 99:    Note:  This is not coded correctly for complex arithmetic!
100:  */

102:   KSP_QCG        *pcgP = (KSP_QCG*)ksp->data;
103:   Mat            Amat,Pmat;
104:   Vec            W,WA,WA2,R,P,ASP,BS,X,B;
105:   PetscScalar    scal,beta,rntrn,step;
106:   PetscReal      q1,q2,xnorm,step1,step2,rnrm,btx,xtax;
107:   PetscReal      ptasp,rtr,wtasp,bstp;
108:   PetscReal      dzero = 0.0,bsnrm;
110:   PetscInt       i,maxit;
111:   PC             pc = ksp->pc;
112:   PCSide         side;
113:   PetscBool      diagonalscale;

116:   PCGetDiagonalScale(ksp->pc,&diagonalscale);
117:   if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
118:   if (ksp->transpose_solve) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Currently does not support transpose solve");

120:   ksp->its = 0;
121:   maxit    = ksp->max_it;
122:   WA       = ksp->work[0];
123:   R        = ksp->work[1];
124:   P        = ksp->work[2];
125:   ASP      = ksp->work[3];
126:   BS       = ksp->work[4];
127:   W        = ksp->work[5];
128:   WA2      = ksp->work[6];
129:   X        = ksp->vec_sol;
130:   B        = ksp->vec_rhs;

132:   if (pcgP->delta <= dzero) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Input error: delta <= 0");
133:   KSPGetPCSide(ksp,&side);
134:   if (side != PC_SYMMETRIC) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Requires symmetric preconditioner!");

136:   /* Initialize variables */
137:   VecSet(W,0.0);  /* W = 0 */
138:   VecSet(X,0.0);  /* X = 0 */
139:   PCGetOperators(pc,&Amat,&Pmat);

141:   /* Compute:  BS = D^{-1} B */
142:   PCApplySymmetricLeft(pc,B,BS);

144:   VecNorm(BS,NORM_2,&bsnrm);
145:   KSPCheckNorm(ksp,bsnrm);
146:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
147:   ksp->its   = 0;
148:   ksp->rnorm = bsnrm;
149:   PetscObjectSAWsGrantAccess((PetscObject)ksp);
150:   KSPLogResidualHistory(ksp,bsnrm);
151:   KSPMonitor(ksp,0,bsnrm);
152:   (*ksp->converged)(ksp,0,bsnrm,&ksp->reason,ksp->cnvP);
153:   if (ksp->reason) return(0);

155:   /* Compute the initial scaled direction and scaled residual */
156:   VecCopy(BS,R);
157:   VecScale(R,-1.0);
158:   VecCopy(R,P);
159:   VecDotRealPart(R,R,&rtr);

161:   for (i=0; i<=maxit; i++) {
162:     PetscObjectSAWsTakeAccess((PetscObject)ksp);
163:     ksp->its++;
164:     PetscObjectSAWsGrantAccess((PetscObject)ksp);

166:     /* Compute:  asp = D^{-T}*A*D^{-1}*p  */
167:     PCApplySymmetricRight(pc,P,WA);
168:     KSP_MatMult(ksp,Amat,WA,WA2);
169:     PCApplySymmetricLeft(pc,WA2,ASP);

171:     /* Check for negative curvature */
172:     VecDotRealPart(P,ASP,&ptasp);
173:     if (ptasp <= dzero) {

175:       /* Scaled negative curvature direction:  Compute a step so that
176:         ||w + step*p|| = delta and QS(w + step*p) is least */

178:       if (!i) {
179:         VecCopy(P,X);
180:         VecNorm(X,NORM_2,&xnorm);
181:         KSPCheckNorm(ksp,xnorm);
182:         scal = pcgP->delta / xnorm;
183:         VecScale(X,scal);
184:       } else {
185:         /* Compute roots of quadratic */
186:         KSPQCGQuadraticRoots(W,P,pcgP->delta,&step1,&step2);
187:         VecDotRealPart(W,ASP,&wtasp);
188:         VecDotRealPart(BS,P,&bstp);
189:         VecCopy(W,X);
190:         q1   = step1*(bstp + wtasp + .5*step1*ptasp);
191:         q2   = step2*(bstp + wtasp + .5*step2*ptasp);
192:         if (q1 <= q2) {
193:           VecAXPY(X,step1,P);
194:         } else {
195:           VecAXPY(X,step2,P);
196:         }
197:       }
198:       pcgP->ltsnrm = pcgP->delta;                       /* convergence in direction of */
199:       ksp->reason  = KSP_CONVERGED_CG_NEG_CURVE;  /* negative curvature */
200:       if (!i) {
201:         PetscInfo1(ksp,"negative curvature: delta=%g\n",(double)pcgP->delta);
202:       } else {
203:         PetscInfo3(ksp,"negative curvature: step1=%g, step2=%g, delta=%g\n",(double)step1,(double)step2,(double)pcgP->delta);
204:       }

206:     } else {
207:       /* Compute step along p */
208:       step = rtr/ptasp;
209:       VecCopy(W,X);        /*  x = w  */
210:       VecAXPY(X,step,P);   /*  x <- step*p + x  */
211:       VecNorm(X,NORM_2,&pcgP->ltsnrm);
212:       KSPCheckNorm(ksp,pcgP->ltsnrm);

214:       if (pcgP->ltsnrm > pcgP->delta) {
215:         /* Since the trial iterate is outside the trust region,
216:             evaluate a constrained step along p so that
217:                     ||w + step*p|| = delta
218:           The positive step is always better in this case. */
219:         if (!i) {
220:           scal = pcgP->delta / pcgP->ltsnrm;
221:           VecScale(X,scal);
222:         } else {
223:           /* Compute roots of quadratic */
224:           KSPQCGQuadraticRoots(W,P,pcgP->delta,&step1,&step2);
225:           VecCopy(W,X);
226:           VecAXPY(X,step1,P);  /*  x <- step1*p + x  */
227:         }
228:         pcgP->ltsnrm = pcgP->delta;
229:         ksp->reason  = KSP_CONVERGED_CG_CONSTRAINED; /* convergence along constrained step */
230:         if (!i) {
231:           PetscInfo1(ksp,"constrained step: delta=%g\n",(double)pcgP->delta);
232:         } else {
233:           PetscInfo3(ksp,"constrained step: step1=%g, step2=%g, delta=%g\n",(double)step1,(double)step2,(double)pcgP->delta);
234:         }

236:       } else {
237:         /* Evaluate the current step */
238:         VecCopy(X,W);  /* update interior iterate */
239:         VecAXPY(R,-step,ASP); /* r <- -step*asp + r */
240:         VecNorm(R,NORM_2,&rnrm);
241:         KSPCheckNorm(ksp,rnrm);
242:         PetscObjectSAWsTakeAccess((PetscObject)ksp);
243:         ksp->rnorm = rnrm;
244:         PetscObjectSAWsGrantAccess((PetscObject)ksp);
245:         KSPLogResidualHistory(ksp,rnrm);
246:         KSPMonitor(ksp,i+1,rnrm);
247:         (*ksp->converged)(ksp,i+1,rnrm,&ksp->reason,ksp->cnvP);
248:         if (ksp->reason) {                 /* convergence for */
249:           PetscInfo3(ksp,"truncated step: step=%g, rnrm=%g, delta=%g\n",(double)PetscRealPart(step),(double)rnrm,(double)pcgP->delta);
250:         }
251:       }
252:     }
253:     if (ksp->reason) break;  /* Convergence has been attained */
254:     else {                   /* Compute a new AS-orthogonal direction */
255:       VecDot(R,R,&rntrn);
256:       beta = rntrn/rtr;
257:       VecAYPX(P,beta,R);  /*  p <- r + beta*p  */
258:       rtr  = PetscRealPart(rntrn);
259:     }
260:   }
261:   if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;

263:   /* Unscale x */
264:   VecCopy(X,WA2);
265:   PCApplySymmetricRight(pc,WA2,X);

267:   KSP_MatMult(ksp,Amat,X,WA);
268:   VecDotRealPart(B,X,&btx);
269:   VecDotRealPart(X,WA,&xtax);

271:   pcgP->quadratic = btx + .5*xtax;
272:   return(0);
273: }

275: PetscErrorCode KSPSetUp_QCG(KSP ksp)
276: {

280:   /* Get work vectors from user code */
281:   KSPSetWorkVecs(ksp,7);
282:   return(0);
283: }

285: PetscErrorCode KSPDestroy_QCG(KSP ksp)
286: {

290:   PetscObjectComposeFunction((PetscObject)ksp,"KSPQCGGetQuadratic_C",NULL);
291:   PetscObjectComposeFunction((PetscObject)ksp,"KSPQCGGetTrialStepNorm_C",NULL);
292:   PetscObjectComposeFunction((PetscObject)ksp,"KSPQCGSetTrustRegionRadius_C",NULL);
293:   KSPDestroyDefault(ksp);
294:   return(0);
295: }

297: static PetscErrorCode  KSPQCGSetTrustRegionRadius_QCG(KSP ksp,PetscReal delta)
298: {
299:   KSP_QCG *cgP = (KSP_QCG*)ksp->data;

302:   cgP->delta = delta;
303:   return(0);
304: }

306: static PetscErrorCode  KSPQCGGetTrialStepNorm_QCG(KSP ksp,PetscReal *ltsnrm)
307: {
308:   KSP_QCG *cgP = (KSP_QCG*)ksp->data;

311:   *ltsnrm = cgP->ltsnrm;
312:   return(0);
313: }

315: static PetscErrorCode  KSPQCGGetQuadratic_QCG(KSP ksp,PetscReal *quadratic)
316: {
317:   KSP_QCG *cgP = (KSP_QCG*)ksp->data;

320:   *quadratic = cgP->quadratic;
321:   return(0);
322: }

324: PetscErrorCode KSPSetFromOptions_QCG(PetscOptionItems *PetscOptionsObject,KSP ksp)
325: {
327:   PetscReal      delta;
328:   KSP_QCG        *cgP = (KSP_QCG*)ksp->data;
329:   PetscBool      flg;

332:   PetscOptionsHead(PetscOptionsObject,"KSP QCG Options");
333:   PetscOptionsReal("-ksp_qcg_trustregionradius","Trust Region Radius","KSPQCGSetTrustRegionRadius",cgP->delta,&delta,&flg);
334:   if (flg) { KSPQCGSetTrustRegionRadius(ksp,delta); }
335:   PetscOptionsTail();
336:   return(0);
337: }

339: /*MC
340:      KSPQCG -   Code to run conjugate gradient method subject to a constraint
341:          on the solution norm. This is used in Trust Region methods for nonlinear equations, SNESNEWTONTR

343:    Options Database Keys:
344: .      -ksp_qcg_trustregionradius <r> - Trust Region Radius

346:    Notes:
347:     This is rarely used directly

349:    Level: developer

351:   Notes:
352:     Use preconditioned conjugate gradient to compute
353:       an approximate minimizer of the quadratic function

355:             q(s) = g^T * s + .5 * s^T * H * s

357:    subject to the Euclidean norm trust region constraint

359:             || D * s || <= delta,

361:    where

363:      delta is the trust region radius,
364:      g is the gradient vector, and
365:      H is Hessian matrix,
366:      D is a scaling matrix.

368:    KSPConvergedReason may be
369: $  KSP_CONVERGED_CG_NEG_CURVE if convergence is reached along a negative curvature direction,
370: $  KSP_CONVERGED_CG_CONSTRAINED if convergence is reached along a constrained step,
371: $  other KSP converged/diverged reasons


374:   Notes:
375:   Currently we allow symmetric preconditioning with the following scaling matrices:
376:       PCNONE:   D = Identity matrix
377:       PCJACOBI: D = diag [d_1, d_2, ...., d_n], where d_i = sqrt(H[i,i])
378:       PCICC:    D = L^T, implemented with forward and backward solves.
379:                 Here L is an incomplete Cholesky factor of H.

381:   References:
382: .  1. - Trond Steihaug, The Conjugate Gradient Method and Trust Regions in Large Scale Optimization,
383:    SIAM Journal on Numerical Analysis, Vol. 20, No. 3 (Jun., 1983).

385: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPQCGSetTrustRegionRadius()
386:            KSPQCGGetTrialStepNorm(), KSPQCGGetQuadratic()
387: M*/

389: PETSC_EXTERN PetscErrorCode KSPCreate_QCG(KSP ksp)
390: {
392:   KSP_QCG        *cgP;

395:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_SYMMETRIC,3);
396:   PetscNewLog(ksp,&cgP);

398:   ksp->data                = (void*)cgP;
399:   ksp->ops->setup          = KSPSetUp_QCG;
400:   ksp->ops->setfromoptions = KSPSetFromOptions_QCG;
401:   ksp->ops->solve          = KSPSolve_QCG;
402:   ksp->ops->destroy        = KSPDestroy_QCG;
403:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
404:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
405:   ksp->ops->view           = 0;

407:   PetscObjectComposeFunction((PetscObject)ksp,"KSPQCGGetQuadratic_C",KSPQCGGetQuadratic_QCG);
408:   PetscObjectComposeFunction((PetscObject)ksp,"KSPQCGGetTrialStepNorm_C",KSPQCGGetTrialStepNorm_QCG);
409:   PetscObjectComposeFunction((PetscObject)ksp,"KSPQCGSetTrustRegionRadius_C",KSPQCGSetTrustRegionRadius_QCG);
410:   cgP->delta = PETSC_MAX_REAL; /* default trust region radius is infinite */
411:   return(0);
412: }

414: /* ---------------------------------------------------------- */
415: /*
416:   KSPQCGQuadraticRoots - Computes the roots of the quadratic,
417:          ||s + step*p|| - delta = 0
418:    such that step1 >= 0 >= step2.
419:    where
420:       delta:
421:         On entry delta must contain scalar delta.
422:         On exit delta is unchanged.
423:       step1:
424:         On entry step1 need not be specified.
425:         On exit step1 contains the non-negative root.
426:       step2:
427:         On entry step2 need not be specified.
428:         On exit step2 contains the non-positive root.
429:    C code is translated from the Fortran version of the MINPACK-2 Project,
430:    Argonne National Laboratory, Brett M. Averick and Richard G. Carter.
431: */
432: static PetscErrorCode KSPQCGQuadraticRoots(Vec s,Vec p,PetscReal delta,PetscReal *step1,PetscReal *step2)
433: {
434:   PetscReal      dsq,ptp,pts,rad,sts;

438:   VecDotRealPart(p,s,&pts);
439:   VecDotRealPart(p,p,&ptp);
440:   VecDotRealPart(s,s,&sts);
441:   dsq  = delta*delta;
442:   rad  = PetscSqrtReal((pts*pts) - ptp*(sts - dsq));
443:   if (pts > 0.0) {
444:     *step2 = -(pts + rad)/ptp;
445:     *step1 = (sts - dsq)/(ptp * *step2);
446:   } else {
447:     *step1 = -(pts - rad)/ptp;
448:     *step2 = (sts - dsq)/(ptp * *step1);
449:   }
450:   return(0);
451: }