Actual source code: qcg.c

  1: #define PETSCKSP_DLL

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

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

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

 13:     Collective on KSP

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

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

 22:     Level: advanced

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

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

 38:   return(0);
 39: }

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

 47:     Collective on KSP

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

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

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

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

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

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

 77:     which satisfies the Euclidian Norm trust region constraint

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

 81:     where

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

 88:     Collective on KSP

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

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

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

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


114: PetscErrorCode KSPSolve_QCG(KSP ksp)
115: {
116: /* 
117:    Correpondence with documentation above:  
118:       B = g = gradient,
119:       X = s = step
120:    Note:  This is not coded correctly for complex arithmetic!
121:  */

123:   KSP_QCG        *pcgP = (KSP_QCG*)ksp->data;
124:   MatStructure   pflag;
125:   Mat            Amat,Pmat;
126:   Vec            W,WA,WA2,R,P,ASP,BS,X,B;
127:   PetscScalar    scal,btx,xtax,beta,rntrn,step;
128:   PetscReal      ptasp,q1,q2,wtasp,bstp,rtr,xnorm,step1,step2,rnrm,p5 = 0.5;
129:   PetscReal      dzero = 0.0,bsnrm;
131:   PetscInt       i,maxit;
132:   PC             pc = ksp->pc;
133:   PCSide         side;
134: #if defined(PETSC_USE_COMPLEX)
135:   PetscScalar    cstep1,cstep2,cbstp,crtr,cwtasp,cptasp;
136: #endif
137:   PetscTruth     diagonalscale;

140:   PCDiagonalScale(ksp->pc,&diagonalscale);
141:   if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
142:   if (ksp->transpose_solve) {
143:     SETERRQ(PETSC_ERR_SUP,"Currently does not support transpose solve");
144:   }

146:   ksp->its = 0;
147:   maxit    = ksp->max_it;
148:   WA       = ksp->work[0];
149:   R        = ksp->work[1];
150:   P        = ksp->work[2];
151:   ASP      = ksp->work[3];
152:   BS       = ksp->work[4];
153:   W        = ksp->work[5];
154:   WA2      = ksp->work[6];
155:   X        = ksp->vec_sol;
156:   B        = ksp->vec_rhs;

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

162:   /* Initialize variables */
163:   VecSet(W,0.0);        /* W = 0 */
164:   VecSet(X,0.0);        /* X = 0 */
165:   PCGetOperators(pc,&Amat,&Pmat,&pflag);

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

170:   VecNorm(BS,NORM_2,&bsnrm);
171:   PetscObjectTakeAccess(ksp);
172:   ksp->its    = 0;
173:   ksp->rnorm  = bsnrm;
174:   PetscObjectGrantAccess(ksp);
175:   KSPLogResidualHistory(ksp,bsnrm);
176:   KSPMonitor(ksp,0,bsnrm);
177:   (*ksp->converged)(ksp,0,bsnrm,&ksp->reason,ksp->cnvP);
178:   if (ksp->reason) return(0);

180:   /* Compute the initial scaled direction and scaled residual */
181:   VecCopy(BS,R);
182:   VecScale(R,-1.0);
183:   VecCopy(R,P);
184: #if defined(PETSC_USE_COMPLEX)
185:   VecDot(R,R,&crtr); rtr = PetscRealPart(crtr);
186: #else
187:   VecDot(R,R,&rtr);
188: #endif

190:   for (i=0; i<=maxit; i++) {
191:     PetscObjectTakeAccess(ksp);
192:     ksp->its++;
193:     PetscObjectGrantAccess(ksp);

195:     /* Compute:  asp = D^{-T}*A*D^{-1}*p  */
196:     PCApplySymmetricRight(pc,P,WA);
197:     MatMult(Amat,WA,WA2);
198:     PCApplySymmetricLeft(pc,WA2,ASP);

200:     /* Check for negative curvature */
201: #if defined(PETSC_USE_COMPLEX)
202:     VecDot(P,ASP,&cptasp);
203:     ptasp = PetscRealPart(cptasp);
204: #else
205:     VecDot(P,ASP,&ptasp);        /* ptasp = p^T asp */
206: #endif
207:     if (ptasp <= dzero) {

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

212:        if (!i) {
213:          VecCopy(P,X);
214:          VecNorm(X,NORM_2,&xnorm);
215:          scal = pcgP->delta / xnorm;
216:          VecScale(X,scal);
217:        } else {
218:          /* Compute roots of quadratic */
219:          QuadraticRoots_Private(W,P,&pcgP->delta,&step1,&step2);
220: #if defined(PETSC_USE_COMPLEX)
221:          VecDot(W,ASP,&cwtasp); wtasp = PetscRealPart(cwtasp);
222:          VecDot(BS,P,&cbstp);   bstp  = PetscRealPart(cbstp);
223: #else
224:          VecDot(W,ASP,&wtasp);
225:          VecDot(BS,P,&bstp);
226: #endif
227:          VecCopy(W,X);
228:          q1 = step1*(bstp + wtasp + p5*step1*ptasp);
229:          q2 = step2*(bstp + wtasp + p5*step2*ptasp);
230: #if defined(PETSC_USE_COMPLEX)
231:          if (q1 <= q2) {
232:            cstep1 = step1; VecAXPY(X,cstep1,P);
233:          } else {
234:            cstep2 = step2; VecAXPY(X,cstep2,P);
235:          }
236: #else
237:          if (q1 <= q2) {VecAXPY(X,step1,P);}
238:          else          {VecAXPY(X,step2,P);}
239: #endif
240:        }
241:        pcgP->ltsnrm = pcgP->delta;                       /* convergence in direction of */
242:        ksp->reason  = KSP_CONVERGED_CG_NEG_CURVE;  /* negative curvature */
243:        if (!i) {
244:          PetscInfo1(ksp,"negative curvature: delta=%G\n",pcgP->delta);
245:        } else {
246:          PetscInfo3(ksp,"negative curvature: step1=%G, step2=%G, delta=%G\n",step1,step2,pcgP->delta);
247:        }
248: 
249:     } else {
250: 
251:        /* Compute step along p */

253:        step = rtr/ptasp;
254:        VecCopy(W,X);           /*  x = w  */
255:        VecAXPY(X,step,P);   /*  x <- step*p + x  */
256:        VecNorm(X,NORM_2,&pcgP->ltsnrm);

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

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

265:          if (!i) {
266:            scal = pcgP->delta / pcgP->ltsnrm;
267:            VecScale(X,scal);
268:          } else {
269:            /* Compute roots of quadratic */
270:            QuadraticRoots_Private(W,P,&pcgP->delta,&step1,&step2);
271:            VecCopy(W,X);
272: #if defined(PETSC_USE_COMPLEX)
273:            cstep1 = step1; VecAXPY(X,cstep1,P);
274: #else
275:            VecAXPY(X,step1,P);  /*  x <- step1*p + x  */
276: #endif
277:          }
278:          pcgP->ltsnrm = pcgP->delta;
279:          ksp->reason  = KSP_CONVERGED_CG_CONSTRAINED;        /* convergence along constrained step */
280:          if (!i) {
281:            PetscInfo1(ksp,"constrained step: delta=%G\n",pcgP->delta);
282:          } else {
283:            PetscInfo3(ksp,"constrained step: step1=%G, step2=%G, delta=%G\n",step1,step2,pcgP->delta);
284:          }

286:        } else {

288:          /* Evaluate the current step */

290:          VecCopy(X,W);        /* update interior iterate */
291:          VecAXPY(R,-step,ASP); /* r <- -step*asp + r */
292:          VecNorm(R,NORM_2,&rnrm);

294:          PetscObjectTakeAccess(ksp);
295:          ksp->rnorm                                    = rnrm;
296:          PetscObjectGrantAccess(ksp);
297:          KSPLogResidualHistory(ksp,rnrm);
298:          KSPMonitor(ksp,i+1,rnrm);
299:          (*ksp->converged)(ksp,i+1,rnrm,&ksp->reason,ksp->cnvP);
300:          if (ksp->reason) {                 /* convergence for */
301: #if defined(PETSC_USE_COMPLEX)               
302:            PetscInfo3(ksp,"truncated step: step=%G, rnrm=%G, delta=%G\n",PetscRealPart(step),rnrm,pcgP->delta);
303: #else
304:            PetscInfo3(ksp,"truncated step: step=%G, rnrm=%G, delta=%G\n",step,rnrm,pcgP->delta);
305: #endif
306:          }
307:       }
308:     }
309:     if (ksp->reason) break;        /* Convergence has been attained */
310:     else {                /* Compute a new AS-orthogonal direction */
311:       VecDot(R,R,&rntrn);
312:       beta = rntrn/rtr;
313:       VecAYPX(P,beta,R);        /*  p <- r + beta*p  */
314: #if defined(PETSC_USE_COMPLEX)
315:       rtr = PetscRealPart(rntrn);
316: #else
317:       rtr = rntrn;
318: #endif
319:     }
320:   }
321:   if (!ksp->reason) {
322:     ksp->reason = KSP_DIVERGED_ITS;
323:   }

325:   /* Unscale x */
326:   VecCopy(X,WA2);
327:   PCApplySymmetricRight(pc,WA2,X);

329:   MatMult(Amat,X,WA);
330:   VecDot(B,X,&btx);
331:   VecDot(X,WA,&xtax);
332: #if defined(PETSC_USE_COMPLEX)
333:   pcgP->quadratic = PetscRealPart(btx) + p5* PetscRealPart(xtax);
334: #else
335:   pcgP->quadratic = btx + p5*xtax;              /* Compute q(x) */
336: #endif
337:   return(0);
338: }

342: PetscErrorCode KSPSetUp_QCG(KSP ksp)
343: {

347:   /* Check user parameters and functions */
348:   if (ksp->pc_side == PC_RIGHT) {
349:     SETERRQ(PETSC_ERR_SUP,"no right preconditioning for QCG");
350:   } else if (ksp->pc_side == PC_LEFT) {
351:     SETERRQ(PETSC_ERR_SUP,"no left preconditioning for QCG");
352:   }

354:   /* Get work vectors from user code */
355:   KSPDefaultGetWork(ksp,7);
356:   return(0);
357: }

361: PetscErrorCode KSPDestroy_QCG(KSP ksp)
362: {

366:   KSPDefaultDestroy(ksp);
367:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPQCGGetQuadratic_C","",PETSC_NULL);
368:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPQCGGetTrialStepNorm_C","",PETSC_NULL);
369:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPQCGSetTrustRegionRadius_C","",PETSC_NULL);
370:   return(0);
371: }

376: PetscErrorCode  KSPQCGSetTrustRegionRadius_QCG(KSP ksp,PetscReal delta)
377: {
378:   KSP_QCG *cgP = (KSP_QCG*)ksp->data;

381:   cgP->delta = delta;
382:   return(0);
383: }

389: PetscErrorCode  KSPQCGGetTrialStepNorm_QCG(KSP ksp,PetscReal *ltsnrm)
390: {
391:   KSP_QCG *cgP = (KSP_QCG*)ksp->data;

394:   *ltsnrm = cgP->ltsnrm;
395:   return(0);
396: }

402: PetscErrorCode  KSPQCGGetQuadratic_QCG(KSP ksp,PetscReal *quadratic)
403: {
404:   KSP_QCG *cgP = (KSP_QCG*)ksp->data;

407:   *quadratic = cgP->quadratic;
408:   return(0);
409: }

414: PetscErrorCode KSPSetFromOptions_QCG(KSP ksp)
415: {
417:   PetscReal      delta;
418:   KSP_QCG        *cgP = (KSP_QCG*)ksp->data;
419:   PetscTruth     flg;

422:   PetscOptionsHead("KSP QCG Options");
423:   PetscOptionsReal("-ksp_qcg_trustregionradius","Trust Region Radius","KSPQCGSetTrustRegionRadius",cgP->delta,&delta,&flg);
424:   if (flg) { KSPQCGSetTrustRegionRadius(ksp,delta); }
425:   PetscOptionsTail();
426:   return(0);
427: }

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

433:    Options Database Keys:
434: .      -ksp_qcg_trustregionradius <r> - Trust Region Radius

436:    Notes: This is rarely used directly

438:    Level: developer

440:   Notes:  Use preconditioned conjugate gradient to compute 
441:       an approximate minimizer of the quadratic function 

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

445:    subject to the Euclidean norm trust region constraint

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

449:    where 

451:      delta is the trust region radius, 
452:      g is the gradient vector, and
453:      H is Hessian matrix,
454:      D is a scaling matrix.

456:    KSPConvergedReason may be 
457: $  KSP_CONVERGED_CG_NEG_CURVE if convergence is reached along a negative curvature direction,
458: $  KSP_CONVERGED_CG_CONSTRAINED if convergence is reached along a constrained step,
459: $  other KSP converged/diverged reasons


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

469:   References:
470:    The Conjugate Gradient Method and Trust Regions in Large Scale Optimization, Trond Steihaug
471:    SIAM Journal on Numerical Analysis, Vol. 20, No. 3 (Jun., 1983), pp. 626-637

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

480: PetscErrorCode  KSPCreate_QCG(KSP ksp)
481: {
483:   KSP_QCG        *cgP;

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

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

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

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