Actual source code: gltr.c

  1: #define PETSCKSP_DLL

  3: #include <math.h>
 4:  #include include/private/kspimpl.h
 5:  #include include/petscblaslapack.h
 6:  #include src/ksp/ksp/impls/cg/gltr/gltr.h

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

 13:     Collective on KSP

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

 19:     Options Database Key:
 20: .   -ksp_gltr_radius <r>

 22:     Level: advanced

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

 32:   if (radius < 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Tolerance must be positive");
 33:   PetscObjectQueryFunction((PetscObject)ksp, "KSPGLTRSetRadius_C", (void (**)(void))&f);
 34:   if (f) {
 35:     (*f)(ksp, radius);
 36:   }
 37:   return(0);
 38: }

 42: /*@
 43:     KSPGLTRGetNormD - Get norm of the direction.

 45:     Collective on KSP

 47:     Input Parameters:
 48: +   ksp    - the iterative context
 49: -   norm_d - the norm of the direction

 51:     Level: advanced

 53: .keywords: KSP, GLTR, get, norm direction
 54: @*/
 55: PetscErrorCode  KSPGLTRGetNormD(KSP ksp,PetscReal *norm_d)
 56: {
 57:   PetscErrorCode ierr, (*f)(KSP, PetscReal *);

 61:   PetscObjectQueryFunction((PetscObject)ksp, "KSPGLTRGetNormD_C", (void (**)(void))&f);
 62:   if (f) {
 63:     (*f)(ksp, norm_d);
 64:   }
 65:   return(0);
 66: }

 70: /*@
 71:     KSPGLTRGetObjFcn - Get objective function value.

 73:     Collective on KSP

 75:     Input Parameters:
 76: +   ksp   - the iterative context
 77: -   o_fcn - the objective function value

 79:     Level: advanced

 81: .keywords: KSP, GLTR, get, objective function
 82: @*/
 83: PetscErrorCode  KSPGLTRGetObjFcn(KSP ksp,PetscReal *o_fcn)
 84: {
 85:   PetscErrorCode ierr, (*f)(KSP, PetscReal *);

 89:   PetscObjectQueryFunction((PetscObject)ksp, "KSPGLTRGetObjFcn_C", (void (**)(void))&f);
 90:   if (f) {
 91:     (*f)(ksp, o_fcn);
 92:   }
 93:   return(0);
 94: }

 98: /*@
 99:     KSPGLTRGetMinEig - Get minimum eigenvalue.

101:     Collective on KSP

103:     Input Parameters:
104: +   ksp   - the iterative context
105: -   e_min - the minimum eigenvalue

107:     Level: advanced

109: .keywords: KSP, GLTR, get, minimum eigenvalue
110: @*/
111: PetscErrorCode  KSPGLTRGetMinEig(KSP ksp,PetscReal *e_min)
112: {
113:   PetscErrorCode ierr, (*f)(KSP, PetscReal *);

117:   PetscObjectQueryFunction((PetscObject)ksp, "KSPGLTRGetMinEig_C", (void (**)(void))&f);
118:   if (f) {
119:     (*f)(ksp, e_min);
120:   }
121:   return(0);
122: }

126: /*@
127:     KSPGLTRGetLambda - Get multiplier on trust-region constraint.

129:     Collective on KSP

131:     Input Parameters:
132: +   ksp    - the iterative context
133: -   lambda - the multiplier

135:     Level: advanced

137: .keywords: KSP, GLTR, get, multiplier
138: @*/
139: PetscErrorCode  KSPGLTRGetLambda(KSP ksp,PetscReal *lambda)
140: {
141:   PetscErrorCode ierr, (*f)(KSP, PetscReal *);

145:   PetscObjectQueryFunction((PetscObject)ksp, "KSPGLTRGetLambda_C", (void (**)(void))&f);
146:   if (f) {
147:     (*f)(ksp, lambda);
148:   }
149:   return(0);
150: }

154: /*
155:   KSPSolve_GLTR - Use preconditioned conjugate gradient to compute
156:   an approximate minimizer of the quadratic function

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

160:    subject to the trust region constraint

162:             || s ||_M <= delta,

164:    where

166:      delta is the trust region radius,
167:      g is the gradient vector, and
168:      H is the Hessian matrix,
169:      M is the positive definite preconditioner matrix.

171:    KSPConvergedReason may be
172: $  KSP_CONVERGED_CG_NEG_CURVE if convergence is reached along a negative curvature direction,
173: $  KSP_CONVERGED_CG_CONSTRAINED if convergence is reached along a constrained step,
174: $  other KSP converged/diverged reasons


177:   Notes:
178:   The preconditioner supplied should be symmetric and positive definite.
179: */

181: PetscErrorCode KSPSolve_GLTR(KSP ksp)
182: {
183: #ifdef PETSC_USE_COMPLEX
184:   SETERRQ(PETSC_ERR_SUP, "GLTR is not available for complex systems");
185: #else
186:   KSP_GLTR *cg = (KSP_GLTR *)ksp->data;
187:   PetscReal *t_soln, *t_diag, *t_offd, *e_valu, *e_vect, *e_rwrk;
188:   PetscInt  *e_iblk, *e_splt, *e_iwrk;

191:   MatStructure  pflag;
192:   Mat Qmat, Mmat;
193:   Vec r, z, p, d;
194:   PC  pc;
195:   PetscReal norm_r, norm_d, norm_dp1, norm_p, dMp;
196:   PetscReal alpha, beta, kappa, rz, rzm1;
197:   PetscReal rr, r2, piv, step;
198:   PetscReal vl, vu;
199:   PetscReal coef1, coef2, coef3, root1, root2, obj1, obj2;
200:   PetscReal norm_t, norm_w, pert;
201:   PetscInt  i, j, max_cg_its, max_lanczos_its, max_newton_its, sigma;
202:   PetscInt  t_size = 0, il, iu, e_valus, e_splts, info;
203:   PetscInt  nrhs, nldb;

205:   KSPConvergedReason reason;
206:   PetscTruth diagonalscale;


210:   /***************************************************************************/
211:   /* Check the arguments and parameters.                                     */
212:   /***************************************************************************/

214:   PCDiagonalScale(ksp->pc, &diagonalscale);
215:   if (diagonalscale) {
216:     SETERRQ1(PETSC_ERR_SUP, "Krylov method %s does not support diagonal scaling", ksp->type_name);
217:   }

219:   if (cg->radius < 0.0) {
220:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Input error: radius < 0");
221:   }

223:   /***************************************************************************/
224:   /* Get the workspace vectors and initialize variables                      */
225:   /***************************************************************************/

227:   r2 = cg->radius * cg->radius;
228:   r  = ksp->work[0];
229:   z  = ksp->work[1];
230:   p  = ksp->work[2];
231:   d  = ksp->vec_sol;
232:   pc = ksp->pc;

234:   PCGetOperators(pc, &Qmat, &Mmat, &pflag);

236:   max_cg_its      = cg->max_cg_its;
237:   max_lanczos_its = cg->max_lanczos_its;
238:   max_newton_its  = cg->max_newton_its;
239:   ksp->its = 0;

241:   /***************************************************************************/
242:   /* Initialize objective function value and minimum eigenvalue.             */
243:   /***************************************************************************/

245:   cg->e_min = 0.0;
246:   cg->o_fcn = 0.0;
247:   cg->lambda = 0.0;

249:   /***************************************************************************/
250:   /* The first phase of GLTR performs a standard conjugate gradient method,  */
251:   /* but stores the values required for the Lanczos matrix.  We switch to    */
252:   /* the Lanczos when the conjugate gradient method breaks down.             */
253:   /***************************************************************************/

255:   VecSet(d, 0.0);                         /* d = 0             */
256:   VecCopy(ksp->vec_rhs, r);         /* r = -grad         */
257:   KSP_PCApply(ksp, r, z);                 /* z = inv(M) r      */
258:   cg->norm_d = 0.0;

260:   /***************************************************************************/
261:   /* Check the preconditioners for numerical problems and for positive       */
262:   /* definiteness.  The check for not-a-number and infinite values need be   */
263:   /* performed only once.                                                    */
264:   /***************************************************************************/

266:   VecDot(r, z, &rz);                 /* rz = r^T inv(M) r */
267:   if ((rz != rz) || (rz && (rz / rz != rz / rz))) {
268:     /*************************************************************************/
269:     /* The preconditioner produced not-a-number or an infinite value.  This  */
270:     /* can appear either due to r having numerical problems or M having      */
271:     /* numerical problems.  Differentiate between the two and then use the   */
272:     /* gradient direction.                                                   */
273:     /*************************************************************************/

275:     ksp->reason = KSP_DIVERGED_NAN;
276:     VecDot(r, r, &rr);                 /* rr = r^T r        */
277:     if ((rr != rr) || (rr && (rr / rr != rr / rr))) {
278:       /***********************************************************************/
279:       /* The right-hand side contains not-a-number or an infinite value.     */
280:       /* The gradient step does not work; return a zero value for the step.  */
281:       /***********************************************************************/

283:       PetscInfo1(ksp, "KSPSolve_GLTR: bad right-hand side: rr=%g\n", rr);
284:     }
285:     else {
286:       /***********************************************************************/
287:       /* The preconditioner contains not-a-number or an infinite value.      */
288:       /***********************************************************************/

290:       PetscInfo1(ksp, "KSPSolve_GLTR: bad preconditioner: rz=%g\n", rz);

292:       if (cg->radius) {
293:         alpha = sqrt(r2 / rr);
294:         VecAXPY(d, alpha, r);         /* d = d + alpha r   */
295:         cg->norm_d = cg->radius;
296: 
297:         /*********************************************************************/
298:         /* Compute objective function.                                       */
299:         /*********************************************************************/

301:         KSP_MatMult(ksp, Qmat, d, z); CHKERRQ(ierr)
302:         VecAYPX(z, -0.5, ksp->vec_rhs);
303:         VecDot(d, z, &cg->o_fcn);
304:         cg->o_fcn = -cg->o_fcn;
305:       }
306:     }
307:     return(0);
308:   }

310:   if (rz < 0.0) {
311:     /*************************************************************************/
312:     /* The preconditioner is indefinite.  Because this is the first          */
313:     /* and we do not have a direction yet, we use the gradient step.  Note   */
314:     /* that we cannot use the preconditioned norm when computing the step    */
315:     /* because the matrix is indefinite.                                     */
316:     /*************************************************************************/

318:     ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
319:     PetscInfo1(ksp, "KSPSolve_GLTR: indefinite preconditioner: rz=%g\n", rz);

321:     if (cg->radius) {
322:       VecDot(r, r, &rr);                 /* rr = r^T r        */
323:       alpha = sqrt(r2 / rr);
324:       VecAXPY(d, alpha, r);         /* d = d + alpha r   */
325:       cg->norm_d = cg->radius;

327:       /***********************************************************************/
328:       /* Compute objective function.                                         */
329:       /***********************************************************************/

331:       KSP_MatMult(ksp, Qmat, d, z); CHKERRQ(ierr)
332:       VecAYPX(z, -0.5, ksp->vec_rhs);
333:       VecDot(d, z, &cg->o_fcn);
334:       cg->o_fcn = -cg->o_fcn;
335:     }
336:     return(0);
337:   }

339:   /***************************************************************************/
340:   /* As far as we know, the preconditioner is positive semidefinite.  Compute*/
341:   /* the residual and check for convergence.  Note that there is no choice   */
342:   /* in the residual that can be used; we must always use the natural        */
343:   /* residual because this is the only one that can be used by the           */
344:   /* preconditioned Lanzcos method.                                          */
345:   /***************************************************************************/

347:   norm_r = sqrt(rz);                                        /* norm_r = |r|_M    */
348:   cg->norm_r[0] = norm_r;

350:   KSPLogResidualHistory(ksp, norm_r);
351:   KSPMonitor(ksp, 0, norm_r);
352:   ksp->rnorm = norm_r;

354:   (*ksp->converged)(ksp, 0, norm_r, &ksp->reason, ksp->cnvP);
355:   if (ksp->reason) {
356:     return(0);
357:   }

359:   /***************************************************************************/
360:   /* We have not converged.  Compute the first direction and check the       */
361:   /* matrix for numerical problems.                                          */
362:   /***************************************************************************/

364:   VecCopy(z, p);                         /* p = z             */
365:   KSP_MatMult(ksp, Qmat, p, z);         /* z = Q * p         */
366:   VecDot(p, z, &kappa);                 /* kappa = p^T Q p   */
367:   if ((kappa != kappa) || (kappa && (kappa / kappa != kappa / kappa))) {
368:     /*************************************************************************/
369:     /* The matrix produced not-a-number or an infinite value.  In this case, */
370:     /* we must stop and use the gradient direction.  This condition need     */
371:     /* only be checked once.                                                 */
372:     /*************************************************************************/

374:     ksp->reason = KSP_DIVERGED_NAN;
375:     PetscInfo1(ksp, "KSPSolve_GLTR: bad matrix: kappa=%g\n", kappa);

377:     if (cg->radius) {
378:       VecDot(r, r, &rr);                 /* rr = r^T r        */
379:       alpha = sqrt(r2 / rr);
380:       VecAXPY(d, alpha, r);         /* d = d + alpha r   */
381:       cg->norm_d = cg->radius;

383:       /***********************************************************************/
384:       /* Compute objective function.                                         */
385:       /***********************************************************************/

387:       KSP_MatMult(ksp, Qmat, d, z); CHKERRQ(ierr)
388:       VecAYPX(z, -0.5, ksp->vec_rhs);
389:       VecDot(d, z, &cg->o_fcn);
390:       cg->o_fcn = -cg->o_fcn;
391:     }
392:     return(0);
393:   }

395:   /***************************************************************************/
396:   /* Check for breakdown of the conjugate gradient method; this occurs when  */
397:   /* kappa is zero.                                                          */
398:   /***************************************************************************/

400:   if (fabs(kappa) <= 0.0) {
401:     /*************************************************************************/
402:     /* The curvature is zero.  In this case, we must stop and use the        */
403:     /* since there the Lanczos matrix is not available.                      */
404:     /*************************************************************************/
405:     ksp->reason = KSP_DIVERGED_BREAKDOWN;
406:     PetscInfo1(ksp, "KSPSolve_GLTR: breakdown: kappa=%g\n", kappa);

408:     if (cg->radius) {
409:       VecDot(r, r, &rr);                 /* rr = r^T r        */
410:       alpha = sqrt(r2 / rr);
411:       VecAXPY(d, alpha, r);         /* d = d + alpha r   */
412:       cg->norm_d = cg->radius;

414:       /***********************************************************************/
415:       /* Compute objective function.                                         */
416:       /***********************************************************************/

418:       KSP_MatMult(ksp, Qmat, d, z); CHKERRQ(ierr)
419:       VecAYPX(z, -0.5, ksp->vec_rhs);
420:       VecDot(d, z, &cg->o_fcn);
421:       cg->o_fcn = -cg->o_fcn;
422:     }
423:     return(0);
424:   }

426:   /***************************************************************************/
427:   /* Initialize variables for calculating the norm of the direction and for  */
428:   /* the Lanczos tridiagonal matrix.  Note that we track the diagonal value  */
429:   /* of the Cholesky factorization of the Lanczos matrix in order to         */
430:   /* determine when negative curvature is encountered.                       */
431:   /***************************************************************************/

433:   dMp = 0.0;
434:   norm_p = rz;
435:   norm_d = 0.0;

437:   cg->diag[t_size] = kappa / rz;
438:   cg->offd[t_size] = 0.0;
439:   ++t_size;

441:   piv = 1.0;

443:   /***************************************************************************/
444:   /* Begin the first part of the GLTR algorithm which runs the conjugate     */
445:   /* gradient method until either the problem is solved, we encounter the    */
446:   /* boundary of the trust region, or the conjugate gradient method breaks   */
447:   /* down.                                                                   */
448:   /***************************************************************************/

450:   for (i = 0; i <= max_cg_its; ++i) {
451:     /*************************************************************************/
452:     /* Know that kappa is nonzero, because we have not broken down, so we    */
453:     /* can compute the steplength.                                           */
454:     /*************************************************************************/

456:     alpha = rz / kappa;
457:     cg->alpha[ksp->its] = alpha;

459:     /*************************************************************************/
460:     /* Compute the diagonal value of the Cholesky factorization of the       */
461:     /* Lanczos matrix and check to see if the Lanczos matrix is indefinite.  */
462:     /* This indicates a direction of negative curvature.                     */
463:     /*************************************************************************/

465:     piv = cg->diag[ksp->its] - cg->offd[ksp->its]*cg->offd[ksp->its] / piv;
466:     if (piv <= 0.0) {
467:       /***********************************************************************/
468:       /* In this case, the matrix is indefinite and we have encountered      */
469:       /* a direction of negative curvature.  Follow the direction to the     */
470:       /* boundary of the trust region.                                       */
471:       /***********************************************************************/

473:       ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
474:       PetscInfo1(ksp, "KSPSolve_GLTR: negative curvature: kappa=%g\n", kappa);

476:       if (cg->radius) {
477:         step = (sqrt(dMp*dMp+norm_p*(r2-norm_d))-dMp)/norm_p;
478:         VecAXPY(d, step, p);         /* d = d + step p    */
479:         cg->norm_d = cg->radius;

481:         /*********************************************************************/
482:         /* Update objective function.                                        */
483:         /*********************************************************************/

485:         cg->o_fcn += step * (0.5 * step * kappa - rz);
486:       }
487:       break;
488:     }

490:     /*************************************************************************/
491:     /* Compute the steplength and check for intersection with the trust      */
492:     /* region.                                                               */
493:     /*************************************************************************/

495:     norm_dp1 = norm_d + alpha*(2.0*dMp + alpha*norm_p);
496:     if (cg->radius && norm_dp1 >= r2) {
497:       /***********************************************************************/
498:       /* In this case, the matrix is positive definite as far as we know.    */
499:       /* However, the direction does beyond the trust region.  Follow the    */
500:       /* direction to the boundary of the trust region.                      */
501:       /***********************************************************************/

503:       ksp->reason = KSP_CONVERGED_CG_CONSTRAINED;
504:       PetscInfo1(ksp, "KSPSolve_GLTR: constrained step: radius=%g\n", cg->radius);

506:       step = (sqrt(dMp*dMp+norm_p*(r2-norm_d))-dMp)/norm_p;
507:       VecAXPY(d, step, p);         /* d = d + step p    */
508:       cg->norm_d = cg->radius;

510:       /***********************************************************************/
511:       /* Update objective function.                                          */
512:       /***********************************************************************/

514:       cg->o_fcn += step * (0.5 * step * kappa - rz);
515:       break;
516:     }

518:     /*************************************************************************/
519:     /* Now we can update the direction, residual, and objective value.       */
520:     /*************************************************************************/

522:     VecAXPY(d, alpha, p);                 /* d = d + alpha p   */
523:     VecAXPY(r, -alpha, z);                        /* r = r - alpha Q p */
524:     KSP_PCApply(ksp, r, z);

526:     norm_d = norm_dp1;
527:     cg->norm_d = sqrt(norm_d);

529:     cg->o_fcn -= 0.5 * alpha * rz;

531:     /*************************************************************************/
532:     /* Check that the preconditioner appears positive definite.              */
533:     /*************************************************************************/

535:     rzm1 = rz;
536:     VecDot(r, z, &rz);                 /* rz = r^T z        */
537:     if (rz < 0.0) {
538:       /***********************************************************************/
539:       /* The preconditioner is indefinite.                                   */
540:       /***********************************************************************/

542:       ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
543:       PetscInfo1(ksp, "KSPSolve_GLTR: cg indefinite preconditioner: rz=%g\n", rz);
544:       break;
545:     }

547:     /*************************************************************************/
548:     /* As far as we know, the preconditioner is positive definite.  Compute  */
549:     /* the residual and check for convergence.  Note that there is no choice */
550:     /* in the residual that can be used; we must always use the natural      */
551:     /* residual because this is the only one that can be used by the         */
552:     /* preconditioned Lanzcos method.                                        */
553:     /*************************************************************************/

555:     norm_r = sqrt(rz);                                        /* norm_r = |r|_M   */
556:     cg->norm_r[ksp->its+1] = norm_r;

558:     KSPLogResidualHistory(ksp, norm_r);
559:     KSPMonitor(ksp, 0, norm_r);
560:     ksp->rnorm = norm_r;
561: 
562:     (*ksp->converged)(ksp, i+1, norm_r, &ksp->reason, ksp->cnvP);
563:     if (ksp->reason) {
564:       PetscInfo2(ksp,"KSPSolve_GLTR: cg truncated step: rnorm=%g, radius=%g\n",norm_r,cg->radius);
565:       break;
566:     }

568:     /*************************************************************************/
569:     /* We have not converged yet, update p and the norms.                    */
570:     /*************************************************************************/

572:     beta = rz / rzm1;
573:     cg->beta[ksp->its] = beta;
574:     VecAYPX(p, beta, z);           /* p = z + beta p    */

576:     dMp = beta*(dMp + alpha*norm_p);
577:     norm_p = beta*(rzm1 + beta*norm_p);

579:     /*************************************************************************/
580:     /* Compute the new direction                                             */
581:     /*************************************************************************/

583:     KSP_MatMult(ksp, Qmat, p, z);         /* z = Q * p         */
584:     VecDot(p, z, &kappa);                 /* kappa = p^T Q p   */

586:     /*************************************************************************/
587:     /* Update the iteration and the Lanczos tridiagonal matrix.              */
588:     /*************************************************************************/

590:     ++ksp->its;

592:     cg->offd[t_size] = sqrt(beta) / fabs(alpha);
593:     cg->diag[t_size] = kappa / rz + beta / alpha;
594:     ++t_size;

596:     /*************************************************************************/
597:     /* Check for breakdown of the conjugate gradient method; this occurs     */
598:     /* when kappa is zero.                                                   */
599:     /*************************************************************************/

601:     if (fabs(kappa) <= 0.0) {
602:       /***********************************************************************/
603:       /* The method breaks down; move along the direction as if the matrix   */
604:       /* were indefinite.                                                    */
605:       /***********************************************************************/

607:       ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
608:       PetscInfo1(ksp, "KSPSolve_GLTR: cg breakdown: kappa=%g\n", kappa);

610:       if (cg->radius) {
611:         step = (sqrt(dMp*dMp+norm_p*(r2-norm_d))-dMp)/norm_p;
612:         VecAXPY(d, step, p);         /* d = d + step p    */
613:         cg->norm_d = cg->radius;
614:       }
615:       break;
616:     }
617:   }

619:   if (!ksp->reason) {
620:     ksp->reason = KSP_DIVERGED_ITS;
621:   }

623:   /***************************************************************************/
624:   /* Check to see if we need to continue with the Lanczos method.            */
625:   /***************************************************************************/

627:   if (!cg->radius) {
628:     /*************************************************************************/
629:     /* There is no radius.  Therefore, we cannot move along the boundary.    */
630:     /*************************************************************************/

632:     return(0);
633:   }

635:   if ((KSP_CONVERGED_CG_CONSTRAINED != ksp->reason) &&
636:       (KSP_CONVERGED_CG_NEG_CURVE != ksp->reason)) {
637:     /*************************************************************************/
638:     /* The method either converged to an interior point or the iteration     */
639:     /* limit was reached.                                                    */
640:     /*************************************************************************/

642:     return(0);
643:   }
644:   reason = ksp->reason;

646:   /***************************************************************************/
647:   /* Switch to contructing the Lanczos basis by way of the conjugate         */
648:   /* directions.                                                             */
649:   /***************************************************************************/
650:   for (i = 0; i < max_lanczos_its; ++i) {
651:     /*************************************************************************/
652:     /* Check for breakdown of the conjugate gradient method; this occurs     */
653:     /* when kappa is zero.                                                   */
654:     /*************************************************************************/

656:     if (fabs(kappa) <= 0.0) {
657:       PetscInfo1(ksp, "KSPSolve_GLTR: lanczos breakdown: kappa=%g\n", kappa);
658:       break;
659:     }

661:     /*************************************************************************/
662:     /* Update the direction and residual.                                    */
663:     /*************************************************************************/
664: 
665:     alpha = rz / kappa;
666:     cg->alpha[ksp->its] = alpha;

668:     VecAXPY(r, -alpha, z);                        /* r = r - alpha Q p */
669:     KSP_PCApply(ksp, r, z);

671:     /*************************************************************************/
672:     /* Check that the preconditioner appears positive definite.              */
673:     /*************************************************************************/

675:     rzm1 = rz;
676:     VecDot(r, z, &rz);                 /* rz = r^T z        */
677:     if (rz < 0.0) {
678:       /***********************************************************************/
679:       /* The preconditioner is indefinite.                                   */
680:       /***********************************************************************/

682:       ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
683:       PetscInfo1(ksp, "KSPSolve_GLTR: lanczos indefinite preconditioner: rz=%g\n", rz);
684:       break;
685:     }

687:     /*************************************************************************/
688:     /* As far as we know, the preconditioner is positive definite.  Compute  */
689:     /* the residual and check for convergence.  Note that there is no choice */
690:     /* in the residual that can be used; we must always use the natural      */
691:     /* residual because this is the only one that can be used by the         */
692:     /* preconditioned Lanzcos method.                                        */
693:     /*************************************************************************/

695:     norm_r = sqrt(rz);                                        /* norm_r = |r|_M   */
696:     cg->norm_r[ksp->its+1] = norm_r;

698:     KSPLogResidualHistory(ksp, norm_r);
699:     KSPMonitor(ksp, 0, norm_r);
700:     ksp->rnorm = norm_r;
701: 
702:     (*ksp->converged)(ksp, i+1, norm_r, &ksp->reason, ksp->cnvP);
703:     if (ksp->reason) {
704:       PetscInfo2(ksp,"KSPSolve_GLTR: lanczos truncated step: rnorm=%g, radius=%g\n",norm_r,cg->radius);
705:       break;
706:     }

708:     /*************************************************************************/
709:     /* Update p and the norms.                                               */
710:     /*************************************************************************/

712:     beta = rz / rzm1;
713:     cg->beta[ksp->its] = beta;
714:     VecAYPX(p, beta, z);           /* p = z + beta p    */

716:     /*************************************************************************/
717:     /* Compute the new direction                                             */
718:     /*************************************************************************/

720:     KSP_MatMult(ksp, Qmat, p, z);         /* z = Q * p         */
721:     VecDot(p, z, &kappa);                 /* kappa = p^T Q p   */

723:     /*************************************************************************/
724:     /* Update the iteration and the Lanczos tridiagonal matrix.              */
725:     /*************************************************************************/

727:     ++ksp->its;

729:     cg->offd[t_size] = sqrt(beta) / fabs(alpha);
730:     cg->diag[t_size] = kappa / rz + beta / alpha;
731:     ++t_size;

733:     norm_r = sqrt(rz);                                        /* norm_r = |r|_M   */
734:   }

736:   /***************************************************************************/
737:   /* We have the Lanczos basis, solve the tridiagonal trust-region problem   */
738:   /* to obtain the Lanczos direction.  We know that the solution lies on     */
739:   /* the boundary of the trust region.  We start by checking that the        */
740:   /* workspace allocated is large enough.                                    */
741:   /***************************************************************************/

743:   if (t_size > cg->alloced) {
744:     if (cg->alloced) {
745:       PetscFree(cg->rwork);
746:       PetscFree(cg->iwork);
747:       cg->alloced += cg->init_alloc;
748:     }
749:     else {
750:       cg->alloced = cg->init_alloc;
751:     }

753:     while (t_size > cg->alloced) {
754:       cg->alloced += cg->init_alloc;
755:     }
756: 
757:     cg->alloced = PetscMin(cg->alloced, t_size);
758:     PetscMalloc(10*sizeof(PetscReal)*cg->alloced, &cg->rwork);
759:     PetscMalloc(5*sizeof(PetscInt)*cg->alloced, &cg->iwork);
760:   }

762:   /***************************************************************************/
763:   /* Set up the required vectors.                                            */
764:   /***************************************************************************/

766:   t_soln = cg->rwork + 0*t_size;                        /* Solution          */
767:   t_diag = cg->rwork + 1*t_size;                        /* Diagonal of T     */
768:   t_offd = cg->rwork + 2*t_size;                        /* Off-diagonal of T */
769:   e_valu = cg->rwork + 3*t_size;                        /* Eigenvalues of T  */
770:   e_vect = cg->rwork + 4*t_size;                        /* Eigenvector of T  */
771:   e_rwrk = cg->rwork + 5*t_size;                        /* Eigen workspace   */
772: 
773:   e_iblk = cg->iwork + 0*t_size;                        /* Eigen blocks      */
774:   e_splt = cg->iwork + 1*t_size;                        /* Eigen splits      */
775:   e_iwrk = cg->iwork + 2*t_size;                        /* Eigen workspace   */

777:   /***************************************************************************/
778:   /* Compute the minimum eigenvalue of T.                                    */
779:   /***************************************************************************/

781:   vl = 0.0;
782:   vu = 0.0;
783:   il = 1;
784:   iu = 1;

786:   LAPACKstebz_("I", "E", &t_size, &vl, &vu, &il, &iu, &cg->eigen_tol,
787:                cg->diag, cg->offd + 1, &e_valus, &e_splts, e_valu,
788:                e_iblk, e_splt, e_rwrk, e_iwrk, &info);

790:   if ((0 != info) || (1 != e_valus)) {
791:     /*************************************************************************/
792:     /* Calculation of the minimum eigenvalue failed.  Return the             */
793:     /* Steihaug-Toint direction.                                             */
794:     /*************************************************************************/

796:     PetscInfo(ksp, "KSPSolve_GLTR: failed to compute eigenvalue.\n");
797:     ksp->reason = reason;
798:     return(0);
799:   }

801:   cg->e_min = e_valu[0];

803:   /***************************************************************************/
804:   /* Compute the initial value of lambda to make (T + lamba I) positive      */
805:   /* definite.                                                               */
806:   /***************************************************************************/

808:   pert = cg->init_pert;
809:   if (e_valu[0] <= 0.0) {
810:     cg->lambda = pert - e_valu[0];
811:   }

813:   while(1) {
814:     for (i = 0; i < t_size; ++i) {
815:       t_diag[i] = cg->diag[i] + cg->lambda;
816:       t_offd[i] = cg->offd[i];
817:     }

819:     LAPACKpttrf_(&t_size, t_diag, t_offd + 1, &info);
820:     if (0 == info) {
821:       break;
822:     }

824:     pert += pert;
825:     cg->lambda = cg->lambda * (1.0 + pert) + pert;
826:   }

828:   /***************************************************************************/
829:   /* Compute the initial step and its norm.                                  */
830:   /***************************************************************************/

832:   nrhs = 1;
833:   nldb = t_size;

835:   t_soln[0] = -cg->norm_r[0];
836:   for (i = 1; i < t_size; ++i) {
837:     t_soln[i] = 0.0;
838:   }

840:   LAPACKpttrs_(&t_size, &nrhs, t_diag, t_offd + 1, t_soln, &nldb, &info);
841:   if (0 != info) {
842:     /*************************************************************************/
843:     /* Calculation of the initial step failed; return the Steihaug-Toint     */
844:     /* direction.                                                            */
845:     /*************************************************************************/

847:     PetscInfo(ksp, "KSPSolve_GLTR: failed to compute step.\n");
848:     ksp->reason = reason;
849:     return(0);
850:   }

852:   norm_t = 0;
853:   for (i = 0; i < t_size; ++i) {
854:     norm_t += t_soln[i] * t_soln[i];
855:   }
856:   norm_t = sqrt(norm_t);

858:   /***************************************************************************/
859:   /* Determine the case we are in.                                           */
860:   /***************************************************************************/

862:   if (norm_t <= cg->radius) {
863:     /*************************************************************************/
864:     /* The step is within the trust region; check if we are in the hard case */
865:     /* and need to move to the boundary by following a direction of negative */
866:     /* curvature.                                                            */
867:     /*************************************************************************/
868: 
869:     if ((e_valu[0] <= 0.0) && (norm_t < cg->radius)) {
870:       /***********************************************************************/
871:       /* This is the hard case; compute the eigenvector associated with the  */
872:       /* minimum eigenvalue and move along this direction to the boundary.   */
873:       /***********************************************************************/

875:       LAPACKstein_(&t_size, cg->diag, cg->offd + 1, &e_valus, e_valu,
876:                    e_iblk, e_splt, e_vect, &nldb,
877:                    e_rwrk, e_iwrk, e_iwrk + t_size, &info);
878:       if (0 != info) {
879:         /*********************************************************************/
880:         /* Calculation of the minimum eigenvalue failed.  Return the         */
881:         /* Steihaug-Toint direction.                                         */
882:         /*********************************************************************/
883: 
884:         PetscInfo(ksp, "KSPSolve_GLTR: failed to compute eigenvector.\n");
885:         ksp->reason = reason;
886:         return(0);
887:       }
888: 
889:       coef1 = 0.0;
890:       coef2 = 0.0;
891:       coef3 = -cg->radius * cg->radius;
892:       for (i = 0; i < t_size; ++i) {
893:         coef1 += e_vect[i] * e_vect[i];
894:         coef2 += e_vect[i] * t_soln[i];
895:         coef3 += t_soln[i] * t_soln[i];
896:       }
897: 
898:       coef3 = sqrt(coef2 * coef2 - coef1 * coef3);
899:       root1 = (-coef2 + coef3) / coef1;
900:       root2 = (-coef2 - coef3) / coef1;

902:       /***********************************************************************/
903:       /* Compute objective value for (t_soln + root1 * e_vect)               */
904:       /***********************************************************************/

906:       for (i = 0; i < t_size; ++i) {
907:         e_rwrk[i] = t_soln[i] + root1 * e_vect[i];
908:       }
909: 
910:       obj1 = e_rwrk[0]*(0.5*(cg->diag[0]*e_rwrk[0]+
911:                              cg->offd[1]*e_rwrk[1])+cg->norm_r[0]);
912:       for (i = 1; i < t_size - 1; ++i) {
913:         obj1 += 0.5*e_rwrk[i]*(cg->offd[i]*e_rwrk[i-1]+
914:                                cg->diag[i]*e_rwrk[i]+
915:                                cg->offd[i+1]*e_rwrk[i+1]);
916:       }
917:       obj1 += 0.5*e_rwrk[i]*(cg->offd[i]*e_rwrk[i-1]+
918:                              cg->diag[i]*e_rwrk[i]);

920:       /***********************************************************************/
921:       /* Compute objective value for (t_soln + root2 * e_vect)               */
922:       /***********************************************************************/

924:       for (i = 0; i < t_size; ++i) {
925:         e_rwrk[i] = t_soln[i] + root2 * e_vect[i];
926:       }
927: 
928:       obj2 = e_rwrk[0]*(0.5*(cg->diag[0]*e_rwrk[0]+
929:                              cg->offd[1]*e_rwrk[1])+cg->norm_r[0]);
930:       for (i = 1; i < t_size - 1; ++i) {
931:         obj2 += 0.5*e_rwrk[i]*(cg->offd[i]*e_rwrk[i-1]+
932:                                cg->diag[i]*e_rwrk[i]+
933:                                cg->offd[i+1]*e_rwrk[i+1]);
934:       }
935:       obj2 += 0.5*e_rwrk[i]*(cg->offd[i]*e_rwrk[i-1]+
936:                              cg->diag[i]*e_rwrk[i]);
937: 
938:       /***********************************************************************/
939:       /* Choose the point with the best objective function value.            */
940:       /***********************************************************************/

942:       if (obj1 <= obj2) {
943:         for (i = 0; i < t_size; ++i) {
944:           t_soln[i] += root1 * e_vect[i];
945:         }
946:       }
947:       else {
948:         for (i = 0; i < t_size; ++i) {
949:           t_soln[i] += root2 * e_vect[i];
950:         }
951:       }
952:     }
953:     else {
954:       /***********************************************************************/
955:       /* The matrix is positive definite or there was no room to move; the   */
956:       /* solution is already contained in t_soln.                            */
957:       /***********************************************************************/
958:     }
959:   }
960:   else {
961:     /*************************************************************************/
962:     /* The step is outside the trust-region.  Compute the correct value for  */
963:     /* lambda by performing Newton's method.                                 */
964:     /*************************************************************************/

966:     for (i = 0; i < max_newton_its; ++i) {
967:       /***********************************************************************/
968:       /* Check for convergence.                                              */
969:       /***********************************************************************/

971:       if (fabs(norm_t - cg->radius) <= cg->newton_tol * cg->radius) {
972:         break;
973:       }

975:       /***********************************************************************/
976:       /* Compute the update.                                                 */
977:       /***********************************************************************/

979:       PetscMemcpy(e_rwrk, t_soln, sizeof(PetscReal)*t_size);
980: 
981:       LAPACKpttrs_(&t_size, &nrhs, t_diag, t_offd + 1, e_rwrk, &nldb, &info);
982:       if (0 != info) {
983:         /*********************************************************************/
984:         /* Calculation of the step failed; return the Steihaug-Toint         */
985:         /* direction.                                                        */
986:         /*********************************************************************/

988:         PetscInfo(ksp, "KSPSolve_GLTR: failed to compute step.\n");
989:         ksp->reason = reason;
990:         return(0);
991:       }

993:       /***********************************************************************/
994:       /* Modify lambda.                                                      */
995:       /***********************************************************************/

997:       norm_w = 0;
998:       for (j = 0; j < t_size; ++j) {
999:         norm_w += t_soln[j] * e_rwrk[j];
1000:       }
1001: 
1002:       cg->lambda += (norm_t - cg->radius)/cg->radius * (norm_t * norm_t) / norm_w;

1004:       /***********************************************************************/
1005:       /* Factor T + lambda I                                                 */
1006:       /***********************************************************************/
1007: 
1008:       for (j = 0; j < t_size; ++j) {
1009:         t_diag[j] = cg->diag[j] + cg->lambda;
1010:         t_offd[j] = cg->offd[j];
1011:       }

1013:       LAPACKpttrf_(&t_size, t_diag, t_offd + 1, &info);
1014:       if (0 != info) {
1015:         /*********************************************************************/
1016:         /* Calculation of factorization failed; return the Steihaug-Toint    */
1017:         /* direction.                                                        */
1018:         /*********************************************************************/

1020:         PetscInfo(ksp, "KSPSolve_GLTR: factorization failed.\n");
1021:         ksp->reason = reason;
1022:         return(0);
1023:       }

1025:       /***********************************************************************/
1026:       /* Compute the new step and its norm.                                  */
1027:       /***********************************************************************/

1029:       t_soln[0] = -cg->norm_r[0];
1030:       for (j = 1; j < t_size; ++j) {
1031:         t_soln[j] = 0.0;
1032:       }

1034:       LAPACKpttrs_(&t_size, &nrhs, t_diag, t_offd + 1, t_soln, &nldb, &info);
1035:       if (0 != info) {
1036:         /*********************************************************************/
1037:         /* Calculation of the step failed; return the Steihaug-Toint         */
1038:         /* direction.                                                        */
1039:         /*********************************************************************/
1040: 
1041:         PetscInfo(ksp, "KSPSolve_GLTR: failed to compute step.\n");
1042:         ksp->reason = reason;
1043:         return(0);
1044:       }

1046:       norm_t = 0;
1047:       for (j = 0; j < t_size; ++j) {
1048:         norm_t += t_soln[j] * t_soln[j];
1049:       }
1050:       norm_t = sqrt(norm_t);
1051:     }

1053:     /*************************************************************************/
1054:     /* Check for convergence.                                                */
1055:     /*************************************************************************/

1057:     if (fabs(norm_t - cg->radius) > cg->newton_tol * cg->radius) {
1058:       /***********************************************************************/
1059:       /* Newton method failed to converge in iteration limit.                */
1060:       /***********************************************************************/

1062:       PetscInfo(ksp, "KSPSolve_GLTR: failed to converge.\n");
1063:       ksp->reason = reason;
1064:       return(0);
1065:     }
1066:   }

1068:   /***************************************************************************/
1069:   /* Recover the norm of the direction and objective function value.         */
1070:   /***************************************************************************/

1072:   cg->norm_d = norm_t;

1074:   cg->o_fcn = t_soln[0]*(0.5*(cg->diag[0]*t_soln[0]+
1075:                               cg->offd[1]*t_soln[1])+cg->norm_r[0]);
1076:   for (i = 1; i < t_size - 1; ++i) {
1077:     cg->o_fcn += 0.5*t_soln[i]*(cg->offd[i]*t_soln[i-1]+
1078:                                  cg->diag[i]*t_soln[i]+
1079:                                 cg->offd[i+1]*t_soln[i+1]);
1080:   }
1081:   cg->o_fcn += 0.5*t_soln[i]*(cg->offd[i]*t_soln[i-1]+
1082:                               cg->diag[i]*t_soln[i]);

1084:   /***************************************************************************/
1085:   /* Recover the direction.                                                  */
1086:   /***************************************************************************/

1088:   sigma = -1;

1090:   /***************************************************************************/
1091:   /* Start conjugate gradient method from the beginning                      */
1092:   /***************************************************************************/

1094:   VecCopy(ksp->vec_rhs, r);         /* r = -grad         */
1095:   KSP_PCApply(ksp, r, z);                 /* z = inv(M) r      */

1097:   /***************************************************************************/
1098:   /* Accumulate Q * s                                                        */
1099:   /***************************************************************************/

1101:   VecCopy(z, d);
1102:   VecScale(d, sigma * t_soln[0] / cg->norm_r[0]);
1103: 
1104:   /***************************************************************************/
1105:   /* Compute the first direction.                                            */
1106:   /***************************************************************************/

1108:   VecCopy(z, p);                         /* p = z             */
1109:   KSP_MatMult(ksp, Qmat, p, z);         /* z = Q * p         */

1111:   for (i = 0; i < ksp->its - 1; ++i) {
1112:     /*************************************************************************/
1113:     /* Update the residual and direction.                                    */
1114:     /*************************************************************************/

1116:     alpha = cg->alpha[i];
1117:     if (alpha >= 0.0) {
1118:       sigma = -sigma;
1119:     }

1121:     VecAXPY(r, -alpha, z);                        /* r = r - alpha Q p */
1122:     KSP_PCApply(ksp, r, z);

1124:     /*************************************************************************/
1125:     /* Accumulate Q * s                                                      */
1126:     /*************************************************************************/

1128:     VecAXPY(d, sigma * t_soln[i+1] / cg->norm_r[i+1], z);

1130:     /*************************************************************************/
1131:     /* Update p.                                                             */
1132:     /*************************************************************************/

1134:     beta = cg->beta[i];
1135:     VecAYPX(p, beta, z);           /* p = z + beta p    */
1136:     KSP_MatMult(ksp, Qmat, p, z);         /* z = Q * p         */
1137:   }

1139:   if (i < ksp->its) {
1140:     /*************************************************************************/
1141:     /* Update the residual and direction.                                    */
1142:     /*************************************************************************/

1144:     alpha = cg->alpha[i];
1145:     if (alpha >= 0.0) {
1146:       sigma = -sigma;
1147:     }

1149:     VecAXPY(r, -alpha, z);                        /* r = r - alpha Q p */
1150:     KSP_PCApply(ksp, r, z);

1152:     /*************************************************************************/
1153:     /* Accumulate Q * s                                                      */
1154:     /*************************************************************************/

1156:     VecAXPY(d, sigma * t_soln[i+1] / cg->norm_r[i+1], z);
1157:   }

1159:   /***************************************************************************/
1160:   /* Set the termination reason.                                             */
1161:   /***************************************************************************/

1163:   ksp->reason = reason;
1164:   return(0);
1165: #endif
1166: }

1170: PetscErrorCode KSPSetUp_GLTR(KSP ksp)
1171: {
1172:   KSP_GLTR *cg = (KSP_GLTR *)ksp->data;
1173:   PetscInt size;

1177:   /* This implementation of CG only handles left preconditioning
1178:    * so generate an error otherwise.
1179:    */
1180:   if (ksp->pc_side == PC_RIGHT) {
1181:     SETERRQ(PETSC_ERR_SUP, "No right preconditioning for KSPGLTR");
1182:   }
1183:   else if (ksp->pc_side == PC_SYMMETRIC) {
1184:     SETERRQ(PETSC_ERR_SUP, "No symmetric preconditioning for KSPGLTR");
1185:   }

1187:   /* get work vectors needed by CG */
1188:   KSPDefaultGetWork(ksp, 3);

1190:   /* allocate workspace for Lanczos matrix */
1191:   cg->max_its = cg->max_cg_its + cg->max_lanczos_its + 1;
1192:   size = 5*cg->max_its*sizeof(PetscReal);

1194:   PetscMalloc(size, &cg->diag);
1195:   PetscMemzero(cg->diag, size);
1196:   PetscLogObjectMemory(ksp, size);
1197:   cg->offd   = cg->diag  + cg->max_its;
1198:   cg->alpha  = cg->offd  + cg->max_its;
1199:   cg->beta   = cg->alpha + cg->max_its;
1200:   cg->norm_r = cg->beta  + cg->max_its;
1201:   return(0);
1202: }

1206: PetscErrorCode KSPDestroy_GLTR(KSP ksp)
1207: {
1208:   KSP_GLTR *cg = (KSP_GLTR *)ksp->data;

1212:   PetscFree(cg->diag);
1213:   if (cg->alloced) {
1214:     PetscFree(cg->rwork);
1215:     PetscFree(cg->iwork);
1216:   }

1218:   KSPDefaultDestroy(ksp);
1219:   /* clear composed functions */
1220:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGLTRSetRadius_C","",PETSC_NULL);
1221:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGLTRGetNormD_C","",PETSC_NULL);
1222:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGLTRGetObjFcn_C","",PETSC_NULL);
1223:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGLTRGetMinEig_C","",PETSC_NULL);
1224:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGLTRGetLambda_C","",PETSC_NULL);
1225:   return(0);
1226: }

1231: PetscErrorCode  KSPGLTRSetRadius_GLTR(KSP ksp,PetscReal radius)
1232: {
1233:   KSP_GLTR *cg = (KSP_GLTR *)ksp->data;

1236:   cg->radius = radius;
1237:   return(0);
1238: }

1242: PetscErrorCode  KSPGLTRGetNormD_GLTR(KSP ksp,PetscReal *norm_d)
1243: {
1244:   KSP_GLTR *cg = (KSP_GLTR *)ksp->data;

1247:   *norm_d = cg->norm_d;
1248:   return(0);
1249: }

1253: PetscErrorCode  KSPGLTRGetObjFcn_GLTR(KSP ksp,PetscReal *o_fcn)
1254: {
1255:   KSP_GLTR *cg = (KSP_GLTR *)ksp->data;

1258:   *o_fcn = cg->o_fcn;
1259:   return(0);
1260: }

1264: PetscErrorCode  KSPGLTRGetMinEig_GLTR(KSP ksp,PetscReal *e_min)
1265: {
1266:   KSP_GLTR *cg = (KSP_GLTR *)ksp->data;

1269:   *e_min = cg->e_min;
1270:   return(0);
1271: }

1275: PetscErrorCode  KSPGLTRGetLambda_GLTR(KSP ksp,PetscReal *lambda)
1276: {
1277:   KSP_GLTR *cg = (KSP_GLTR *)ksp->data;

1280:   *lambda = cg->lambda;
1281:   return(0);
1282: }

1287: PetscErrorCode KSPSetFromOptions_GLTR(KSP ksp)
1288: {
1290:   KSP_GLTR *cg = (KSP_GLTR *)ksp->data;

1293:   PetscOptionsHead("KSP GLTR options");

1295:   PetscOptionsReal("-ksp_gltr_radius", "Trust Region Radius", "KSPGLTRSetRadius", cg->radius, &cg->radius, PETSC_NULL);

1297:   PetscOptionsReal("-ksp_gltr_init_pert", "Initial perturbation", "", cg->init_pert, &cg->init_pert, PETSC_NULL);
1298:   PetscOptionsReal("-ksp_gltr_eigen_tol", "Eigenvalue tolerance", "", cg->eigen_tol, &cg->eigen_tol, PETSC_NULL);
1299:   PetscOptionsReal("-ksp_gltr_newton_tol", "Newton tolerance", "", cg->newton_tol, &cg->newton_tol, PETSC_NULL);

1301:   PetscOptionsInt("-ksp_gltr_max_cg_its", "Maximum Conjugate Gradient Iters", "", cg->max_cg_its, &cg->max_cg_its, PETSC_NULL);
1302:   PetscOptionsInt("-ksp_gltr_max_lanczos_its", "Maximum Lanczos Iters", "", cg->max_lanczos_its, &cg->max_lanczos_its, PETSC_NULL);
1303:   PetscOptionsInt("-ksp_gltr_max_newton_its", "Maximum Newton Iters", "", cg->max_newton_its, &cg->max_newton_its, PETSC_NULL);
1304:   PetscOptionsTail();
1305:   return(0);
1306: }

1308: /*MC
1309:      KSPGLTR -   Code to run conjugate gradient method subject to a constraint
1310:          on the solution norm. This is used in Trust Region methods for
1311:          nonlinear equations, SNESTR

1313:    Options Database Keys:
1314: .      -ksp_gltr_radius <r> - Trust Region Radius

1316:    Notes: This is rarely used directly

1318:    Level: developer

1320: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPGLTRSetRadius(), KSPGLTRGetNormD(), KSPGLTRGetObjFcn(), KSPGLTRGetMinEig(), KSPGLTRGetLambda()
1321: M*/

1326: PetscErrorCode  KSPCreate_GLTR(KSP ksp)
1327: {
1329:   KSP_GLTR *cg;


1333:   PetscNew(KSP_GLTR, &cg);
1334:   PetscLogObjectMemory(ksp, sizeof(KSP_GLTR));

1336:   cg->radius = PETSC_MAX;

1338:   cg->init_pert = 1.0e-8;
1339:   cg->eigen_tol = 1.0e-10;
1340:   cg->newton_tol = 1.0e-6;

1342:   cg->alloced = 0;
1343:   cg->init_alloc = 1024;

1345:   cg->max_cg_its = 10000;
1346:   cg->max_lanczos_its = 20;
1347:   cg->max_newton_its  = 10;

1349:   cg->max_its = cg->max_cg_its + cg->max_lanczos_its + 1;

1351:   ksp->data = (void *) cg;
1352:   ksp->pc_side = PC_LEFT;

1354:   /* Sets the functions that are associated with this data structure
1355:    * (in C++ this is the same as defining virtual functions)
1356:    */
1357:   ksp->ops->setup          = KSPSetUp_GLTR;
1358:   ksp->ops->solve          = KSPSolve_GLTR;
1359:   ksp->ops->destroy        = KSPDestroy_GLTR;
1360:   ksp->ops->setfromoptions = KSPSetFromOptions_GLTR;
1361:   ksp->ops->buildsolution  = KSPDefaultBuildSolution;
1362:   ksp->ops->buildresidual  = KSPDefaultBuildResidual;
1363:   ksp->ops->view           = 0;

1365:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,
1366:                                     "KSPGLTRSetRadius_C",
1367:                                     "KSPGLTRSetRadius_GLTR",
1368:                                      KSPGLTRSetRadius_GLTR);
1369:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,
1370:                                     "KSPGLTRGetNormD_C",
1371:                                     "KSPGLTRGetNormD_GLTR",
1372:                                      KSPGLTRGetNormD_GLTR);
1373:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,
1374:                                     "KSPGLTRGetObjFcn_C",
1375:                                     "KSPGLTRGetObjFcn_GLTR",
1376:                                      KSPGLTRGetObjFcn_GLTR);
1377:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,
1378:                                     "KSPGLTRGetMinEig_C",
1379:                                     "KSPGLTRGetMinEig_GLTR",
1380:                                      KSPGLTRGetMinEig_GLTR);
1381:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,
1382:                                     "KSPGLTRGetLambda_C",
1383:                                     "KSPGLTRGetLambda_GLTR",
1384:                                      KSPGLTRGetLambda_GLTR);
1385:   return(0);
1386: }