Actual source code: gltr.c
1: #define PETSCKSP_DLL
3: #include private/kspimpl.h
4: #include ../src/ksp/ksp/impls/cg/gltr/gltrimpl.h
5: #include petscblaslapack.h
7: #define GLTR_PRECONDITIONED_DIRECTION 0
8: #define GLTR_UNPRECONDITIONED_DIRECTION 1
9: #define GLTR_DIRECTION_TYPES 2
11: static const char *DType_Table[64] = {"preconditioned", "unpreconditioned"};
15: /*@
16: KSPGLTRSetRadius - Sets the radius of the trust region.
18: Collective on KSP
20: Input Parameters:
21: + ksp - the iterative context
22: - radius - the trust region radius (Infinity is the default)
24: Options Database Key:
25: . -ksp_gltr_radius <r>
27: Level: advanced
29: .keywords: KSP, GLTR, set, trust region radius
30: @*/
31: PetscErrorCode KSPGLTRSetRadius(KSP ksp, PetscReal radius)
32: {
33: PetscErrorCode ierr, (*f)(KSP, PetscReal);
37: if (radius < 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Radius negative");
38: PetscObjectQueryFunction((PetscObject)ksp, "KSPGLTRSetRadius_C", (void (**)(void))&f);
39: if (f) {
40: (*f)(ksp, radius);
41: }
42: return(0);
43: }
47: /*@
48: KSPGLTRGetNormD - Get norm of the direction.
50: Collective on KSP
52: Input Parameters:
53: + ksp - the iterative context
54: - norm_d - the norm of the direction
56: Level: advanced
58: .keywords: KSP, GLTR, get, norm direction
59: @*/
60: PetscErrorCode KSPGLTRGetNormD(KSP ksp, PetscReal *norm_d)
61: {
62: PetscErrorCode ierr, (*f)(KSP, PetscReal *);
66: PetscObjectQueryFunction((PetscObject)ksp, "KSPGLTRGetNormD_C", (void (**)(void))&f);
67: if (f) {
68: (*f)(ksp, norm_d);
69: }
70: return(0);
71: }
75: /*@
76: KSPGLTRGetObjFcn - Get objective function value.
78: Collective on KSP
80: Input Parameters:
81: + ksp - the iterative context
82: - o_fcn - the objective function value
84: Level: advanced
86: .keywords: KSP, GLTR, get, objective function
87: @*/
88: PetscErrorCode KSPGLTRGetObjFcn(KSP ksp, PetscReal *o_fcn)
89: {
90: PetscErrorCode ierr, (*f)(KSP, PetscReal *);
94: PetscObjectQueryFunction((PetscObject)ksp, "KSPGLTRGetObjFcn_C", (void (**)(void))&f);
95: if (f) {
96: (*f)(ksp, o_fcn);
97: }
98: return(0);
99: }
103: /*@
104: KSPGLTRGetMinEig - Get minimum eigenvalue.
106: Collective on KSP
108: Input Parameters:
109: + ksp - the iterative context
110: - e_min - the minimum eigenvalue
112: Level: advanced
114: .keywords: KSP, GLTR, get, minimum eigenvalue
115: @*/
116: PetscErrorCode KSPGLTRGetMinEig(KSP ksp, PetscReal *e_min)
117: {
118: PetscErrorCode ierr, (*f)(KSP, PetscReal *);
122: PetscObjectQueryFunction((PetscObject)ksp, "KSPGLTRGetMinEig_C", (void (**)(void))&f);
123: if (f) {
124: (*f)(ksp, e_min);
125: }
126: return(0);
127: }
131: /*@
132: KSPGLTRGetLambda - Get multiplier on trust-region constraint.
134: Collective on KSP
136: Input Parameters:
137: + ksp - the iterative context
138: - lambda - the multiplier
140: Level: advanced
142: .keywords: KSP, GLTR, get, multiplier
143: @*/
144: PetscErrorCode KSPGLTRGetLambda(KSP ksp, PetscReal *lambda)
145: {
146: PetscErrorCode ierr, (*f)(KSP, PetscReal *);
150: PetscObjectQueryFunction((PetscObject)ksp, "KSPGLTRGetLambda_C", (void (**)(void))&f);
151: if (f) {
152: (*f)(ksp, lambda);
153: }
154: return(0);
155: }
159: PetscErrorCode KSPSolve_GLTR(KSP ksp)
160: {
161: #ifdef PETSC_USE_COMPLEX
162: SETERRQ(PETSC_ERR_SUP, "GLTR is not available for complex systems");
163: #else
164: KSP_GLTR *cg = (KSP_GLTR *)ksp->data;
165: PetscReal *t_soln, *t_diag, *t_offd, *e_valu, *e_vect, *e_rwrk;
166: PetscBLASInt *e_iblk, *e_splt, *e_iwrk;
169: MatStructure pflag;
170: Mat Qmat, Mmat;
171: Vec r, z, p, d;
172: PC pc;
174: PetscReal norm_r, norm_d, norm_dp1, norm_p, dMp;
175: PetscReal alpha, beta, kappa, rz, rzm1;
176: PetscReal rr, r2, piv, step;
177: PetscReal vl, vu;
178: PetscReal coef1, coef2, coef3, root1, root2, obj1, obj2;
179: PetscReal norm_t, norm_w, pert;
181: PetscInt i, j, max_cg_its, max_lanczos_its, max_newton_its, sigma;
182: PetscBLASInt t_size = 0, l_size = 0, il, iu, e_valus, info;
183: PetscBLASInt nrhs, nldb;
185: #if !defined(PETSC_MISSING_LAPACK_STEBZ)
186: PetscBLASInt e_splts;
187: #endif
188: PetscTruth diagonalscale;
191: /***************************************************************************/
192: /* Check the arguments and parameters. */
193: /***************************************************************************/
195: PCDiagonalScale(ksp->pc, &diagonalscale);
196: if (diagonalscale) {
197: SETERRQ1(PETSC_ERR_SUP, "Krylov method %s does not support diagonal scaling", ((PetscObject)ksp)->type_name);
198: }
200: if (cg->radius < 0.0) {
201: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Input error: radius < 0");
202: }
204: /***************************************************************************/
205: /* Get the workspace vectors and initialize variables */
206: /***************************************************************************/
208: r2 = cg->radius * cg->radius;
209: r = ksp->work[0];
210: z = ksp->work[1];
211: p = ksp->work[2];
212: d = ksp->vec_sol;
213: pc = ksp->pc;
215: PCGetOperators(pc, &Qmat, &Mmat, &pflag);
217: VecGetSize(d, &max_cg_its);
218: max_cg_its = PetscMin(max_cg_its, ksp->max_it);
219: max_lanczos_its = cg->max_lanczos_its;
220: max_newton_its = cg->max_newton_its;
221: ksp->its = 0;
223: /***************************************************************************/
224: /* Initialize objective function direction, and minimum eigenvalue. */
225: /***************************************************************************/
227: cg->o_fcn = 0.0;
229: VecSet(d, 0.0); /* d = 0 */
230: cg->norm_d = 0.0;
232: cg->e_min = 0.0;
233: cg->lambda = 0.0;
235: /***************************************************************************/
236: /* The first phase of GLTR performs a standard conjugate gradient method, */
237: /* but stores the values required for the Lanczos matrix. We switch to */
238: /* the Lanczos when the conjugate gradient method breaks down. Check the */
239: /* right-hand side for numerical problems. The check for not-a-number and */
240: /* infinite values need be performed only once. */
241: /***************************************************************************/
243: VecCopy(ksp->vec_rhs, r); /* r = -grad */
244: VecDot(r, r, &rr); /* rr = r^T r */
245: if PetscIsInfOrNanScalar(rr) {
246: /*************************************************************************/
247: /* The right-hand side contains not-a-number or an infinite value. */
248: /* The gradient step does not work; return a zero value for the step. */
249: /*************************************************************************/
251: ksp->reason = KSP_DIVERGED_NAN;
252: PetscInfo1(ksp, "KSPSolve_GLTR: bad right-hand side: rr=%g\n", rr);
253: return(0);
254: }
256: /***************************************************************************/
257: /* Check the preconditioner for numerical problems and for positive */
258: /* definiteness. The check for not-a-number and infinite values need be */
259: /* performed only once. */
260: /***************************************************************************/
262: KSP_PCApply(ksp, r, z); /* z = inv(M) r */
263: VecDot(r, z, &rz); /* rz = r^T inv(M) r */
264: if PetscIsInfOrNanScalar(rz) {
265: /*************************************************************************/
266: /* The preconditioner contains not-a-number or an infinite value. */
267: /* Return the gradient direction intersected with the trust region. */
268: /*************************************************************************/
270: ksp->reason = KSP_DIVERGED_NAN;
271: PetscInfo1(ksp, "KSPSolve_GLTR: bad preconditioner: rz=%g\n", rz);
273: if (cg->radius) {
274: if (r2 >= rr) {
275: alpha = 1.0;
276: cg->norm_d = sqrt(rr);
277: }
278: else {
279: alpha = sqrt(r2 / rr);
280: cg->norm_d = cg->radius;
281: }
283: VecAXPY(d, alpha, r); /* d = d + alpha r */
284:
285: /***********************************************************************/
286: /* Compute objective function. */
287: /***********************************************************************/
289: KSP_MatMult(ksp, Qmat, d, z);CHKERRQ(ierr)
290: VecAYPX(z, -0.5, ksp->vec_rhs);
291: VecDot(d, z, &cg->o_fcn);
292: cg->o_fcn = -cg->o_fcn;
293: ++ksp->its;
294: }
295: return(0);
296: }
298: if (rz < 0.0) {
299: /*************************************************************************/
300: /* The preconditioner is indefinite. Because this is the first */
301: /* and we do not have a direction yet, we use the gradient step. Note */
302: /* that we cannot use the preconditioned norm when computing the step */
303: /* because the matrix is indefinite. */
304: /*************************************************************************/
306: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
307: PetscInfo1(ksp, "KSPSolve_GLTR: indefinite preconditioner: rz=%g\n", rz);
309: if (cg->radius) {
310: if (r2 >= rr) {
311: alpha = 1.0;
312: cg->norm_d = sqrt(rr);
313: }
314: else {
315: alpha = sqrt(r2 / rr);
316: cg->norm_d = cg->radius;
317: }
319: VecAXPY(d, alpha, r); /* d = d + alpha r */
321: /***********************************************************************/
322: /* Compute objective function. */
323: /***********************************************************************/
325: KSP_MatMult(ksp, Qmat, d, z);CHKERRQ(ierr)
326: VecAYPX(z, -0.5, ksp->vec_rhs);
327: VecDot(d, z, &cg->o_fcn);
328: cg->o_fcn = -cg->o_fcn;
329: ++ksp->its;
330: }
331: return(0);
332: }
334: /***************************************************************************/
335: /* As far as we know, the preconditioner is positive semidefinite. */
336: /* Compute and log the residual. Check convergence because this */
337: /* initializes things, but do not terminate until at least one conjugate */
338: /* gradient iteration has been performed. */
339: /***************************************************************************/
341: cg->norm_r[0] = sqrt(rz); /* norm_r = |r|_M */
343: switch(ksp->normtype) {
344: case KSP_NORM_PRECONDITIONED:
345: VecNorm(z, NORM_2, &norm_r); /* norm_r = |z| */
346: break;
348: case KSP_NORM_UNPRECONDITIONED:
349: norm_r = sqrt(rr); /* norm_r = |r| */
350: break;
352: case KSP_NORM_NATURAL:
353: norm_r = cg->norm_r[0]; /* norm_r = |r|_M */
354: break;
356: default:
357: norm_r = 0.0;
358: break;
359: }
361: KSPLogResidualHistory(ksp, norm_r);
362: KSPMonitor(ksp, ksp->its, norm_r);
363: ksp->rnorm = norm_r;
365: (*ksp->converged)(ksp, ksp->its, norm_r, &ksp->reason, ksp->cnvP);
367: /***************************************************************************/
368: /* Compute the first direction and update the iteration. */
369: /***************************************************************************/
371: VecCopy(z, p); /* p = z */
372: KSP_MatMult(ksp, Qmat, p, z); /* z = Q * p */
373: ++ksp->its;
375: /***************************************************************************/
376: /* Check the matrix for numerical problems. */
377: /***************************************************************************/
379: VecDot(p, z, &kappa); /* kappa = p^T Q p */
380: if PetscIsInfOrNanScalar(kappa) {
381: /*************************************************************************/
382: /* The matrix produced not-a-number or an infinite value. In this case, */
383: /* we must stop and use the gradient direction. This condition need */
384: /* only be checked once. */
385: /*************************************************************************/
387: ksp->reason = KSP_DIVERGED_NAN;
388: PetscInfo1(ksp, "KSPSolve_GLTR: bad matrix: kappa=%g\n", kappa);
390: if (cg->radius) {
391: if (r2 >= rr) {
392: alpha = 1.0;
393: cg->norm_d = sqrt(rr);
394: }
395: else {
396: alpha = sqrt(r2 / rr);
397: cg->norm_d = cg->radius;
398: }
400: VecAXPY(d, alpha, r); /* d = d + alpha r */
402: /***********************************************************************/
403: /* Compute objective function. */
404: /***********************************************************************/
406: KSP_MatMult(ksp, Qmat, d, z);CHKERRQ(ierr)
407: VecAYPX(z, -0.5, ksp->vec_rhs);
408: VecDot(d, z, &cg->o_fcn);
409: cg->o_fcn = -cg->o_fcn;
410: ++ksp->its;
411: }
412: return(0);
413: }
415: /***************************************************************************/
416: /* Initialize variables for calculating the norm of the direction and for */
417: /* the Lanczos tridiagonal matrix. Note that we track the diagonal value */
418: /* of the Cholesky factorization of the Lanczos matrix in order to */
419: /* determine when negative curvature is encountered. */
420: /***************************************************************************/
422: dMp = 0.0;
423: norm_d = 0.0;
424: switch(cg->dtype) {
425: case GLTR_PRECONDITIONED_DIRECTION:
426: norm_p = rz;
427: break;
429: default:
430: VecDot(p, p, &norm_p);
431: break;
432: }
434: cg->diag[t_size] = kappa / rz;
435: cg->offd[t_size] = 0.0;
436: ++t_size;
438: piv = 1.0;
440: /***************************************************************************/
441: /* Check for breakdown of the conjugate gradient method; this occurs when */
442: /* kappa is zero. */
443: /***************************************************************************/
445: if (fabs(kappa) <= 0.0) {
446: /*************************************************************************/
447: /* The curvature is zero. In this case, we must stop and use follow */
448: /* the direction of negative curvature since the Lanczos matrix is zero. */
449: /*************************************************************************/
451: ksp->reason = KSP_DIVERGED_BREAKDOWN;
452: PetscInfo1(ksp, "KSPSolve_GLTR: breakdown: kappa=%g\n", kappa);
454: if (cg->radius && norm_p > 0.0) {
455: /***********************************************************************/
456: /* Follow direction of negative curvature to the boundary of the */
457: /* trust region. */
458: /***********************************************************************/
460: step = sqrt(r2 / norm_p);
461: cg->norm_d = cg->radius;
463: VecAXPY(d, step, p); /* d = d + step p */
465: /***********************************************************************/
466: /* Update objective function. */
467: /***********************************************************************/
469: cg->o_fcn += step * (0.5 * step * kappa - rz);
470: }
471: else if (cg->radius) {
472: /***********************************************************************/
473: /* The norm of the preconditioned direction is zero; use the gradient */
474: /* step. */
475: /***********************************************************************/
477: if (r2 >= rr) {
478: alpha = 1.0;
479: cg->norm_d = sqrt(rr);
480: }
481: else {
482: alpha = sqrt(r2 / rr);
483: cg->norm_d = cg->radius;
484: }
486: VecAXPY(d, alpha, r); /* d = d + alpha r */
488: /***********************************************************************/
489: /* Compute objective function. */
490: /***********************************************************************/
492: KSP_MatMult(ksp, Qmat, d, z);CHKERRQ(ierr)
493: VecAYPX(z, -0.5, ksp->vec_rhs);
494: VecDot(d, z, &cg->o_fcn);
495: cg->o_fcn = -cg->o_fcn;
496: ++ksp->its;
497: }
498: return(0);
499: }
501: /***************************************************************************/
502: /* Begin the first part of the GLTR algorithm which runs the conjugate */
503: /* gradient method until either the problem is solved, we encounter the */
504: /* boundary of the trust region, or the conjugate gradient method breaks */
505: /* down. */
506: /***************************************************************************/
508: while(1) {
509: /*************************************************************************/
510: /* Know that kappa is nonzero, because we have not broken down, so we */
511: /* can compute the steplength. */
512: /*************************************************************************/
514: alpha = rz / kappa;
515: cg->alpha[l_size] = alpha;
517: /*************************************************************************/
518: /* Compute the diagonal value of the Cholesky factorization of the */
519: /* Lanczos matrix and check to see if the Lanczos matrix is indefinite. */
520: /* This indicates a direction of negative curvature. */
521: /*************************************************************************/
523: piv = cg->diag[l_size] - cg->offd[l_size]*cg->offd[l_size] / piv;
524: if (piv <= 0.0) {
525: /***********************************************************************/
526: /* In this case, the matrix is indefinite and we have encountered */
527: /* a direction of negative curvature. Follow the direction to the */
528: /* boundary of the trust region. */
529: /***********************************************************************/
531: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
532: PetscInfo1(ksp, "KSPSolve_GLTR: negative curvature: kappa=%g\n", kappa);
534: if (cg->radius && norm_p > 0.0) {
535: /*********************************************************************/
536: /* Follow direction of negative curvature to boundary. */
537: /*********************************************************************/
539: step = (sqrt(dMp*dMp+norm_p*(r2-norm_d))-dMp)/norm_p;
540: cg->norm_d = cg->radius;
542: VecAXPY(d, step, p); /* d = d + step p */
544: /*********************************************************************/
545: /* Update objective function. */
546: /*********************************************************************/
548: cg->o_fcn += step * (0.5 * step * kappa - rz);
549: }
550: else if (cg->radius) {
551: /*********************************************************************/
552: /* The norm of the direction is zero; there is nothing to follow. */
553: /*********************************************************************/
554: }
555: break;
556: }
558: /*************************************************************************/
559: /* Compute the steplength and check for intersection with the trust */
560: /* region. */
561: /*************************************************************************/
563: norm_dp1 = norm_d + alpha*(2.0*dMp + alpha*norm_p);
564: if (cg->radius && norm_dp1 >= r2) {
565: /***********************************************************************/
566: /* In this case, the matrix is positive definite as far as we know. */
567: /* However, the full step goes beyond the trust region. */
568: /***********************************************************************/
570: ksp->reason = KSP_CONVERGED_CG_CONSTRAINED;
571: PetscInfo1(ksp, "KSPSolve_GLTR: constrained step: radius=%g\n", cg->radius);
573: if (norm_p > 0.0) {
574: /*********************************************************************/
575: /* Follow the direction to the boundary of the trust region. */
576: /*********************************************************************/
578: step = (sqrt(dMp*dMp+norm_p*(r2-norm_d))-dMp)/norm_p;
579: cg->norm_d = cg->radius;
581: VecAXPY(d, step, p); /* d = d + step p */
583: /*********************************************************************/
584: /* Update objective function. */
585: /*********************************************************************/
587: cg->o_fcn += step * (0.5 * step * kappa - rz);
588: }
589: else {
590: /*********************************************************************/
591: /* The norm of the direction is zero; there is nothing to follow. */
592: /*********************************************************************/
593: }
594: break;
595: }
597: /*************************************************************************/
598: /* Now we can update the direction and residual. */
599: /*************************************************************************/
601: VecAXPY(d, alpha, p); /* d = d + alpha p */
602: VecAXPY(r, -alpha, z); /* r = r - alpha Q p */
603: KSP_PCApply(ksp, r, z); /* z = inv(M) r */
605: switch(cg->dtype) {
606: case GLTR_PRECONDITIONED_DIRECTION:
607: norm_d = norm_dp1;
608: break;
610: default:
611: VecDot(d, d, &norm_d);
612: break;
613: }
614: cg->norm_d = sqrt(norm_d);
616: /*************************************************************************/
617: /* Update objective function. */
618: /*************************************************************************/
620: cg->o_fcn -= 0.5 * alpha * rz;
622: /*************************************************************************/
623: /* Check that the preconditioner appears positive semidefinite. */
624: /*************************************************************************/
626: rzm1 = rz;
627: VecDot(r, z, &rz); /* rz = r^T z */
628: if (rz < 0.0) {
629: /***********************************************************************/
630: /* The preconditioner is indefinite. */
631: /***********************************************************************/
633: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
634: PetscInfo1(ksp, "KSPSolve_GLTR: cg indefinite preconditioner: rz=%g\n", rz);
635: break;
636: }
638: /*************************************************************************/
639: /* As far as we know, the preconditioner is positive semidefinite. */
640: /* Compute the residual and check for convergence. */
641: /*************************************************************************/
643: cg->norm_r[l_size+1] = sqrt(rz); /* norm_r = |r|_M */
645: switch(ksp->normtype) {
646: case KSP_NORM_PRECONDITIONED:
647: VecNorm(z, NORM_2, &norm_r);/* norm_r = |z| */
648: break;
650: case KSP_NORM_UNPRECONDITIONED:
651: VecNorm(r, NORM_2, &norm_r);/* norm_r = |r| */
652: break;
654: case KSP_NORM_NATURAL:
655: norm_r = cg->norm_r[l_size+1]; /* norm_r = |r|_M */
656: break;
658: default:
659: norm_r = 0.0;
660: break;
661: }
663: KSPLogResidualHistory(ksp, norm_r);
664: KSPMonitor(ksp, ksp->its, norm_r);
665: ksp->rnorm = norm_r;
666:
667: (*ksp->converged)(ksp, ksp->its, norm_r, &ksp->reason, ksp->cnvP);
668: if (ksp->reason) {
669: /***********************************************************************/
670: /* The method has converged. */
671: /***********************************************************************/
673: PetscInfo2(ksp, "KSPSolve_GLTR: cg truncated step: rnorm=%g, radius=%g\n", norm_r, cg->radius);
674: break;
675: }
677: /*************************************************************************/
678: /* We have not converged yet. Check for breakdown. */
679: /*************************************************************************/
681: beta = rz / rzm1;
682: if (fabs(beta) <= 0.0) {
683: /***********************************************************************/
684: /* Conjugate gradients has broken down. */
685: /***********************************************************************/
687: ksp->reason = KSP_DIVERGED_BREAKDOWN;
688: PetscInfo1(ksp, "KSPSolve_GLTR: breakdown: beta=%g\n", beta);
689: break;
690: }
692: /*************************************************************************/
693: /* Check iteration limit. */
694: /*************************************************************************/
696: if (ksp->its >= max_cg_its) {
697: ksp->reason = KSP_DIVERGED_ITS;
698: PetscInfo1(ksp, "KSPSolve_GLTR: iterlim: its=%d\n", ksp->its);
699: break;
700: }
702: /*************************************************************************/
703: /* Update p and the norms. */
704: /*************************************************************************/
706: cg->beta[l_size] = beta;
707: VecAYPX(p, beta, z); /* p = z + beta p */
709: switch(cg->dtype) {
710: case GLTR_PRECONDITIONED_DIRECTION:
711: dMp = beta*(dMp + alpha*norm_p);
712: norm_p = beta*(rzm1 + beta*norm_p);
713: break;
715: default:
716: VecDot(d, p, &dMp);
717: VecDot(p, p, &norm_p);
718: break;
719: }
721: /*************************************************************************/
722: /* Compute the new direction and update the iteration. */
723: /*************************************************************************/
725: KSP_MatMult(ksp, Qmat, p, z); /* z = Q * p */
726: VecDot(p, z, &kappa); /* kappa = p^T Q p */
727: ++ksp->its;
729: /*************************************************************************/
730: /* Update the the Lanczos tridiagonal matrix. */
731: /*************************************************************************/
733: ++l_size;
734: cg->offd[t_size] = sqrt(beta) / fabs(alpha);
735: cg->diag[t_size] = kappa / rz + beta / alpha;
736: ++t_size;
738: /*************************************************************************/
739: /* Check for breakdown of the conjugate gradient method; this occurs */
740: /* when kappa is zero. */
741: /*************************************************************************/
743: if (fabs(kappa) <= 0.0) {
744: /***********************************************************************/
745: /* The method breaks down; move along the direction as if the matrix */
746: /* were indefinite. */
747: /***********************************************************************/
749: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
750: PetscInfo1(ksp, "KSPSolve_GLTR: cg breakdown: kappa=%g\n", kappa);
752: if (cg->radius && norm_p > 0.0) {
753: /*********************************************************************/
754: /* Follow direction to boundary. */
755: /*********************************************************************/
757: step = (sqrt(dMp*dMp+norm_p*(r2-norm_d))-dMp)/norm_p;
758: cg->norm_d = cg->radius;
760: VecAXPY(d, step, p); /* d = d + step p */
762: /*********************************************************************/
763: /* Update objective function. */
764: /*********************************************************************/
766: cg->o_fcn += step * (0.5 * step * kappa - rz);
767: }
768: else if (cg->radius) {
769: /*********************************************************************/
770: /* The norm of the direction is zero; there is nothing to follow. */
771: /*********************************************************************/
772: }
773: break;
774: }
775: }
777: /***************************************************************************/
778: /* Check to see if we need to continue with the Lanczos method. */
779: /***************************************************************************/
781: if (!cg->radius) {
782: /*************************************************************************/
783: /* There is no radius. Therefore, we cannot move along the boundary. */
784: /*************************************************************************/
786: return(0);
787: }
789: if (KSP_CONVERGED_CG_NEG_CURVE != ksp->reason) {
790: /*************************************************************************/
791: /* The method either converged to an interior point, hit the boundary of */
792: /* the trust-region without encountering a direction of negative */
793: /* curvature or the iteration limit was reached. */
794: /*************************************************************************/
796: return(0);
797: }
799: /***************************************************************************/
800: /* Switch to contructing the Lanczos basis by way of the conjugate */
801: /* directions. */
802: /***************************************************************************/
804: for (i = 0; i < max_lanczos_its; ++i) {
805: /*************************************************************************/
806: /* Check for breakdown of the conjugate gradient method; this occurs */
807: /* when kappa is zero. */
808: /*************************************************************************/
810: if (fabs(kappa) <= 0.0) {
811: ksp->reason = KSP_DIVERGED_BREAKDOWN;
812: PetscInfo1(ksp, "KSPSolve_GLTR: lanczos breakdown: kappa=%g\n", kappa);
813: break;
814: }
816: /*************************************************************************/
817: /* Update the direction and residual. */
818: /*************************************************************************/
819:
820: alpha = rz / kappa;
821: cg->alpha[l_size] = alpha;
823: VecAXPY(r, -alpha, z); /* r = r - alpha Q p */
824: KSP_PCApply(ksp, r, z); /* z = inv(M) r */
826: /*************************************************************************/
827: /* Check that the preconditioner appears positive semidefinite. */
828: /*************************************************************************/
830: rzm1 = rz;
831: VecDot(r, z, &rz); /* rz = r^T z */
832: if (rz < 0.0) {
833: /***********************************************************************/
834: /* The preconditioner is indefinite. */
835: /***********************************************************************/
837: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
838: PetscInfo1(ksp, "KSPSolve_GLTR: lanczos indefinite preconditioner: rz=%g\n", rz);
839: break;
840: }
842: /*************************************************************************/
843: /* As far as we know, the preconditioner is positive definite. Compute */
844: /* the residual. Do NOT check for convergence. */
845: /*************************************************************************/
847: cg->norm_r[l_size+1] = sqrt(rz); /* norm_r = |r|_M */
849: switch(ksp->normtype) {
850: case KSP_NORM_PRECONDITIONED:
851: VecNorm(z, NORM_2, &norm_r);/* norm_r = |z| */
852: break;
854: case KSP_NORM_UNPRECONDITIONED:
855: VecNorm(r, NORM_2, &norm_r);/* norm_r = |r| */
856: break;
858: case KSP_NORM_NATURAL:
859: norm_r = cg->norm_r[l_size+1]; /* norm_r = |r|_M */
860: break;
862: default:
863: norm_r = 0.0;
864: break;
865: }
867: KSPLogResidualHistory(ksp, norm_r);
868: KSPMonitor(ksp, ksp->its, norm_r);
869: ksp->rnorm = norm_r;
870:
871: /*************************************************************************/
872: /* Check for breakdown. */
873: /*************************************************************************/
875: beta = rz / rzm1;
876: if (fabs(beta) <= 0.0) {
877: /***********************************************************************/
878: /* Conjugate gradients has broken down. */
879: /***********************************************************************/
881: ksp->reason = KSP_DIVERGED_BREAKDOWN;
882: PetscInfo1(ksp, "KSPSolve_GLTR: breakdown: beta=%g\n", beta);
883: break;
884: }
886: /*************************************************************************/
887: /* Update p and the norms. */
888: /*************************************************************************/
890: cg->beta[l_size] = beta;
891: VecAYPX(p, beta, z); /* p = z + beta p */
893: /*************************************************************************/
894: /* Compute the new direction and update the iteration. */
895: /*************************************************************************/
897: KSP_MatMult(ksp, Qmat, p, z); /* z = Q * p */
898: VecDot(p, z, &kappa); /* kappa = p^T Q p */
899: ++ksp->its;
901: /*************************************************************************/
902: /* Update the Lanczos tridiagonal matrix. */
903: /*************************************************************************/
905: ++l_size;
906: cg->offd[t_size] = sqrt(beta) / fabs(alpha);
907: cg->diag[t_size] = kappa / rz + beta / alpha;
908: ++t_size;
909: }
911: /***************************************************************************/
912: /* We have the Lanczos basis, solve the tridiagonal trust-region problem */
913: /* to obtain the Lanczos direction. We know that the solution lies on */
914: /* the boundary of the trust region. We start by checking that the */
915: /* workspace allocated is large enough. */
916: /***************************************************************************/
917: /* Note that the current version only computes the solution by using the */
918: /* preconditioned direction. Need to think about how to do the */
919: /* unpreconditioned direction calculation. */
920: /***************************************************************************/
922: if (t_size > cg->alloced) {
923: if (cg->alloced) {
924: PetscFree(cg->rwork);
925: PetscFree(cg->iwork);
926: cg->alloced += cg->init_alloc;
927: }
928: else {
929: cg->alloced = cg->init_alloc;
930: }
932: while (t_size > cg->alloced) {
933: cg->alloced += cg->init_alloc;
934: }
935:
936: cg->alloced = PetscMin(cg->alloced, t_size);
937: PetscMalloc2(10*cg->alloced,PetscReal, &cg->rwork,5*cg->alloced,PetscBLASInt, &cg->iwork);
938: }
940: /***************************************************************************/
941: /* Set up the required vectors. */
942: /***************************************************************************/
944: t_soln = cg->rwork + 0*t_size; /* Solution */
945: t_diag = cg->rwork + 1*t_size; /* Diagonal of T */
946: t_offd = cg->rwork + 2*t_size; /* Off-diagonal of T */
947: e_valu = cg->rwork + 3*t_size; /* Eigenvalues of T */
948: e_vect = cg->rwork + 4*t_size; /* Eigenvector of T */
949: e_rwrk = cg->rwork + 5*t_size; /* Eigen workspace */
950:
951: e_iblk = cg->iwork + 0*t_size; /* Eigen blocks */
952: e_splt = cg->iwork + 1*t_size; /* Eigen splits */
953: e_iwrk = cg->iwork + 2*t_size; /* Eigen workspace */
955: /***************************************************************************/
956: /* Compute the minimum eigenvalue of T. */
957: /***************************************************************************/
959: vl = 0.0;
960: vu = 0.0;
961: il = 1;
962: iu = 1;
964: #if defined(PETSC_MISSING_LAPACK_STEBZ)
965: SETERRQ(PETSC_ERR_SUP,"STEBZ - Lapack routine is unavailable.");
966: #else
967: LAPACKstebz_("I", "E", &t_size, &vl, &vu, &il, &iu, &cg->eigen_tol,
968: cg->diag, cg->offd + 1, &e_valus, &e_splts, e_valu,
969: e_iblk, e_splt, e_rwrk, e_iwrk, &info);
971: if ((0 != info) || (1 != e_valus)) {
972: /*************************************************************************/
973: /* Calculation of the minimum eigenvalue failed. Return the */
974: /* Steihaug-Toint direction. */
975: /*************************************************************************/
977: PetscInfo(ksp, "KSPSolve_GLTR: failed to compute eigenvalue.\n");
978: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
979: return(0);
980: }
982: cg->e_min = e_valu[0];
984: /***************************************************************************/
985: /* Compute the initial value of lambda to make (T + lamba I) positive */
986: /* definite. */
987: /***************************************************************************/
989: pert = cg->init_pert;
990: if (e_valu[0] < 0.0) {
991: cg->lambda = pert - e_valu[0];
992: }
993: #endif
995: while(1) {
996: for (i = 0; i < t_size; ++i) {
997: t_diag[i] = cg->diag[i] + cg->lambda;
998: t_offd[i] = cg->offd[i];
999: }
1001: #if defined(PETSC_MISSING_LAPACK_PTTRF)
1002: SETERRQ(PETSC_ERR_SUP,"PTTRF - Lapack routine is unavailable.");
1003: #else
1004: LAPACKpttrf_(&t_size, t_diag, t_offd + 1, &info);
1006: if (0 == info) {
1007: break;
1008: }
1010: pert += pert;
1011: cg->lambda = cg->lambda * (1.0 + pert) + pert;
1012: #endif
1013: }
1015: /***************************************************************************/
1016: /* Compute the initial step and its norm. */
1017: /***************************************************************************/
1019: nrhs = 1;
1020: nldb = t_size;
1022: t_soln[0] = -cg->norm_r[0];
1023: for (i = 1; i < t_size; ++i) {
1024: t_soln[i] = 0.0;
1025: }
1027: #if defined(PETSC_MISSING_LAPACK_PTTRS)
1028: SETERRQ(PETSC_ERR_SUP,"PTTRS - Lapack routine is unavailable.");
1029: #else
1030: LAPACKpttrs_(&t_size, &nrhs, t_diag, t_offd + 1, t_soln, &nldb, &info);
1031: #endif
1033: if (0 != info) {
1034: /*************************************************************************/
1035: /* Calculation of the initial step failed; return the Steihaug-Toint */
1036: /* direction. */
1037: /*************************************************************************/
1039: PetscInfo(ksp, "KSPSolve_GLTR: failed to compute step.\n");
1040: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
1041: return(0);
1042: }
1044: norm_t = 0.;
1045: for (i = 0; i < t_size; ++i) {
1046: norm_t += t_soln[i] * t_soln[i];
1047: }
1048: norm_t = sqrt(norm_t);
1050: /***************************************************************************/
1051: /* Determine the case we are in. */
1052: /***************************************************************************/
1054: if (norm_t <= cg->radius) {
1055: /*************************************************************************/
1056: /* The step is within the trust region; check if we are in the hard case */
1057: /* and need to move to the boundary by following a direction of negative */
1058: /* curvature. */
1059: /*************************************************************************/
1060:
1061: if ((e_valu[0] <= 0.0) && (norm_t < cg->radius)) {
1062: /***********************************************************************/
1063: /* This is the hard case; compute the eigenvector associated with the */
1064: /* minimum eigenvalue and move along this direction to the boundary. */
1065: /***********************************************************************/
1067: #if defined(PETSC_MISSING_LAPACK_STEIN)
1068: SETERRQ(PETSC_ERR_SUP,"STEIN - Lapack routine is unavailable.");
1069: #else
1070: LAPACKstein_(&t_size, cg->diag, cg->offd + 1, &e_valus, e_valu,
1071: e_iblk, e_splt, e_vect, &nldb,
1072: e_rwrk, e_iwrk, e_iwrk + t_size, &info);
1073: #endif
1075: if (0 != info) {
1076: /*********************************************************************/
1077: /* Calculation of the minimum eigenvalue failed. Return the */
1078: /* Steihaug-Toint direction. */
1079: /*********************************************************************/
1080:
1081: PetscInfo(ksp, "KSPSolve_GLTR: failed to compute eigenvector.\n");
1082: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
1083: return(0);
1084: }
1085:
1086: coef1 = 0.0;
1087: coef2 = 0.0;
1088: coef3 = -cg->radius * cg->radius;
1089: for (i = 0; i < t_size; ++i) {
1090: coef1 += e_vect[i] * e_vect[i];
1091: coef2 += e_vect[i] * t_soln[i];
1092: coef3 += t_soln[i] * t_soln[i];
1093: }
1094:
1095: coef3 = sqrt(coef2 * coef2 - coef1 * coef3);
1096: root1 = (-coef2 + coef3) / coef1;
1097: root2 = (-coef2 - coef3) / coef1;
1099: /***********************************************************************/
1100: /* Compute objective value for (t_soln + root1 * e_vect) */
1101: /***********************************************************************/
1103: for (i = 0; i < t_size; ++i) {
1104: e_rwrk[i] = t_soln[i] + root1 * e_vect[i];
1105: }
1106:
1107: obj1 = e_rwrk[0]*(0.5*(cg->diag[0]*e_rwrk[0]+
1108: cg->offd[1]*e_rwrk[1])+cg->norm_r[0]);
1109: for (i = 1; i < t_size - 1; ++i) {
1110: obj1 += 0.5*e_rwrk[i]*(cg->offd[i]*e_rwrk[i-1]+
1111: cg->diag[i]*e_rwrk[i]+
1112: cg->offd[i+1]*e_rwrk[i+1]);
1113: }
1114: obj1 += 0.5*e_rwrk[i]*(cg->offd[i]*e_rwrk[i-1]+
1115: cg->diag[i]*e_rwrk[i]);
1117: /***********************************************************************/
1118: /* Compute objective value for (t_soln + root2 * e_vect) */
1119: /***********************************************************************/
1121: for (i = 0; i < t_size; ++i) {
1122: e_rwrk[i] = t_soln[i] + root2 * e_vect[i];
1123: }
1124:
1125: obj2 = e_rwrk[0]*(0.5*(cg->diag[0]*e_rwrk[0]+
1126: cg->offd[1]*e_rwrk[1])+cg->norm_r[0]);
1127: for (i = 1; i < t_size - 1; ++i) {
1128: obj2 += 0.5*e_rwrk[i]*(cg->offd[i]*e_rwrk[i-1]+
1129: cg->diag[i]*e_rwrk[i]+
1130: cg->offd[i+1]*e_rwrk[i+1]);
1131: }
1132: obj2 += 0.5*e_rwrk[i]*(cg->offd[i]*e_rwrk[i-1]+
1133: cg->diag[i]*e_rwrk[i]);
1134:
1135: /***********************************************************************/
1136: /* Choose the point with the best objective function value. */
1137: /***********************************************************************/
1139: if (obj1 <= obj2) {
1140: for (i = 0; i < t_size; ++i) {
1141: t_soln[i] += root1 * e_vect[i];
1142: }
1143: }
1144: else {
1145: for (i = 0; i < t_size; ++i) {
1146: t_soln[i] += root2 * e_vect[i];
1147: }
1148: }
1149: }
1150: else {
1151: /***********************************************************************/
1152: /* The matrix is positive definite or there was no room to move; the */
1153: /* solution is already contained in t_soln. */
1154: /***********************************************************************/
1155: }
1156: }
1157: else {
1158: /*************************************************************************/
1159: /* The step is outside the trust-region. Compute the correct value for */
1160: /* lambda by performing Newton's method. */
1161: /*************************************************************************/
1163: for (i = 0; i < max_newton_its; ++i) {
1164: /***********************************************************************/
1165: /* Check for convergence. */
1166: /***********************************************************************/
1168: if (fabs(norm_t - cg->radius) <= cg->newton_tol * cg->radius) {
1169: break;
1170: }
1172: /***********************************************************************/
1173: /* Compute the update. */
1174: /***********************************************************************/
1176: PetscMemcpy(e_rwrk, t_soln, sizeof(PetscReal)*t_size);
1178: #if defined(PETSC_MISSING_LAPACK_PTTRS)
1179: SETERRQ(PETSC_ERR_SUP,"PTTRS - Lapack routine is unavailable.");
1180: #else
1181: LAPACKpttrs_(&t_size, &nrhs, t_diag, t_offd + 1, e_rwrk, &nldb, &info);
1182: #endif
1184: if (0 != info) {
1185: /*********************************************************************/
1186: /* Calculation of the step failed; return the Steihaug-Toint */
1187: /* direction. */
1188: /*********************************************************************/
1190: PetscInfo(ksp, "KSPSolve_GLTR: failed to compute step.\n");
1191: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
1192: return(0);
1193: }
1195: /***********************************************************************/
1196: /* Modify lambda. */
1197: /***********************************************************************/
1199: norm_w = 0.;
1200: for (j = 0; j < t_size; ++j) {
1201: norm_w += t_soln[j] * e_rwrk[j];
1202: }
1203:
1204: cg->lambda += (norm_t - cg->radius)/cg->radius * (norm_t * norm_t) / norm_w;
1206: /***********************************************************************/
1207: /* Factor T + lambda I */
1208: /***********************************************************************/
1209:
1210: for (j = 0; j < t_size; ++j) {
1211: t_diag[j] = cg->diag[j] + cg->lambda;
1212: t_offd[j] = cg->offd[j];
1213: }
1215: #if defined(PETSC_MISSING_LAPACK_PTTRF)
1216: SETERRQ(PETSC_ERR_SUP,"PTTRF - Lapack routine is unavailable.");
1217: #else
1218: LAPACKpttrf_(&t_size, t_diag, t_offd + 1, &info);
1219: #endif
1221: if (0 != info) {
1222: /*********************************************************************/
1223: /* Calculation of factorization failed; return the Steihaug-Toint */
1224: /* direction. */
1225: /*********************************************************************/
1227: PetscInfo(ksp, "KSPSolve_GLTR: factorization failed.\n");
1228: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
1229: return(0);
1230: }
1232: /***********************************************************************/
1233: /* Compute the new step and its norm. */
1234: /***********************************************************************/
1236: t_soln[0] = -cg->norm_r[0];
1237: for (j = 1; j < t_size; ++j) {
1238: t_soln[j] = 0.0;
1239: }
1241: #if defined(PETSC_MISSING_LAPACK_PTTRS)
1242: SETERRQ(PETSC_ERR_SUP,"PTTRS - Lapack routine is unavailable.");
1243: #else
1244: LAPACKpttrs_(&t_size, &nrhs, t_diag, t_offd + 1, t_soln, &nldb, &info);
1245: #endif
1247: if (0 != info) {
1248: /*********************************************************************/
1249: /* Calculation of the step failed; return the Steihaug-Toint */
1250: /* direction. */
1251: /*********************************************************************/
1252:
1253: PetscInfo(ksp, "KSPSolve_GLTR: failed to compute step.\n");
1254: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
1255: return(0);
1256: }
1258: norm_t = 0.;
1259: for (j = 0; j < t_size; ++j) {
1260: norm_t += t_soln[j] * t_soln[j];
1261: }
1262: norm_t = sqrt(norm_t);
1263: }
1265: /*************************************************************************/
1266: /* Check for convergence. */
1267: /*************************************************************************/
1269: if (fabs(norm_t - cg->radius) > cg->newton_tol * cg->radius) {
1270: /***********************************************************************/
1271: /* Newton method failed to converge in iteration limit. */
1272: /***********************************************************************/
1274: PetscInfo(ksp, "KSPSolve_GLTR: failed to converge.\n");
1275: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
1276: return(0);
1277: }
1278: }
1280: /***************************************************************************/
1281: /* Recover the norm of the direction and objective function value. */
1282: /***************************************************************************/
1284: cg->norm_d = norm_t;
1286: cg->o_fcn = t_soln[0]*(0.5*(cg->diag[0]*t_soln[0]+cg->offd[1]*t_soln[1])+cg->norm_r[0]);
1287: for (i = 1; i < t_size - 1; ++i) {
1288: cg->o_fcn += 0.5*t_soln[i]*(cg->offd[i]*t_soln[i-1]+cg->diag[i]*t_soln[i]+cg->offd[i+1]*t_soln[i+1]);
1289: }
1290: cg->o_fcn += 0.5*t_soln[i]*(cg->offd[i]*t_soln[i-1]+cg->diag[i]*t_soln[i]);
1292: /***************************************************************************/
1293: /* Recover the direction. */
1294: /***************************************************************************/
1296: sigma = -1;
1298: /***************************************************************************/
1299: /* Start conjugate gradient method from the beginning */
1300: /***************************************************************************/
1302: VecCopy(ksp->vec_rhs, r); /* r = -grad */
1303: KSP_PCApply(ksp, r, z); /* z = inv(M) r */
1305: /***************************************************************************/
1306: /* Accumulate Q * s */
1307: /***************************************************************************/
1309: VecCopy(z, d);
1310: VecScale(d, sigma * t_soln[0] / cg->norm_r[0]);
1311:
1312: /***************************************************************************/
1313: /* Compute the first direction. */
1314: /***************************************************************************/
1316: VecCopy(z, p); /* p = z */
1317: KSP_MatMult(ksp, Qmat, p, z); /* z = Q * p */
1318: ++ksp->its;
1320: for (i = 0; i < l_size - 1; ++i) {
1321: /*************************************************************************/
1322: /* Update the residual and direction. */
1323: /*************************************************************************/
1325: alpha = cg->alpha[i];
1326: if (alpha >= 0.0) {
1327: sigma = -sigma;
1328: }
1330: VecAXPY(r, -alpha, z); /* r = r - alpha Q p */
1331: KSP_PCApply(ksp, r, z); /* z = inv(M) r */
1333: /*************************************************************************/
1334: /* Accumulate Q * s */
1335: /*************************************************************************/
1337: VecAXPY(d, sigma * t_soln[i+1] / cg->norm_r[i+1], z);
1339: /*************************************************************************/
1340: /* Update p. */
1341: /*************************************************************************/
1343: beta = cg->beta[i];
1344: VecAYPX(p, beta, z); /* p = z + beta p */
1345: KSP_MatMult(ksp, Qmat, p, z); /* z = Q * p */
1346: ++ksp->its;
1347: }
1349: /***************************************************************************/
1350: /* Update the residual and direction. */
1351: /***************************************************************************/
1353: alpha = cg->alpha[i];
1354: if (alpha >= 0.0) {
1355: sigma = -sigma;
1356: }
1358: VecAXPY(r, -alpha, z); /* r = r - alpha Q p */
1359: KSP_PCApply(ksp, r, z); /* z = inv(M) r */
1361: /***************************************************************************/
1362: /* Accumulate Q * s */
1363: /***************************************************************************/
1365: VecAXPY(d, sigma * t_soln[i+1] / cg->norm_r[i+1], z);
1367: /***************************************************************************/
1368: /* Set the termination reason. */
1369: /***************************************************************************/
1371: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
1372: return(0);
1373: #endif
1374: }
1378: PetscErrorCode KSPSetUp_GLTR(KSP ksp)
1379: {
1380: KSP_GLTR *cg = (KSP_GLTR *)ksp->data;
1381: PetscInt max_its;
1386: /***************************************************************************/
1387: /* This implementation of CG only handles left preconditioning so generate */
1388: /* an error otherwise. */
1389: /***************************************************************************/
1391: if (ksp->pc_side == PC_RIGHT) {
1392: SETERRQ(PETSC_ERR_SUP, "No right preconditioning for KSPGLTR");
1393: } else if (ksp->pc_side == PC_SYMMETRIC) {
1394: SETERRQ(PETSC_ERR_SUP, "No symmetric preconditioning for KSPGLTR");
1395: }
1397: /***************************************************************************/
1398: /* Determine the total maximum number of iterations. */
1399: /***************************************************************************/
1401: max_its = ksp->max_it + cg->max_lanczos_its + 1;
1403: /***************************************************************************/
1404: /* Set work vectors needed by conjugate gradient method and allocate */
1405: /* workspace for Lanczos matrix. */
1406: /***************************************************************************/
1408: KSPDefaultGetWork(ksp, 3);
1410: PetscMalloc5(max_its,PetscReal,&cg->diag,max_its,PetscReal,&cg->offd,max_its,PetscReal,&cg->alpha,max_its,PetscReal,&cg->beta,max_its,PetscReal,&cg->norm_r);
1411: PetscMemzero(cg->diag, max_its*sizeof(PetscReal));
1412: PetscMemzero(cg->offd, max_its*sizeof(PetscReal));
1413: PetscMemzero(cg->alpha, max_its*sizeof(PetscReal));
1414: PetscMemzero(cg->beta, max_its*sizeof(PetscReal));
1415: PetscMemzero(cg->norm_r, max_its*sizeof(PetscReal));
1416: PetscLogObjectMemory(ksp, 5*max_its*sizeof(PetscReal));
1417: return(0);
1418: }
1422: PetscErrorCode KSPDestroy_GLTR(KSP ksp)
1423: {
1424: KSP_GLTR *cg = (KSP_GLTR *)ksp->data;
1429: /***************************************************************************/
1430: /* Free memory allocated for the data. */
1431: /***************************************************************************/
1433: PetscFree5(cg->diag,cg->offd,cg->alpha,cg->beta,cg->norm_r);
1434: if (cg->alloced) {
1435: PetscFree2(cg->rwork,cg->iwork);
1436: }
1438: /***************************************************************************/
1439: /* Clear composed functions */
1440: /***************************************************************************/
1442: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGLTRSetRadius_C","",PETSC_NULL);
1443: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGLTRGetNormD_C","",PETSC_NULL);
1444: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGLTRGetObjFcn_C","",PETSC_NULL);
1445: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGLTRGetMinEig_C","",PETSC_NULL);
1446: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGLTRGetLambda_C","",PETSC_NULL);
1448: /***************************************************************************/
1449: /* Destroy KSP object. */
1450: /***************************************************************************/
1451: KSPDefaultDestroy(ksp);
1452: return(0);
1453: }
1458: PetscErrorCode KSPGLTRSetRadius_GLTR(KSP ksp, PetscReal radius)
1459: {
1460: KSP_GLTR *cg = (KSP_GLTR *)ksp->data;
1463: cg->radius = radius;
1464: return(0);
1465: }
1469: PetscErrorCode KSPGLTRGetNormD_GLTR(KSP ksp, PetscReal *norm_d)
1470: {
1471: KSP_GLTR *cg = (KSP_GLTR *)ksp->data;
1474: *norm_d = cg->norm_d;
1475: return(0);
1476: }
1480: PetscErrorCode KSPGLTRGetObjFcn_GLTR(KSP ksp, PetscReal *o_fcn)
1481: {
1482: KSP_GLTR *cg = (KSP_GLTR *)ksp->data;
1485: *o_fcn = cg->o_fcn;
1486: return(0);
1487: }
1491: PetscErrorCode KSPGLTRGetMinEig_GLTR(KSP ksp, PetscReal *e_min)
1492: {
1493: KSP_GLTR *cg = (KSP_GLTR *)ksp->data;
1496: *e_min = cg->e_min;
1497: return(0);
1498: }
1502: PetscErrorCode KSPGLTRGetLambda_GLTR(KSP ksp, PetscReal *lambda)
1503: {
1504: KSP_GLTR *cg = (KSP_GLTR *)ksp->data;
1507: *lambda = cg->lambda;
1508: return(0);
1509: }
1514: PetscErrorCode KSPSetFromOptions_GLTR(KSP ksp)
1515: {
1517: KSP_GLTR *cg = (KSP_GLTR *)ksp->data;
1520: PetscOptionsHead("KSP GLTR options");
1522: PetscOptionsReal("-ksp_gltr_radius", "Trust Region Radius", "KSPGLTRSetRadius", cg->radius, &cg->radius, PETSC_NULL);
1524: PetscOptionsReal("-ksp_gltr_init_pert", "Initial perturbation", "", cg->init_pert, &cg->init_pert, PETSC_NULL);
1525: PetscOptionsReal("-ksp_gltr_eigen_tol", "Eigenvalue tolerance", "", cg->eigen_tol, &cg->eigen_tol, PETSC_NULL);
1526: PetscOptionsReal("-ksp_gltr_newton_tol", "Newton tolerance", "", cg->newton_tol, &cg->newton_tol, PETSC_NULL);
1528: PetscOptionsInt("-ksp_gltr_max_lanczos_its", "Maximum Lanczos Iters", "", cg->max_lanczos_its, &cg->max_lanczos_its, PETSC_NULL);
1529: PetscOptionsInt("-ksp_gltr_max_newton_its", "Maximum Newton Iters", "", cg->max_newton_its, &cg->max_newton_its, PETSC_NULL);
1531: PetscOptionsEList("-ksp_gltr_dtype", "Norm used for direction", "", DType_Table, GLTR_DIRECTION_TYPES, DType_Table[cg->dtype], &cg->dtype, PETSC_NULL);
1533: PetscOptionsTail();
1534: return(0);
1535: }
1537: /*MC
1538: KSPGLTR - Code to run conjugate gradient method subject to a constraint
1539: on the solution norm. This is used in Trust Region methods for
1540: nonlinear equations, SNESTR
1542: Options Database Keys:
1543: . -ksp_gltr_radius <r> - Trust Region Radius
1545: Notes: This is rarely used directly
1547: Use preconditioned conjugate gradient to compute
1548: an approximate minimizer of the quadratic function
1550: q(s) = g^T * s + .5 * s^T * H * s
1552: subject to the trust region constraint
1554: || s || <= delta,
1556: where
1558: delta is the trust region radius,
1559: g is the gradient vector,
1560: H is the Hessian approximation,
1561: M is the positive definite preconditioner matrix.
1563: KSPConvergedReason may be
1564: $ KSP_CONVERGED_CG_NEG_CURVE if convergence is reached along a negative curvature direction,
1565: $ KSP_CONVERGED_CG_CONSTRAINED if convergence is reached along a constrained step,
1566: $ other KSP converged/diverged reasons
1568: Notes:
1569: The preconditioner supplied should be symmetric and positive definite.
1571: Level: developer
1573: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPGLTRSetRadius(), KSPGLTRGetNormD(), KSPGLTRGetObjFcn(), KSPGLTRGetMinEig(), KSPGLTRGetLambda()
1574: M*/
1579: PetscErrorCode KSPCreate_GLTR(KSP ksp)
1580: {
1582: KSP_GLTR *cg;
1585: PetscNewLog(ksp, KSP_GLTR, &cg);
1586: cg->radius = 0.0;
1587: cg->dtype = GLTR_UNPRECONDITIONED_DIRECTION;
1589: cg->init_pert = 1.0e-8;
1590: cg->eigen_tol = 1.0e-10;
1591: cg->newton_tol = 1.0e-6;
1593: cg->alloced = 0;
1594: cg->init_alloc = 1024;
1596: cg->max_lanczos_its = 20;
1597: cg->max_newton_its = 10;
1599: ksp->data = (void *) cg;
1600: ksp->pc_side = PC_LEFT;
1601: ksp->normtype = KSP_NORM_UNPRECONDITIONED;
1603: /***************************************************************************/
1604: /* Sets the functions that are associated with this data structure */
1605: /* (in C++ this is the same as defining virtual functions). */
1606: /***************************************************************************/
1608: ksp->ops->setup = KSPSetUp_GLTR;
1609: ksp->ops->solve = KSPSolve_GLTR;
1610: ksp->ops->destroy = KSPDestroy_GLTR;
1611: ksp->ops->setfromoptions = KSPSetFromOptions_GLTR;
1612: ksp->ops->buildsolution = KSPDefaultBuildSolution;
1613: ksp->ops->buildresidual = KSPDefaultBuildResidual;
1614: ksp->ops->view = 0;
1616: PetscObjectComposeFunctionDynamic((PetscObject)ksp,
1617: "KSPGLTRSetRadius_C",
1618: "KSPGLTRSetRadius_GLTR",
1619: KSPGLTRSetRadius_GLTR);
1620: PetscObjectComposeFunctionDynamic((PetscObject)ksp,
1621: "KSPGLTRGetNormD_C",
1622: "KSPGLTRGetNormD_GLTR",
1623: KSPGLTRGetNormD_GLTR);
1624: PetscObjectComposeFunctionDynamic((PetscObject)ksp,
1625: "KSPGLTRGetObjFcn_C",
1626: "KSPGLTRGetObjFcn_GLTR",
1627: KSPGLTRGetObjFcn_GLTR);
1628: PetscObjectComposeFunctionDynamic((PetscObject)ksp,
1629: "KSPGLTRGetMinEig_C",
1630: "KSPGLTRGetMinEig_GLTR",
1631: KSPGLTRGetMinEig_GLTR);
1632: PetscObjectComposeFunctionDynamic((PetscObject)ksp,
1633: "KSPGLTRGetLambda_C",
1634: "KSPGLTRGetLambda_GLTR",
1635: KSPGLTRGetLambda_GLTR);
1636: return(0);
1637: }