Actual source code: qcg.c

 2:  #include src/ksp/ksp/kspimpl.h
 3:  #include src/ksp/ksp/impls/qcg/qcg.h

  5: static PetscErrorCode QuadraticRoots_Private(Vec,Vec,PetscReal*,PetscReal*,PetscReal*);

  9: /*@
 10:     KSPQCGSetTrustRegionRadius - Sets the radius of the trust region.

 12:     Collective on KSP

 14:     Input Parameters:
 15: +   ksp   - the iterative context
 16: -   delta - the trust region radius (Infinity is the default)

 18:     Options Database Key:
 19: .   -ksp_qcg_trustregionradius <delta>

 21:     Level: advanced

 23: .keywords: KSP, QCG, set, trust region radius
 24: @*/
 25: PetscErrorCode KSPQCGSetTrustRegionRadius(KSP ksp,PetscReal delta)
 26: {
 27:   PetscErrorCode ierr,(*f)(KSP,PetscReal);

 31:   if (delta < 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Tolerance must be non-negative");
 32:   PetscObjectQueryFunction((PetscObject)ksp,"KSPQCGSetTrustRegionRadius_C",(void (**)(void))&f);
 33:   if (f) {
 34:     (*f)(ksp,delta);
 35:   }

 37:   return(0);
 38: }

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

 46:     Collective on KSP

 48:     Input Parameter:
 49: .   ksp - the iterative context

 51:     Output Parameter:
 52: .   tsnorm - the norm

 54:     Level: advanced
 55: @*/
 56: PetscErrorCode KSPQCGGetTrialStepNorm(KSP ksp,PetscReal *tsnorm)
 57: {
 58:   PetscErrorCode ierr,(*f)(KSP,PetscReal*);

 62:   PetscObjectQueryFunction((PetscObject)ksp,"KSPQCGGetTrialStepNorm_C",(void (**)(void))&f);
 63:   if (f) {
 64:     (*f)(ksp,tsnorm);
 65:   }
 66:   return(0);
 67: }

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

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

 76:     which satisfies the Euclidian Norm trust region constraint

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

 80:     where

 82:      delta is the trust region radius, 
 83:      g is the gradient vector, and
 84:      H is Hessian matrix,
 85:      D is a scaling matrix.

 87:     Collective on KSP

 89:     Input Parameter:
 90: .   ksp - the iterative context

 92:     Output Parameter:
 93: .   quadratic - the quadratic function evaluated at the new iterate

 95:     Level: advanced
 96: @*/
 97: PetscErrorCode KSPQCGGetQuadratic(KSP ksp,PetscReal *quadratic)
 98: {
 99:   PetscErrorCode ierr,(*f)(KSP,PetscReal*);

103:   PetscObjectQueryFunction((PetscObject)ksp,"KSPQCGGetQuadratic_C",(void (**)(void))&f);
104:   if (f) {
105:     (*f)(ksp,quadratic);
106:   }
107:   return(0);
108: }

112: /* 
113:   KSPSolve_QCG - Use preconditioned conjugate gradient to compute 
114:   an approximate minimizer of the quadratic function 

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

118:    subject to the Euclidean norm trust region constraint

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

122:    where 

124:      delta is the trust region radius, 
125:      g is the gradient vector, and
126:      H is Hessian matrix,
127:      D is a scaling matrix.

129:    KSPConvergedReason may be 
130: $  KSP_CONVERGED_QCG_NEG_CURVE if convergence is reached along a negative curvature direction,
131: $  KSP_CONVERGED_QCG_CONSTRAINED if convergence is reached along a constrained step,
132: $  other KSP converged/diverged reasons

134:   This method is intended for use in conjunction with the TAO trust region method
135:   for unconstrained minimization (see www.mcs.anl.gov/tao).

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

144:  We should perhaps rewrite using PCApplyBAorAB().
145:  */
146: PetscErrorCode KSPSolve_QCG(KSP ksp)
147: {
148: /* 
149:    Correpondence with documentation above:  
150:       B = g = gradient,
151:       X = s = step
152:    Note:  This is not coded correctly for complex arithmetic!
153:  */

155:   KSP_QCG        *pcgP = (KSP_QCG*)ksp->data;
156:   MatStructure   pflag;
157:   Mat            Amat,Pmat;
158:   Vec            W,WA,WA2,R,P,ASP,BS,X,B;
159:   PetscScalar    zero = 0.0,negone = -1.0,scal,nstep,btx,xtax,beta,rntrn,step;
160:   PetscReal      ptasp,q1,q2,wtasp,bstp,rtr,xnorm,step1,step2,rnrm,p5 = 0.5;
161:   PetscReal      dzero = 0.0,bsnrm;
163:   PetscInt       i,maxit;
164:   PC             pc = ksp->pc;
165:   PCSide         side;
166: #if defined(PETSC_USE_COMPLEX)
167:   PetscScalar    cstep1,cstep2,cbstp,crtr,cwtasp,cptasp;
168: #endif
169:   PetscTruth     diagonalscale;

172:   PCDiagonalScale(ksp->pc,&diagonalscale);
173:   if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",ksp->type_name);
174:   if (ksp->transpose_solve) {
175:     SETERRQ(PETSC_ERR_SUP,"Currently does not support transpose solve");
176:   }

178:   ksp->its = 0;
179:   maxit    = ksp->max_it;
180:   WA       = ksp->work[0];
181:   R        = ksp->work[1];
182:   P        = ksp->work[2];
183:   ASP      = ksp->work[3];
184:   BS       = ksp->work[4];
185:   W        = ksp->work[5];
186:   WA2      = ksp->work[6];
187:   X        = ksp->vec_sol;
188:   B        = ksp->vec_rhs;

190:   if (pcgP->delta <= dzero) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Input error: delta <= 0");
191:   KSPGetPreconditionerSide(ksp,&side);
192:   if (side != PC_SYMMETRIC) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Requires symmetric preconditioner!");

194:   /* Initialize variables */
195:   VecSet(&zero,W);        /* W = 0 */
196:   VecSet(&zero,X);        /* X = 0 */
197:   PCGetOperators(pc,&Amat,&Pmat,&pflag);

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

202:   VecNorm(BS,NORM_2,&bsnrm);
203:   PetscObjectTakeAccess(ksp);
204:   ksp->its    = 0;
205:   ksp->rnorm  = bsnrm;
206:   PetscObjectGrantAccess(ksp);
207:   KSPLogResidualHistory(ksp,bsnrm);
208:   KSPMonitor(ksp,0,bsnrm);
209:   (*ksp->converged)(ksp,0,bsnrm,&ksp->reason,ksp->cnvP);
210:   if (ksp->reason) return(0);

212:   /* Compute the initial scaled direction and scaled residual */
213:   VecCopy(BS,R);
214:   VecScale(&negone,R);
215:   VecCopy(R,P);
216: #if defined(PETSC_USE_COMPLEX)
217:   VecDot(R,R,&crtr); rtr = PetscRealPart(crtr);
218: #else
219:   VecDot(R,R,&rtr);
220: #endif

222:   for (i=0; i<=maxit; i++) {
223:     PetscObjectTakeAccess(ksp);
224:     ksp->its++;
225:     PetscObjectGrantAccess(ksp);

227:     /* Compute:  asp = D^{-T}*A*D^{-1}*p  */
228:     PCApplySymmetricRight(pc,P,WA);
229:     MatMult(Amat,WA,WA2);
230:     PCApplySymmetricLeft(pc,WA2,ASP);

232:     /* Check for negative curvature */
233: #if defined(PETSC_USE_COMPLEX)
234:     VecDot(P,ASP,&cptasp);
235:     ptasp = PetscRealPart(cptasp);
236: #else
237:     VecDot(P,ASP,&ptasp);        /* ptasp = p^T asp */
238: #endif
239:     if (ptasp <= dzero) {

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

244:        if (!i) {
245:          VecCopy(P,X);
246:          VecNorm(X,NORM_2,&xnorm);
247:          scal = pcgP->delta / xnorm;
248:          VecScale(&scal,X);
249:        } else {
250:          /* Compute roots of quadratic */
251:          QuadraticRoots_Private(W,P,&pcgP->delta,&step1,&step2);
252: #if defined(PETSC_USE_COMPLEX)
253:          VecDot(W,ASP,&cwtasp); wtasp = PetscRealPart(cwtasp);
254:          VecDot(BS,P,&cbstp);   bstp  = PetscRealPart(cbstp);
255: #else
256:          VecDot(W,ASP,&wtasp);
257:          VecDot(BS,P,&bstp);
258: #endif
259:          VecCopy(W,X);
260:          q1 = step1*(bstp + wtasp + p5*step1*ptasp);
261:          q2 = step2*(bstp + wtasp + p5*step2*ptasp);
262: #if defined(PETSC_USE_COMPLEX)
263:          if (q1 <= q2) {
264:            cstep1 = step1; VecAXPY(&cstep1,P,X);
265:          } else {
266:            cstep2 = step2; VecAXPY(&cstep2,P,X);
267:          }
268: #else
269:          if (q1 <= q2) {VecAXPY(&step1,P,X);}
270:          else          {VecAXPY(&step2,P,X);}
271: #endif
272:        }
273:        pcgP->ltsnrm = pcgP->delta;                       /* convergence in direction of */
274:        ksp->reason  = KSP_CONVERGED_QCG_NEG_CURVE;  /* negative curvature */
275:        if (!i) {
276:          PetscLogInfo(ksp,"KSPSolve_QCG: negative curvature: delta=%g\n",pcgP->delta);
277:        } else {
278:          PetscLogInfo(ksp,"KSPSolve_QCG: negative curvature: step1=%g, step2=%g, delta=%g\n",step1,step2,pcgP->delta);
279:        }
280: 
281:     } else {
282: 
283:        /* Compute step along p */

285:        step = rtr/ptasp;
286:        VecCopy(W,X);           /*  x = w  */
287:        VecAXPY(&step,P,X);   /*  x <- step*p + x  */
288:        VecNorm(X,NORM_2,&pcgP->ltsnrm);

290:        if (pcgP->ltsnrm > pcgP->delta) {

292:          /* Since the trial iterate is outside the trust region, 
293:              evaluate a constrained step along p so that 
294:                       ||w + step*p|| = delta 
295:             The positive step is always better in this case. */

297:          if (!i) {
298:            scal = pcgP->delta / pcgP->ltsnrm;
299:            VecScale(&scal,X);
300:          } else {
301:            /* Compute roots of quadratic */
302:            QuadraticRoots_Private(W,P,&pcgP->delta,&step1,&step2);
303:            VecCopy(W,X);
304: #if defined(PETSC_USE_COMPLEX)
305:            cstep1 = step1; VecAXPY(&cstep1,P,X);
306: #else
307:            VecAXPY(&step1,P,X);  /*  x <- step1*p + x  */
308: #endif
309:          }
310:          pcgP->ltsnrm = pcgP->delta;
311:          ksp->reason  = KSP_CONVERGED_QCG_CONSTRAINED;        /* convergence along constrained step */
312:          if (!i) {
313:            PetscLogInfo(ksp,"KSPSolve_QCG: constrained step: delta=%g\n",pcgP->delta);
314:          } else {
315:            PetscLogInfo(ksp,"KSPSolve_QCG: constrained step: step1=%g, step2=%g, delta=%g\n",step1,step2,pcgP->delta);
316:          }

318:        } else {

320:          /* Evaluate the current step */

322:          VecCopy(X,W);        /* update interior iterate */
323:          nstep = -step;
324:          VecAXPY(&nstep,ASP,R); /* r <- -step*asp + r */
325:          VecNorm(R,NORM_2,&rnrm);

327:          PetscObjectTakeAccess(ksp);
328:          ksp->rnorm                                    = rnrm;
329:          PetscObjectGrantAccess(ksp);
330:          KSPLogResidualHistory(ksp,rnrm);
331:          KSPMonitor(ksp,i+1,rnrm);
332:          (*ksp->converged)(ksp,i+1,rnrm,&ksp->reason,ksp->cnvP);
333:          if (ksp->reason) {                 /* convergence for */
334: #if defined(PETSC_USE_COMPLEX)               
335:            PetscLogInfo(ksp,"KSPSolve_QCG: truncated step: step=%g, rnrm=%g, delta=%g\n",PetscRealPart(step),rnrm,pcgP->delta);
336: #else
337:            PetscLogInfo(ksp,"KSPSolve_QCG: truncated step: step=%g, rnrm=%g, delta=%g\n",step,rnrm,pcgP->delta);
338: #endif
339:          }
340:       }
341:     }
342:     if (ksp->reason) break;        /* Convergence has been attained */
343:     else {                /* Compute a new AS-orthogonal direction */
344:       VecDot(R,R,&rntrn);
345:       beta = rntrn/rtr;
346:       VecAYPX(&beta,R,P);        /*  p <- r + beta*p  */
347: #if defined(PETSC_USE_COMPLEX)
348:       rtr = PetscRealPart(rntrn);
349: #else
350:       rtr = rntrn;
351: #endif
352:     }
353:   }
354:   if (!ksp->reason) {
355:     ksp->reason = KSP_DIVERGED_ITS;
356:   }

358:   /* Unscale x */
359:   VecCopy(X,WA2);
360:   PCApplySymmetricRight(pc,WA2,X);

362:   MatMult(Amat,X,WA);
363:   VecDot(B,X,&btx);
364:   VecDot(X,WA,&xtax);
365: #if defined(PETSC_USE_COMPLEX)
366:   pcgP->quadratic = PetscRealPart(btx) + p5* PetscRealPart(xtax);
367: #else
368:   pcgP->quadratic = btx + p5*xtax;              /* Compute q(x) */
369: #endif
370:   return(0);
371: }

375: PetscErrorCode KSPSetUp_QCG(KSP ksp)
376: {

380:   /* Check user parameters and functions */
381:   if (ksp->pc_side == PC_RIGHT) {
382:     SETERRQ(2,"no right preconditioning for QCG");
383:   } else if (ksp->pc_side == PC_LEFT) {
384:     SETERRQ(2,"no left preconditioning for QCG");
385:   }

387:   /* Get work vectors from user code */
388:   KSPDefaultGetWork(ksp,7);
389:   return(0);
390: }

394: PetscErrorCode KSPDestroy_QCG(KSP ksp)
395: {
396:   KSP_QCG        *cgP = (KSP_QCG*)ksp->data;

400:   KSPDefaultFreeWork(ksp);
401: 
402:   /* Free the context variable */
403:   PetscFree(cgP);
404:   return(0);
405: }

410: PetscErrorCode KSPQCGSetTrustRegionRadius_QCG(KSP ksp,PetscReal delta)
411: {
412:   KSP_QCG *cgP = (KSP_QCG*)ksp->data;

415:   cgP->delta = delta;
416:   return(0);
417: }

423: PetscErrorCode KSPQCGGetTrialStepNorm_QCG(KSP ksp,PetscReal *ltsnrm)
424: {
425:   KSP_QCG *cgP = (KSP_QCG*)ksp->data;

428:   *ltsnrm = cgP->ltsnrm;
429:   return(0);
430: }

436: PetscErrorCode KSPQCGGetQuadratic_QCG(KSP ksp,PetscReal *quadratic)
437: {
438:   KSP_QCG *cgP = (KSP_QCG*)ksp->data;

441:   *quadratic = cgP->quadratic;
442:   return(0);
443: }

448: PetscErrorCode KSPSetFromOptions_QCG(KSP ksp)
449: {
451:   PetscReal      delta;
452:   KSP_QCG        *cgP = (KSP_QCG*)ksp->data;
453:   PetscTruth     flg;

456:   PetscOptionsHead("KSP QCG Options");
457:   PetscOptionsReal("-ksp_qcg_trustregionradius","Trust Region Radius","KSPQCGSetTrustRegionRadius",cgP->delta,&delta,&flg);
458:   if (flg) { KSPQCGSetTrustRegionRadius(ksp,delta); }
459:   PetscOptionsTail();
460:   return(0);
461: }

463: /*MC
464:      KSPQCG -   Code to run conjugate gradient method subject to a constraint
465:          on the solution norm. This is used in Trust Region methods for nonlinear equations, SNESTR

467:    Options Database Keys:
468: .      -ksp_qcg_trustregionradius <r> - Trust Region Radius

470:    Notes: This is rarely used directly

472:    Level: developer

474: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPQCGSetTrustRegionRadius()
475:            KSPQCGGetTrialStepNorm(), KSPQCGGetQuadratic()
476: M*/

481: PetscErrorCode KSPCreate_QCG(KSP ksp)
482: {
484:   KSP_QCG        *cgP;

487:   PetscMalloc(sizeof(KSP_QCG),&cgP);
488:   PetscMemzero(cgP,sizeof(KSP_QCG));
489:   PetscLogObjectMemory(ksp,sizeof(KSP_QCG));
490:   ksp->data                      = (void*)cgP;
491:   ksp->pc_side                   = PC_SYMMETRIC;
492:   ksp->ops->setup                = KSPSetUp_QCG;
493:   ksp->ops->setfromoptions       = KSPSetFromOptions_QCG;
494:   ksp->ops->solve                = KSPSolve_QCG;
495:   ksp->ops->destroy              = KSPDestroy_QCG;
496:   ksp->ops->buildsolution        = KSPDefaultBuildSolution;
497:   ksp->ops->buildresidual        = KSPDefaultBuildResidual;
498:   ksp->ops->setfromoptions       = 0;
499:   ksp->ops->view                 = 0;

501:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPQCGGetQuadratic_C",
502:                                     "KSPQCGGetQuadratic_QCG",
503:                                      KSPQCGGetQuadratic_QCG);
504:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPQCGGetTrialStepNorm_C",
505:                                     "KSPQCGGetTrialStepNorm_QCG",
506:                                      KSPQCGGetTrialStepNorm_QCG);
507:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPQCGSetTrustRegionRadius_C",
508:                                     "KSPQCGSetTrustRegionRadius_QCG",
509:                                      KSPQCGSetTrustRegionRadius_QCG);
510:   cgP->delta = PETSC_MAX; /* default trust region radius is infinite */
511:   return(0);
512: }

515: /* ---------------------------------------------------------- */
518: /* 
519:   QuadraticRoots_Private - Computes the roots of the quadratic,
520:          ||s + step*p|| - delta = 0 
521:    such that step1 >= 0 >= step2.
522:    where
523:       delta:
524:         On entry delta must contain scalar delta.
525:         On exit delta is unchanged.
526:       step1:
527:         On entry step1 need not be specified.
528:         On exit step1 contains the non-negative root.
529:       step2:
530:         On entry step2 need not be specified.
531:         On exit step2 contains the non-positive root.
532:    C code is translated from the Fortran version of the MINPACK-2 Project,
533:    Argonne National Laboratory, Brett M. Averick and Richard G. Carter.
534: */
535: static PetscErrorCode QuadraticRoots_Private(Vec s,Vec p,PetscReal *delta,PetscReal *step1,PetscReal *step2)
536: {
537:   PetscReal      zero = 0.0,dsq,ptp,pts,rad,sts;
539: #if defined(PETSC_USE_COMPLEX)
540:   PetscScalar    cptp,cpts,csts;
541: #endif

544: #if defined(PETSC_USE_COMPLEX)
545:   VecDot(p,s,&cpts); pts = PetscRealPart(cpts);
546:   VecDot(p,p,&cptp); ptp = PetscRealPart(cptp);
547:   VecDot(s,s,&csts); sts = PetscRealPart(csts);
548: #else
549:   VecDot(p,s,&pts);
550:   VecDot(p,p,&ptp);
551:   VecDot(s,s,&sts);
552: #endif
553:   dsq  = (*delta)*(*delta);
554:   rad  = sqrt((pts*pts) - ptp*(sts - dsq));
555:   if (pts > zero) {
556:     *step2 = -(pts + rad)/ptp;
557:     *step1 = (sts - dsq)/(ptp * *step2);
558:   } else {
559:     *step1 = -(pts - rad)/ptp;
560:     *step2 = (sts - dsq)/(ptp * *step1);
561:   }
562:   return(0);
563: }