Actual source code: gltr.c
1: #define PETSCKSP_DLL
3: #include private/kspimpl.h
4: #include petscblaslapack.h
5: #include ../src/ksp/ksp/impls/cg/gltr/gltr.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] = {
12: "preconditioned", "unpreconditioned"
13: };
17: /*@
18: KSPGLTRSetRadius - Sets the radius of the trust region.
20: Collective on KSP
22: Input Parameters:
23: + ksp - the iterative context
24: - radius - the trust region radius (Infinity is the default)
26: Options Database Key:
27: . -ksp_gltr_radius <r>
29: Level: advanced
31: .keywords: KSP, GLTR, set, trust region radius
32: @*/
33: PetscErrorCode KSPGLTRSetRadius(KSP ksp, PetscReal radius)
34: {
35: PetscErrorCode ierr, (*f)(KSP, PetscReal);
39: if (radius < 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Radius negative");
40: PetscObjectQueryFunction((PetscObject)ksp, "KSPGLTRSetRadius_C", (void (**)(void))&f);
41: if (f) {
42: (*f)(ksp, radius);
43: }
44: return(0);
45: }
49: /*@
50: KSPGLTRGetNormD - Get norm of the direction.
52: Collective on KSP
54: Input Parameters:
55: + ksp - the iterative context
56: - norm_d - the norm of the direction
58: Level: advanced
60: .keywords: KSP, GLTR, get, norm direction
61: @*/
62: PetscErrorCode KSPGLTRGetNormD(KSP ksp, PetscReal *norm_d)
63: {
64: PetscErrorCode ierr, (*f)(KSP, PetscReal *);
68: PetscObjectQueryFunction((PetscObject)ksp, "KSPGLTRGetNormD_C", (void (**)(void))&f);
69: if (f) {
70: (*f)(ksp, norm_d);
71: }
72: return(0);
73: }
77: /*@
78: KSPGLTRGetObjFcn - Get objective function value.
80: Collective on KSP
82: Input Parameters:
83: + ksp - the iterative context
84: - o_fcn - the objective function value
86: Level: advanced
88: .keywords: KSP, GLTR, get, objective function
89: @*/
90: PetscErrorCode KSPGLTRGetObjFcn(KSP ksp, PetscReal *o_fcn)
91: {
92: PetscErrorCode ierr, (*f)(KSP, PetscReal *);
96: PetscObjectQueryFunction((PetscObject)ksp, "KSPGLTRGetObjFcn_C", (void (**)(void))&f);
97: if (f) {
98: (*f)(ksp, o_fcn);
99: }
100: return(0);
101: }
105: /*@
106: KSPGLTRGetMinEig - Get minimum eigenvalue.
108: Collective on KSP
110: Input Parameters:
111: + ksp - the iterative context
112: - e_min - the minimum eigenvalue
114: Level: advanced
116: .keywords: KSP, GLTR, get, minimum eigenvalue
117: @*/
118: PetscErrorCode KSPGLTRGetMinEig(KSP ksp, PetscReal *e_min)
119: {
120: PetscErrorCode ierr, (*f)(KSP, PetscReal *);
124: PetscObjectQueryFunction((PetscObject)ksp, "KSPGLTRGetMinEig_C", (void (**)(void))&f);
125: if (f) {
126: (*f)(ksp, e_min);
127: }
128: return(0);
129: }
133: /*@
134: KSPGLTRGetLambda - Get multiplier on trust-region constraint.
136: Collective on KSP
138: Input Parameters:
139: + ksp - the iterative context
140: - lambda - the multiplier
142: Level: advanced
144: .keywords: KSP, GLTR, get, multiplier
145: @*/
146: PetscErrorCode KSPGLTRGetLambda(KSP ksp, PetscReal *lambda)
147: {
148: PetscErrorCode ierr, (*f)(KSP, PetscReal *);
152: PetscObjectQueryFunction((PetscObject)ksp, "KSPGLTRGetLambda_C", (void (**)(void))&f);
153: if (f) {
154: (*f)(ksp, lambda);
155: }
156: return(0);
157: }
161: /*
162: KSPSolve_GLTR - Use preconditioned conjugate gradient to compute
163: an approximate minimizer of the quadratic function
165: q(s) = g^T * s + .5 * s^T * H * s
167: subject to the trust region constraint
169: || s || <= delta,
171: where
173: delta is the trust region radius,
174: g is the gradient vector,
175: H is the Hessian approximation,
176: M is the positive definite preconditioner matrix.
178: KSPConvergedReason may be
179: $ KSP_CONVERGED_CG_NEG_CURVE if convergence is reached along a negative curvature direction,
180: $ KSP_CONVERGED_CG_CONSTRAINED if convergence is reached along a constrained step,
181: $ other KSP converged/diverged reasons
183: Notes:
184: The preconditioner supplied should be symmetric and positive definite.
185: */
186: PetscErrorCode KSPSolve_GLTR(KSP ksp)
187: {
188: #ifdef PETSC_USE_COMPLEX
189: SETERRQ(PETSC_ERR_SUP, "GLTR is not available for complex systems");
190: #else
191: KSP_GLTR *cg = (KSP_GLTR *)ksp->data;
192: PetscReal *t_soln, *t_diag, *t_offd, *e_valu, *e_vect, *e_rwrk;
193: PetscBLASInt *e_iblk, *e_splt, *e_iwrk;
196: MatStructure pflag;
197: Mat Qmat, Mmat;
198: Vec r, z, p, d;
199: PC pc;
201: PetscReal norm_r, norm_d, norm_dp1, norm_p, dMp;
202: PetscReal alpha, beta, kappa, rz, rzm1;
203: PetscReal rr, r2, piv, step;
204: PetscReal vl, vu;
205: PetscReal coef1, coef2, coef3, root1, root2, obj1, obj2;
206: PetscReal norm_t, norm_w, pert;
208: PetscInt i, j, max_cg_its, max_lanczos_its, max_newton_its, sigma;
209: PetscBLASInt t_size = 0, l_size = 0, il, iu, e_valus, info;
210: PetscBLASInt nrhs, nldb;
212: #if !defined(PETSC_MISSING_LAPACK_STEBZ)
213: PetscBLASInt e_splts;
214: #endif
216: PetscTruth diagonalscale;
220: /***************************************************************************/
221: /* Check the arguments and parameters. */
222: /***************************************************************************/
224: PCDiagonalScale(ksp->pc, &diagonalscale);
225: if (diagonalscale) {
226: SETERRQ1(PETSC_ERR_SUP, "Krylov method %s does not support diagonal scaling", ((PetscObject)ksp)->type_name);
227: }
229: if (cg->radius < 0.0) {
230: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Input error: radius < 0");
231: }
233: /***************************************************************************/
234: /* Get the workspace vectors and initialize variables */
235: /***************************************************************************/
237: r2 = cg->radius * cg->radius;
238: r = ksp->work[0];
239: z = ksp->work[1];
240: p = ksp->work[2];
241: d = ksp->vec_sol;
242: pc = ksp->pc;
244: PCGetOperators(pc, &Qmat, &Mmat, &pflag);
246: VecGetSize(d, &max_cg_its);
247: max_cg_its = PetscMin(max_cg_its, ksp->max_it);
248: max_lanczos_its = cg->max_lanczos_its;
249: max_newton_its = cg->max_newton_its;
250: ksp->its = 0;
252: /***************************************************************************/
253: /* Initialize objective function direction, and minimum eigenvalue. */
254: /***************************************************************************/
256: cg->o_fcn = 0.0;
258: VecSet(d, 0.0); /* d = 0 */
259: cg->norm_d = 0.0;
261: cg->e_min = 0.0;
262: cg->lambda = 0.0;
264: /***************************************************************************/
265: /* The first phase of GLTR performs a standard conjugate gradient method, */
266: /* but stores the values required for the Lanczos matrix. We switch to */
267: /* the Lanczos when the conjugate gradient method breaks down. Check the */
268: /* right-hand side for numerical problems. The check for not-a-number and */
269: /* infinite values need be performed only once. */
270: /***************************************************************************/
272: VecCopy(ksp->vec_rhs, r); /* r = -grad */
273: VecDot(r, r, &rr); /* rr = r^T r */
274: if PetscIsInfOrNanScalar(rr) {
275: /*************************************************************************/
276: /* The right-hand side contains not-a-number or an infinite value. */
277: /* The gradient step does not work; return a zero value for the step. */
278: /*************************************************************************/
280: ksp->reason = KSP_DIVERGED_NAN;
281: PetscInfo1(ksp, "KSPSolve_GLTR: bad right-hand side: rr=%g\n", rr);
282: return(0);
283: }
285: /***************************************************************************/
286: /* Check the preconditioner for numerical problems and for positive */
287: /* definiteness. The check for not-a-number and infinite values need be */
288: /* performed only once. */
289: /***************************************************************************/
291: KSP_PCApply(ksp, r, z); /* z = inv(M) r */
292: VecDot(r, z, &rz); /* rz = r^T inv(M) r */
293: if PetscIsInfOrNanScalar(rz) {
294: /*************************************************************************/
295: /* The preconditioner contains not-a-number or an infinite value. */
296: /* Return the gradient direction intersected with the trust region. */
297: /*************************************************************************/
299: ksp->reason = KSP_DIVERGED_NAN;
300: PetscInfo1(ksp, "KSPSolve_GLTR: bad preconditioner: rz=%g\n", rz);
302: if (cg->radius) {
303: if (r2 >= rr) {
304: alpha = 1.0;
305: cg->norm_d = sqrt(rr);
306: }
307: else {
308: alpha = sqrt(r2 / rr);
309: cg->norm_d = cg->radius;
310: }
312: VecAXPY(d, alpha, r); /* d = d + alpha r */
313:
314: /***********************************************************************/
315: /* Compute objective function. */
316: /***********************************************************************/
318: KSP_MatMult(ksp, Qmat, d, z); CHKERRQ(ierr)
319: VecAYPX(z, -0.5, ksp->vec_rhs);
320: VecDot(d, z, &cg->o_fcn);
321: cg->o_fcn = -cg->o_fcn;
322: ++ksp->its;
323: }
324: return(0);
325: }
327: if (rz < 0.0) {
328: /*************************************************************************/
329: /* The preconditioner is indefinite. Because this is the first */
330: /* and we do not have a direction yet, we use the gradient step. Note */
331: /* that we cannot use the preconditioned norm when computing the step */
332: /* because the matrix is indefinite. */
333: /*************************************************************************/
335: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
336: PetscInfo1(ksp, "KSPSolve_GLTR: indefinite preconditioner: rz=%g\n", rz);
338: if (cg->radius) {
339: if (r2 >= rr) {
340: alpha = 1.0;
341: cg->norm_d = sqrt(rr);
342: }
343: else {
344: alpha = sqrt(r2 / rr);
345: cg->norm_d = cg->radius;
346: }
348: VecAXPY(d, alpha, r); /* d = d + alpha r */
350: /***********************************************************************/
351: /* Compute objective function. */
352: /***********************************************************************/
354: KSP_MatMult(ksp, Qmat, d, z); CHKERRQ(ierr)
355: VecAYPX(z, -0.5, ksp->vec_rhs);
356: VecDot(d, z, &cg->o_fcn);
357: cg->o_fcn = -cg->o_fcn;
358: ++ksp->its;
359: }
360: return(0);
361: }
363: /***************************************************************************/
364: /* As far as we know, the preconditioner is positive semidefinite. */
365: /* Compute and log the residual. Check convergence because this */
366: /* initializes things, but do not terminate until at least one conjugate */
367: /* gradient iteration has been performed. */
368: /***************************************************************************/
370: cg->norm_r[0] = sqrt(rz); /* norm_r = |r|_M */
372: switch(ksp->normtype) {
373: case KSP_NORM_PRECONDITIONED:
374: VecNorm(z, NORM_2, &norm_r); /* norm_r = |z| */
375: break;
377: case KSP_NORM_UNPRECONDITIONED:
378: norm_r = sqrt(rr); /* norm_r = |r| */
379: break;
381: case KSP_NORM_NATURAL:
382: norm_r = cg->norm_r[0]; /* norm_r = |r|_M */
383: break;
385: default:
386: norm_r = 0.0;
387: break;
388: }
390: KSPLogResidualHistory(ksp, norm_r);
391: KSPMonitor(ksp, ksp->its, norm_r);
392: ksp->rnorm = norm_r;
394: (*ksp->converged)(ksp, ksp->its, norm_r, &ksp->reason, ksp->cnvP);
396: /***************************************************************************/
397: /* Compute the first direction and update the iteration. */
398: /***************************************************************************/
400: VecCopy(z, p); /* p = z */
401: KSP_MatMult(ksp, Qmat, p, z); /* z = Q * p */
402: ++ksp->its;
404: /***************************************************************************/
405: /* Check the matrix for numerical problems. */
406: /***************************************************************************/
408: VecDot(p, z, &kappa); /* kappa = p^T Q p */
409: if PetscIsInfOrNanScalar(kappa) {
410: /*************************************************************************/
411: /* The matrix produced not-a-number or an infinite value. In this case, */
412: /* we must stop and use the gradient direction. This condition need */
413: /* only be checked once. */
414: /*************************************************************************/
416: ksp->reason = KSP_DIVERGED_NAN;
417: PetscInfo1(ksp, "KSPSolve_GLTR: bad matrix: kappa=%g\n", kappa);
419: if (cg->radius) {
420: if (r2 >= rr) {
421: alpha = 1.0;
422: cg->norm_d = sqrt(rr);
423: }
424: else {
425: alpha = sqrt(r2 / rr);
426: cg->norm_d = cg->radius;
427: }
429: VecAXPY(d, alpha, r); /* d = d + alpha r */
431: /***********************************************************************/
432: /* Compute objective function. */
433: /***********************************************************************/
435: KSP_MatMult(ksp, Qmat, d, z); CHKERRQ(ierr)
436: VecAYPX(z, -0.5, ksp->vec_rhs);
437: VecDot(d, z, &cg->o_fcn);
438: cg->o_fcn = -cg->o_fcn;
439: ++ksp->its;
440: }
441: return(0);
442: }
444: /***************************************************************************/
445: /* Initialize variables for calculating the norm of the direction and for */
446: /* the Lanczos tridiagonal matrix. Note that we track the diagonal value */
447: /* of the Cholesky factorization of the Lanczos matrix in order to */
448: /* determine when negative curvature is encountered. */
449: /***************************************************************************/
451: dMp = 0.0;
452: norm_d = 0.0;
453: switch(cg->dtype) {
454: case GLTR_PRECONDITIONED_DIRECTION:
455: norm_p = rz;
456: break;
458: default:
459: VecDot(p, p, &norm_p);
460: break;
461: }
463: cg->diag[t_size] = kappa / rz;
464: cg->offd[t_size] = 0.0;
465: ++t_size;
467: piv = 1.0;
469: /***************************************************************************/
470: /* Check for breakdown of the conjugate gradient method; this occurs when */
471: /* kappa is zero. */
472: /***************************************************************************/
474: if (fabs(kappa) <= 0.0) {
475: /*************************************************************************/
476: /* The curvature is zero. In this case, we must stop and use follow */
477: /* the direction of negative curvature since the Lanczos matrix is zero. */
478: /*************************************************************************/
480: ksp->reason = KSP_DIVERGED_BREAKDOWN;
481: PetscInfo1(ksp, "KSPSolve_GLTR: breakdown: kappa=%g\n", kappa);
483: if (cg->radius && norm_p > 0.0) {
484: /***********************************************************************/
485: /* Follow direction of negative curvature to the boundary of the */
486: /* trust region. */
487: /***********************************************************************/
489: step = sqrt(r2 / norm_p);
490: cg->norm_d = cg->radius;
492: VecAXPY(d, step, p); /* d = d + step p */
494: /***********************************************************************/
495: /* Update objective function. */
496: /***********************************************************************/
498: cg->o_fcn += step * (0.5 * step * kappa - rz);
499: }
500: else if (cg->radius) {
501: /***********************************************************************/
502: /* The norm of the preconditioned direction is zero; use the gradient */
503: /* step. */
504: /***********************************************************************/
506: if (r2 >= rr) {
507: alpha = 1.0;
508: cg->norm_d = sqrt(rr);
509: }
510: else {
511: alpha = sqrt(r2 / rr);
512: cg->norm_d = cg->radius;
513: }
515: VecAXPY(d, alpha, r); /* d = d + alpha r */
517: /***********************************************************************/
518: /* Compute objective function. */
519: /***********************************************************************/
521: KSP_MatMult(ksp, Qmat, d, z); CHKERRQ(ierr)
522: VecAYPX(z, -0.5, ksp->vec_rhs);
523: VecDot(d, z, &cg->o_fcn);
524: cg->o_fcn = -cg->o_fcn;
525: ++ksp->its;
526: }
527: return(0);
528: }
530: /***************************************************************************/
531: /* Begin the first part of the GLTR algorithm which runs the conjugate */
532: /* gradient method until either the problem is solved, we encounter the */
533: /* boundary of the trust region, or the conjugate gradient method breaks */
534: /* down. */
535: /***************************************************************************/
537: while(1) {
538: /*************************************************************************/
539: /* Know that kappa is nonzero, because we have not broken down, so we */
540: /* can compute the steplength. */
541: /*************************************************************************/
543: alpha = rz / kappa;
544: cg->alpha[l_size] = alpha;
546: /*************************************************************************/
547: /* Compute the diagonal value of the Cholesky factorization of the */
548: /* Lanczos matrix and check to see if the Lanczos matrix is indefinite. */
549: /* This indicates a direction of negative curvature. */
550: /*************************************************************************/
552: piv = cg->diag[l_size] - cg->offd[l_size]*cg->offd[l_size] / piv;
553: if (piv <= 0.0) {
554: /***********************************************************************/
555: /* In this case, the matrix is indefinite and we have encountered */
556: /* a direction of negative curvature. Follow the direction to the */
557: /* boundary of the trust region. */
558: /***********************************************************************/
560: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
561: PetscInfo1(ksp, "KSPSolve_GLTR: negative curvature: kappa=%g\n", kappa);
563: if (cg->radius && norm_p > 0.0) {
564: /*********************************************************************/
565: /* Follow direction of negative curvature to boundary. */
566: /*********************************************************************/
568: step = (sqrt(dMp*dMp+norm_p*(r2-norm_d))-dMp)/norm_p;
569: cg->norm_d = cg->radius;
571: VecAXPY(d, step, p); /* d = d + step p */
573: /*********************************************************************/
574: /* Update objective function. */
575: /*********************************************************************/
577: cg->o_fcn += step * (0.5 * step * kappa - rz);
578: }
579: else if (cg->radius) {
580: /*********************************************************************/
581: /* The norm of the direction is zero; there is nothing to follow. */
582: /*********************************************************************/
583: }
584: break;
585: }
587: /*************************************************************************/
588: /* Compute the steplength and check for intersection with the trust */
589: /* region. */
590: /*************************************************************************/
592: norm_dp1 = norm_d + alpha*(2.0*dMp + alpha*norm_p);
593: if (cg->radius && norm_dp1 >= r2) {
594: /***********************************************************************/
595: /* In this case, the matrix is positive definite as far as we know. */
596: /* However, the full step goes beyond the trust region. */
597: /***********************************************************************/
599: ksp->reason = KSP_CONVERGED_CG_CONSTRAINED;
600: PetscInfo1(ksp, "KSPSolve_GLTR: constrained step: radius=%g\n", cg->radius);
602: if (norm_p > 0.0) {
603: /*********************************************************************/
604: /* Follow the direction to the boundary of the trust region. */
605: /*********************************************************************/
607: step = (sqrt(dMp*dMp+norm_p*(r2-norm_d))-dMp)/norm_p;
608: cg->norm_d = cg->radius;
610: VecAXPY(d, step, p); /* d = d + step p */
612: /*********************************************************************/
613: /* Update objective function. */
614: /*********************************************************************/
616: cg->o_fcn += step * (0.5 * step * kappa - rz);
617: }
618: else {
619: /*********************************************************************/
620: /* The norm of the direction is zero; there is nothing to follow. */
621: /*********************************************************************/
622: }
623: break;
624: }
626: /*************************************************************************/
627: /* Now we can update the direction and residual. */
628: /*************************************************************************/
630: VecAXPY(d, alpha, p); /* d = d + alpha p */
631: VecAXPY(r, -alpha, z); /* r = r - alpha Q p */
632: KSP_PCApply(ksp, r, z); /* z = inv(M) r */
634: switch(cg->dtype) {
635: case GLTR_PRECONDITIONED_DIRECTION:
636: norm_d = norm_dp1;
637: break;
639: default:
640: VecDot(d, d, &norm_d);
641: break;
642: }
643: cg->norm_d = sqrt(norm_d);
645: /*************************************************************************/
646: /* Update objective function. */
647: /*************************************************************************/
649: cg->o_fcn -= 0.5 * alpha * rz;
651: /*************************************************************************/
652: /* Check that the preconditioner appears positive semidefinite. */
653: /*************************************************************************/
655: rzm1 = rz;
656: VecDot(r, z, &rz); /* rz = r^T z */
657: if (rz < 0.0) {
658: /***********************************************************************/
659: /* The preconditioner is indefinite. */
660: /***********************************************************************/
662: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
663: PetscInfo1(ksp, "KSPSolve_GLTR: cg indefinite preconditioner: rz=%g\n", rz);
664: break;
665: }
667: /*************************************************************************/
668: /* As far as we know, the preconditioner is positive semidefinite. */
669: /* Compute the residual and check for convergence. */
670: /*************************************************************************/
672: cg->norm_r[l_size+1] = sqrt(rz); /* norm_r = |r|_M */
674: switch(ksp->normtype) {
675: case KSP_NORM_PRECONDITIONED:
676: VecNorm(z, NORM_2, &norm_r); /* norm_r = |z| */
677: break;
679: case KSP_NORM_UNPRECONDITIONED:
680: VecNorm(r, NORM_2, &norm_r); /* norm_r = |r| */
681: break;
683: case KSP_NORM_NATURAL:
684: norm_r = cg->norm_r[l_size+1]; /* norm_r = |r|_M */
685: break;
687: default:
688: norm_r = 0.0;
689: break;
690: }
692: KSPLogResidualHistory(ksp, norm_r);
693: KSPMonitor(ksp, ksp->its, norm_r);
694: ksp->rnorm = norm_r;
695:
696: (*ksp->converged)(ksp, ksp->its, norm_r, &ksp->reason, ksp->cnvP);
697: if (ksp->reason) {
698: /***********************************************************************/
699: /* The method has converged. */
700: /***********************************************************************/
702: PetscInfo2(ksp, "KSPSolve_GLTR: cg truncated step: rnorm=%g, radius=%g\n", norm_r, cg->radius);
703: break;
704: }
706: /*************************************************************************/
707: /* We have not converged yet. Check for breakdown. */
708: /*************************************************************************/
710: beta = rz / rzm1;
711: if (fabs(beta) <= 0.0) {
712: /***********************************************************************/
713: /* Conjugate gradients has broken down. */
714: /***********************************************************************/
716: ksp->reason = KSP_DIVERGED_BREAKDOWN;
717: PetscInfo1(ksp, "KSPSolve_GLTR: breakdown: beta=%g\n", beta);
718: break;
719: }
721: /*************************************************************************/
722: /* Check iteration limit. */
723: /*************************************************************************/
725: if (ksp->its >= max_cg_its) {
726: ksp->reason = KSP_DIVERGED_ITS;
727: PetscInfo1(ksp, "KSPSolve_GLTR: iterlim: its=%d\n", ksp->its);
728: break;
729: }
731: /*************************************************************************/
732: /* Update p and the norms. */
733: /*************************************************************************/
735: cg->beta[l_size] = beta;
736: VecAYPX(p, beta, z); /* p = z + beta p */
738: switch(cg->dtype) {
739: case GLTR_PRECONDITIONED_DIRECTION:
740: dMp = beta*(dMp + alpha*norm_p);
741: norm_p = beta*(rzm1 + beta*norm_p);
742: break;
744: default:
745: VecDot(d, p, &dMp);
746: VecDot(p, p, &norm_p);
747: break;
748: }
750: /*************************************************************************/
751: /* Compute the new direction and update the iteration. */
752: /*************************************************************************/
754: KSP_MatMult(ksp, Qmat, p, z); /* z = Q * p */
755: VecDot(p, z, &kappa); /* kappa = p^T Q p */
756: ++ksp->its;
758: /*************************************************************************/
759: /* Update the the Lanczos tridiagonal matrix. */
760: /*************************************************************************/
762: ++l_size;
763: cg->offd[t_size] = sqrt(beta) / fabs(alpha);
764: cg->diag[t_size] = kappa / rz + beta / alpha;
765: ++t_size;
767: /*************************************************************************/
768: /* Check for breakdown of the conjugate gradient method; this occurs */
769: /* when kappa is zero. */
770: /*************************************************************************/
772: if (fabs(kappa) <= 0.0) {
773: /***********************************************************************/
774: /* The method breaks down; move along the direction as if the matrix */
775: /* were indefinite. */
776: /***********************************************************************/
778: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
779: PetscInfo1(ksp, "KSPSolve_GLTR: cg breakdown: kappa=%g\n", kappa);
781: if (cg->radius && norm_p > 0.0) {
782: /*********************************************************************/
783: /* Follow direction to boundary. */
784: /*********************************************************************/
786: step = (sqrt(dMp*dMp+norm_p*(r2-norm_d))-dMp)/norm_p;
787: cg->norm_d = cg->radius;
789: VecAXPY(d, step, p); /* d = d + step p */
791: /*********************************************************************/
792: /* Update objective function. */
793: /*********************************************************************/
795: cg->o_fcn += step * (0.5 * step * kappa - rz);
796: }
797: else if (cg->radius) {
798: /*********************************************************************/
799: /* The norm of the direction is zero; there is nothing to follow. */
800: /*********************************************************************/
801: }
802: break;
803: }
804: }
806: /***************************************************************************/
807: /* Check to see if we need to continue with the Lanczos method. */
808: /***************************************************************************/
810: if (!cg->radius) {
811: /*************************************************************************/
812: /* There is no radius. Therefore, we cannot move along the boundary. */
813: /*************************************************************************/
815: return(0);
816: }
818: if (KSP_CONVERGED_CG_NEG_CURVE != ksp->reason) {
819: /*************************************************************************/
820: /* The method either converged to an interior point, hit the boundary of */
821: /* the trust-region without encountering a direction of negative */
822: /* curvature or the iteration limit was reached. */
823: /*************************************************************************/
825: return(0);
826: }
828: /***************************************************************************/
829: /* Switch to contructing the Lanczos basis by way of the conjugate */
830: /* directions. */
831: /***************************************************************************/
833: for (i = 0; i < max_lanczos_its; ++i) {
834: /*************************************************************************/
835: /* Check for breakdown of the conjugate gradient method; this occurs */
836: /* when kappa is zero. */
837: /*************************************************************************/
839: if (fabs(kappa) <= 0.0) {
840: ksp->reason = KSP_DIVERGED_BREAKDOWN;
841: PetscInfo1(ksp, "KSPSolve_GLTR: lanczos breakdown: kappa=%g\n", kappa);
842: break;
843: }
845: /*************************************************************************/
846: /* Update the direction and residual. */
847: /*************************************************************************/
848:
849: alpha = rz / kappa;
850: cg->alpha[l_size] = alpha;
852: VecAXPY(r, -alpha, z); /* r = r - alpha Q p */
853: KSP_PCApply(ksp, r, z); /* z = inv(M) r */
855: /*************************************************************************/
856: /* Check that the preconditioner appears positive semidefinite. */
857: /*************************************************************************/
859: rzm1 = rz;
860: VecDot(r, z, &rz); /* rz = r^T z */
861: if (rz < 0.0) {
862: /***********************************************************************/
863: /* The preconditioner is indefinite. */
864: /***********************************************************************/
866: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
867: PetscInfo1(ksp, "KSPSolve_GLTR: lanczos indefinite preconditioner: rz=%g\n", rz);
868: break;
869: }
871: /*************************************************************************/
872: /* As far as we know, the preconditioner is positive definite. Compute */
873: /* the residual. Do NOT check for convergence. */
874: /*************************************************************************/
876: cg->norm_r[l_size+1] = sqrt(rz); /* norm_r = |r|_M */
878: switch(ksp->normtype) {
879: case KSP_NORM_PRECONDITIONED:
880: VecNorm(z, NORM_2, &norm_r); /* norm_r = |z| */
881: break;
883: case KSP_NORM_UNPRECONDITIONED:
884: VecNorm(r, NORM_2, &norm_r); /* norm_r = |r| */
885: break;
887: case KSP_NORM_NATURAL:
888: norm_r = cg->norm_r[l_size+1]; /* norm_r = |r|_M */
889: break;
891: default:
892: norm_r = 0.0;
893: break;
894: }
896: KSPLogResidualHistory(ksp, norm_r);
897: KSPMonitor(ksp, ksp->its, norm_r);
898: ksp->rnorm = norm_r;
899:
900: /*************************************************************************/
901: /* Check for breakdown. */
902: /*************************************************************************/
904: beta = rz / rzm1;
905: if (fabs(beta) <= 0.0) {
906: /***********************************************************************/
907: /* Conjugate gradients has broken down. */
908: /***********************************************************************/
910: ksp->reason = KSP_DIVERGED_BREAKDOWN;
911: PetscInfo1(ksp, "KSPSolve_GLTR: breakdown: beta=%g\n", beta);
912: break;
913: }
915: /*************************************************************************/
916: /* Update p and the norms. */
917: /*************************************************************************/
919: cg->beta[l_size] = beta;
920: VecAYPX(p, beta, z); /* p = z + beta p */
922: /*************************************************************************/
923: /* Compute the new direction and update the iteration. */
924: /*************************************************************************/
926: KSP_MatMult(ksp, Qmat, p, z); /* z = Q * p */
927: VecDot(p, z, &kappa); /* kappa = p^T Q p */
928: ++ksp->its;
930: /*************************************************************************/
931: /* Update the Lanczos tridiagonal matrix. */
932: /*************************************************************************/
934: ++l_size;
935: cg->offd[t_size] = sqrt(beta) / fabs(alpha);
936: cg->diag[t_size] = kappa / rz + beta / alpha;
937: ++t_size;
938: }
940: /***************************************************************************/
941: /* We have the Lanczos basis, solve the tridiagonal trust-region problem */
942: /* to obtain the Lanczos direction. We know that the solution lies on */
943: /* the boundary of the trust region. We start by checking that the */
944: /* workspace allocated is large enough. */
945: /***************************************************************************/
946: /* Note that the current version only computes the solution by using the */
947: /* preconditioned direction. Need to think about how to do the */
948: /* unpreconditioned direction calculation. */
949: /***************************************************************************/
951: if (t_size > cg->alloced) {
952: if (cg->alloced) {
953: PetscFree(cg->rwork);
954: PetscFree(cg->iwork);
955: cg->alloced += cg->init_alloc;
956: }
957: else {
958: cg->alloced = cg->init_alloc;
959: }
961: while (t_size > cg->alloced) {
962: cg->alloced += cg->init_alloc;
963: }
964:
965: cg->alloced = PetscMin(cg->alloced, t_size);
966: PetscMalloc(10*sizeof(PetscReal)*cg->alloced, &cg->rwork);
967: PetscMalloc(5*sizeof(PetscBLASInt)*cg->alloced, &cg->iwork);
968: }
970: /***************************************************************************/
971: /* Set up the required vectors. */
972: /***************************************************************************/
974: t_soln = cg->rwork + 0*t_size; /* Solution */
975: t_diag = cg->rwork + 1*t_size; /* Diagonal of T */
976: t_offd = cg->rwork + 2*t_size; /* Off-diagonal of T */
977: e_valu = cg->rwork + 3*t_size; /* Eigenvalues of T */
978: e_vect = cg->rwork + 4*t_size; /* Eigenvector of T */
979: e_rwrk = cg->rwork + 5*t_size; /* Eigen workspace */
980:
981: e_iblk = cg->iwork + 0*t_size; /* Eigen blocks */
982: e_splt = cg->iwork + 1*t_size; /* Eigen splits */
983: e_iwrk = cg->iwork + 2*t_size; /* Eigen workspace */
985: /***************************************************************************/
986: /* Compute the minimum eigenvalue of T. */
987: /***************************************************************************/
989: vl = 0.0;
990: vu = 0.0;
991: il = 1;
992: iu = 1;
994: #if defined(PETSC_MISSING_LAPACK_STEBZ)
995: SETERRQ(PETSC_ERR_SUP,"STEBZ - Lapack routine is unavailable.");
996: #else
997: LAPACKstebz_("I", "E", &t_size, &vl, &vu, &il, &iu, &cg->eigen_tol,
998: cg->diag, cg->offd + 1, &e_valus, &e_splts, e_valu,
999: e_iblk, e_splt, e_rwrk, e_iwrk, &info);
1001: if ((0 != info) || (1 != e_valus)) {
1002: /*************************************************************************/
1003: /* Calculation of the minimum eigenvalue failed. Return the */
1004: /* Steihaug-Toint direction. */
1005: /*************************************************************************/
1007: PetscInfo(ksp, "KSPSolve_GLTR: failed to compute eigenvalue.\n");
1008: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
1009: return(0);
1010: }
1012: cg->e_min = e_valu[0];
1014: /***************************************************************************/
1015: /* Compute the initial value of lambda to make (T + lamba I) positive */
1016: /* definite. */
1017: /***************************************************************************/
1019: pert = cg->init_pert;
1020: if (e_valu[0] < 0.0) {
1021: cg->lambda = pert - e_valu[0];
1022: }
1023: #endif
1025: while(1) {
1026: for (i = 0; i < t_size; ++i) {
1027: t_diag[i] = cg->diag[i] + cg->lambda;
1028: t_offd[i] = cg->offd[i];
1029: }
1031: #if defined(PETSC_MISSING_LAPACK_PTTRF)
1032: SETERRQ(PETSC_ERR_SUP,"PTTRF - Lapack routine is unavailable.");
1033: #else
1034: LAPACKpttrf_(&t_size, t_diag, t_offd + 1, &info);
1036: if (0 == info) {
1037: break;
1038: }
1040: pert += pert;
1041: cg->lambda = cg->lambda * (1.0 + pert) + pert;
1042: #endif
1043: }
1045: /***************************************************************************/
1046: /* Compute the initial step and its norm. */
1047: /***************************************************************************/
1049: nrhs = 1;
1050: nldb = t_size;
1052: t_soln[0] = -cg->norm_r[0];
1053: for (i = 1; i < t_size; ++i) {
1054: t_soln[i] = 0.0;
1055: }
1057: #if defined(PETSC_MISSING_LAPACK_PTTRS)
1058: SETERRQ(PETSC_ERR_SUP,"PTTRS - Lapack routine is unavailable.");
1059: #else
1060: LAPACKpttrs_(&t_size, &nrhs, t_diag, t_offd + 1, t_soln, &nldb, &info);
1061: #endif
1063: if (0 != info) {
1064: /*************************************************************************/
1065: /* Calculation of the initial step failed; return the Steihaug-Toint */
1066: /* direction. */
1067: /*************************************************************************/
1069: PetscInfo(ksp, "KSPSolve_GLTR: failed to compute step.\n");
1070: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
1071: return(0);
1072: }
1074: norm_t = 0;
1075: for (i = 0; i < t_size; ++i) {
1076: norm_t += t_soln[i] * t_soln[i];
1077: }
1078: norm_t = sqrt(norm_t);
1080: /***************************************************************************/
1081: /* Determine the case we are in. */
1082: /***************************************************************************/
1084: if (norm_t <= cg->radius) {
1085: /*************************************************************************/
1086: /* The step is within the trust region; check if we are in the hard case */
1087: /* and need to move to the boundary by following a direction of negative */
1088: /* curvature. */
1089: /*************************************************************************/
1090:
1091: if ((e_valu[0] <= 0.0) && (norm_t < cg->radius)) {
1092: /***********************************************************************/
1093: /* This is the hard case; compute the eigenvector associated with the */
1094: /* minimum eigenvalue and move along this direction to the boundary. */
1095: /***********************************************************************/
1097: #if defined(PETSC_MISSING_LAPACK_STEIN)
1098: SETERRQ(PETSC_ERR_SUP,"STEIN - Lapack routine is unavailable.");
1099: #else
1100: LAPACKstein_(&t_size, cg->diag, cg->offd + 1, &e_valus, e_valu,
1101: e_iblk, e_splt, e_vect, &nldb,
1102: e_rwrk, e_iwrk, e_iwrk + t_size, &info);
1103: #endif
1105: if (0 != info) {
1106: /*********************************************************************/
1107: /* Calculation of the minimum eigenvalue failed. Return the */
1108: /* Steihaug-Toint direction. */
1109: /*********************************************************************/
1110:
1111: PetscInfo(ksp, "KSPSolve_GLTR: failed to compute eigenvector.\n");
1112: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
1113: return(0);
1114: }
1115:
1116: coef1 = 0.0;
1117: coef2 = 0.0;
1118: coef3 = -cg->radius * cg->radius;
1119: for (i = 0; i < t_size; ++i) {
1120: coef1 += e_vect[i] * e_vect[i];
1121: coef2 += e_vect[i] * t_soln[i];
1122: coef3 += t_soln[i] * t_soln[i];
1123: }
1124:
1125: coef3 = sqrt(coef2 * coef2 - coef1 * coef3);
1126: root1 = (-coef2 + coef3) / coef1;
1127: root2 = (-coef2 - coef3) / coef1;
1129: /***********************************************************************/
1130: /* Compute objective value for (t_soln + root1 * e_vect) */
1131: /***********************************************************************/
1133: for (i = 0; i < t_size; ++i) {
1134: e_rwrk[i] = t_soln[i] + root1 * e_vect[i];
1135: }
1136:
1137: obj1 = e_rwrk[0]*(0.5*(cg->diag[0]*e_rwrk[0]+
1138: cg->offd[1]*e_rwrk[1])+cg->norm_r[0]);
1139: for (i = 1; i < t_size - 1; ++i) {
1140: obj1 += 0.5*e_rwrk[i]*(cg->offd[i]*e_rwrk[i-1]+
1141: cg->diag[i]*e_rwrk[i]+
1142: cg->offd[i+1]*e_rwrk[i+1]);
1143: }
1144: obj1 += 0.5*e_rwrk[i]*(cg->offd[i]*e_rwrk[i-1]+
1145: cg->diag[i]*e_rwrk[i]);
1147: /***********************************************************************/
1148: /* Compute objective value for (t_soln + root2 * e_vect) */
1149: /***********************************************************************/
1151: for (i = 0; i < t_size; ++i) {
1152: e_rwrk[i] = t_soln[i] + root2 * e_vect[i];
1153: }
1154:
1155: obj2 = e_rwrk[0]*(0.5*(cg->diag[0]*e_rwrk[0]+
1156: cg->offd[1]*e_rwrk[1])+cg->norm_r[0]);
1157: for (i = 1; i < t_size - 1; ++i) {
1158: obj2 += 0.5*e_rwrk[i]*(cg->offd[i]*e_rwrk[i-1]+
1159: cg->diag[i]*e_rwrk[i]+
1160: cg->offd[i+1]*e_rwrk[i+1]);
1161: }
1162: obj2 += 0.5*e_rwrk[i]*(cg->offd[i]*e_rwrk[i-1]+
1163: cg->diag[i]*e_rwrk[i]);
1164:
1165: /***********************************************************************/
1166: /* Choose the point with the best objective function value. */
1167: /***********************************************************************/
1169: if (obj1 <= obj2) {
1170: for (i = 0; i < t_size; ++i) {
1171: t_soln[i] += root1 * e_vect[i];
1172: }
1173: }
1174: else {
1175: for (i = 0; i < t_size; ++i) {
1176: t_soln[i] += root2 * e_vect[i];
1177: }
1178: }
1179: }
1180: else {
1181: /***********************************************************************/
1182: /* The matrix is positive definite or there was no room to move; the */
1183: /* solution is already contained in t_soln. */
1184: /***********************************************************************/
1185: }
1186: }
1187: else {
1188: /*************************************************************************/
1189: /* The step is outside the trust-region. Compute the correct value for */
1190: /* lambda by performing Newton's method. */
1191: /*************************************************************************/
1193: for (i = 0; i < max_newton_its; ++i) {
1194: /***********************************************************************/
1195: /* Check for convergence. */
1196: /***********************************************************************/
1198: if (fabs(norm_t - cg->radius) <= cg->newton_tol * cg->radius) {
1199: break;
1200: }
1202: /***********************************************************************/
1203: /* Compute the update. */
1204: /***********************************************************************/
1206: PetscMemcpy(e_rwrk, t_soln, sizeof(PetscReal)*t_size);
1208: #if defined(PETSC_MISSING_LAPACK_PTTRS)
1209: SETERRQ(PETSC_ERR_SUP,"PTTRS - Lapack routine is unavailable.");
1210: #else
1211: LAPACKpttrs_(&t_size, &nrhs, t_diag, t_offd + 1, e_rwrk, &nldb, &info);
1212: #endif
1214: if (0 != info) {
1215: /*********************************************************************/
1216: /* Calculation of the step failed; return the Steihaug-Toint */
1217: /* direction. */
1218: /*********************************************************************/
1220: PetscInfo(ksp, "KSPSolve_GLTR: failed to compute step.\n");
1221: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
1222: return(0);
1223: }
1225: /***********************************************************************/
1226: /* Modify lambda. */
1227: /***********************************************************************/
1229: norm_w = 0;
1230: for (j = 0; j < t_size; ++j) {
1231: norm_w += t_soln[j] * e_rwrk[j];
1232: }
1233:
1234: cg->lambda += (norm_t - cg->radius)/cg->radius * (norm_t * norm_t) / norm_w;
1236: /***********************************************************************/
1237: /* Factor T + lambda I */
1238: /***********************************************************************/
1239:
1240: for (j = 0; j < t_size; ++j) {
1241: t_diag[j] = cg->diag[j] + cg->lambda;
1242: t_offd[j] = cg->offd[j];
1243: }
1245: #if defined(PETSC_MISSING_LAPACK_PTTRF)
1246: SETERRQ(PETSC_ERR_SUP,"PTTRF - Lapack routine is unavailable.");
1247: #else
1248: LAPACKpttrf_(&t_size, t_diag, t_offd + 1, &info);
1249: #endif
1251: if (0 != info) {
1252: /*********************************************************************/
1253: /* Calculation of factorization failed; return the Steihaug-Toint */
1254: /* direction. */
1255: /*********************************************************************/
1257: PetscInfo(ksp, "KSPSolve_GLTR: factorization failed.\n");
1258: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
1259: return(0);
1260: }
1262: /***********************************************************************/
1263: /* Compute the new step and its norm. */
1264: /***********************************************************************/
1266: t_soln[0] = -cg->norm_r[0];
1267: for (j = 1; j < t_size; ++j) {
1268: t_soln[j] = 0.0;
1269: }
1271: #if defined(PETSC_MISSING_LAPACK_PTTRS)
1272: SETERRQ(PETSC_ERR_SUP,"PTTRS - Lapack routine is unavailable.");
1273: #else
1274: LAPACKpttrs_(&t_size, &nrhs, t_diag, t_offd + 1, t_soln, &nldb, &info);
1275: #endif
1277: if (0 != info) {
1278: /*********************************************************************/
1279: /* Calculation of the step failed; return the Steihaug-Toint */
1280: /* direction. */
1281: /*********************************************************************/
1282:
1283: PetscInfo(ksp, "KSPSolve_GLTR: failed to compute step.\n");
1284: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
1285: return(0);
1286: }
1288: norm_t = 0;
1289: for (j = 0; j < t_size; ++j) {
1290: norm_t += t_soln[j] * t_soln[j];
1291: }
1292: norm_t = sqrt(norm_t);
1293: }
1295: /*************************************************************************/
1296: /* Check for convergence. */
1297: /*************************************************************************/
1299: if (fabs(norm_t - cg->radius) > cg->newton_tol * cg->radius) {
1300: /***********************************************************************/
1301: /* Newton method failed to converge in iteration limit. */
1302: /***********************************************************************/
1304: PetscInfo(ksp, "KSPSolve_GLTR: failed to converge.\n");
1305: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
1306: return(0);
1307: }
1308: }
1310: /***************************************************************************/
1311: /* Recover the norm of the direction and objective function value. */
1312: /***************************************************************************/
1314: cg->norm_d = norm_t;
1316: cg->o_fcn = t_soln[0]*(0.5*(cg->diag[0]*t_soln[0]+
1317: cg->offd[1]*t_soln[1])+cg->norm_r[0]);
1318: for (i = 1; i < t_size - 1; ++i) {
1319: cg->o_fcn += 0.5*t_soln[i]*(cg->offd[i]*t_soln[i-1]+
1320: cg->diag[i]*t_soln[i]+
1321: cg->offd[i+1]*t_soln[i+1]);
1322: }
1323: cg->o_fcn += 0.5*t_soln[i]*(cg->offd[i]*t_soln[i-1]+
1324: cg->diag[i]*t_soln[i]);
1326: /***************************************************************************/
1327: /* Recover the direction. */
1328: /***************************************************************************/
1330: sigma = -1;
1332: /***************************************************************************/
1333: /* Start conjugate gradient method from the beginning */
1334: /***************************************************************************/
1336: VecCopy(ksp->vec_rhs, r); /* r = -grad */
1337: KSP_PCApply(ksp, r, z); /* z = inv(M) r */
1339: /***************************************************************************/
1340: /* Accumulate Q * s */
1341: /***************************************************************************/
1343: VecCopy(z, d);
1344: VecScale(d, sigma * t_soln[0] / cg->norm_r[0]);
1345:
1346: /***************************************************************************/
1347: /* Compute the first direction. */
1348: /***************************************************************************/
1350: VecCopy(z, p); /* p = z */
1351: KSP_MatMult(ksp, Qmat, p, z); /* z = Q * p */
1352: ++ksp->its;
1354: for (i = 0; i < l_size - 1; ++i) {
1355: /*************************************************************************/
1356: /* Update the residual and direction. */
1357: /*************************************************************************/
1359: alpha = cg->alpha[i];
1360: if (alpha >= 0.0) {
1361: sigma = -sigma;
1362: }
1364: VecAXPY(r, -alpha, z); /* r = r - alpha Q p */
1365: KSP_PCApply(ksp, r, z); /* z = inv(M) r */
1367: /*************************************************************************/
1368: /* Accumulate Q * s */
1369: /*************************************************************************/
1371: VecAXPY(d, sigma * t_soln[i+1] / cg->norm_r[i+1], z);
1373: /*************************************************************************/
1374: /* Update p. */
1375: /*************************************************************************/
1377: beta = cg->beta[i];
1378: VecAYPX(p, beta, z); /* p = z + beta p */
1379: KSP_MatMult(ksp, Qmat, p, z); /* z = Q * p */
1380: ++ksp->its;
1381: }
1383: /***************************************************************************/
1384: /* Update the residual and direction. */
1385: /***************************************************************************/
1387: alpha = cg->alpha[i];
1388: if (alpha >= 0.0) {
1389: sigma = -sigma;
1390: }
1392: VecAXPY(r, -alpha, z); /* r = r - alpha Q p */
1393: KSP_PCApply(ksp, r, z); /* z = inv(M) r */
1395: /***************************************************************************/
1396: /* Accumulate Q * s */
1397: /***************************************************************************/
1399: VecAXPY(d, sigma * t_soln[i+1] / cg->norm_r[i+1], z);
1401: /***************************************************************************/
1402: /* Set the termination reason. */
1403: /***************************************************************************/
1405: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE;
1406: return(0);
1407: #endif
1408: }
1412: PetscErrorCode KSPSetUp_GLTR(KSP ksp)
1413: {
1414: KSP_GLTR *cg = (KSP_GLTR *)ksp->data;
1415: PetscInt max_its, size;
1420: /***************************************************************************/
1421: /* This implementation of CG only handles left preconditioning so generate */
1422: /* an error otherwise. */
1423: /***************************************************************************/
1425: if (ksp->pc_side == PC_RIGHT) {
1426: SETERRQ(PETSC_ERR_SUP, "No right preconditioning for KSPGLTR");
1427: } else if (ksp->pc_side == PC_SYMMETRIC) {
1428: SETERRQ(PETSC_ERR_SUP, "No symmetric preconditioning for KSPGLTR");
1429: }
1431: /***************************************************************************/
1432: /* Determine the total maximum number of iterations. */
1433: /***************************************************************************/
1435: max_its = ksp->max_it + cg->max_lanczos_its + 1;
1437: /***************************************************************************/
1438: /* Set work vectors needed by conjugate gradient method and allocate */
1439: /* workspace for Lanczos matrix. */
1440: /***************************************************************************/
1442: KSPDefaultGetWork(ksp, 3);
1444: size = 5*max_its*sizeof(PetscReal);
1445: PetscMalloc(size, &cg->diag);
1446: PetscMemzero(cg->diag, size);
1447: PetscLogObjectMemory(ksp, size);
1449: cg->offd = cg->diag + max_its;
1450: cg->alpha = cg->offd + max_its;
1451: cg->beta = cg->alpha + max_its;
1452: cg->norm_r = cg->beta + max_its;
1453: return(0);
1454: }
1458: PetscErrorCode KSPDestroy_GLTR(KSP ksp)
1459: {
1460: KSP_GLTR *cg = (KSP_GLTR *)ksp->data;
1465: /***************************************************************************/
1466: /* Free memory allocated for the data. */
1467: /***************************************************************************/
1469: PetscFree(cg->diag);
1470: if (cg->alloced) {
1471: PetscFree(cg->rwork);
1472: PetscFree(cg->iwork);
1473: }
1475: /***************************************************************************/
1476: /* Clear composed functions */
1477: /***************************************************************************/
1479: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGLTRSetRadius_C","",PETSC_NULL);
1480: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGLTRGetNormD_C","",PETSC_NULL);
1481: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGLTRGetObjFcn_C","",PETSC_NULL);
1482: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGLTRGetMinEig_C","",PETSC_NULL);
1483: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGLTRGetLambda_C","",PETSC_NULL);
1485: /***************************************************************************/
1486: /* Destroy KSP object. */
1487: /***************************************************************************/
1488: KSPDefaultDestroy(ksp);
1489: return(0);
1490: }
1495: PetscErrorCode KSPGLTRSetRadius_GLTR(KSP ksp, PetscReal radius)
1496: {
1497: KSP_GLTR *cg = (KSP_GLTR *)ksp->data;
1500: cg->radius = radius;
1501: return(0);
1502: }
1506: PetscErrorCode KSPGLTRGetNormD_GLTR(KSP ksp, PetscReal *norm_d)
1507: {
1508: KSP_GLTR *cg = (KSP_GLTR *)ksp->data;
1511: *norm_d = cg->norm_d;
1512: return(0);
1513: }
1517: PetscErrorCode KSPGLTRGetObjFcn_GLTR(KSP ksp, PetscReal *o_fcn)
1518: {
1519: KSP_GLTR *cg = (KSP_GLTR *)ksp->data;
1522: *o_fcn = cg->o_fcn;
1523: return(0);
1524: }
1528: PetscErrorCode KSPGLTRGetMinEig_GLTR(KSP ksp, PetscReal *e_min)
1529: {
1530: KSP_GLTR *cg = (KSP_GLTR *)ksp->data;
1533: *e_min = cg->e_min;
1534: return(0);
1535: }
1539: PetscErrorCode KSPGLTRGetLambda_GLTR(KSP ksp, PetscReal *lambda)
1540: {
1541: KSP_GLTR *cg = (KSP_GLTR *)ksp->data;
1544: *lambda = cg->lambda;
1545: return(0);
1546: }
1551: PetscErrorCode KSPSetFromOptions_GLTR(KSP ksp)
1552: {
1554: KSP_GLTR *cg = (KSP_GLTR *)ksp->data;
1557: PetscOptionsHead("KSP GLTR options");
1559: PetscOptionsReal("-ksp_gltr_radius", "Trust Region Radius", "KSPGLTRSetRadius", cg->radius, &cg->radius, PETSC_NULL);
1561: PetscOptionsReal("-ksp_gltr_init_pert", "Initial perturbation", "", cg->init_pert, &cg->init_pert, PETSC_NULL);
1562: PetscOptionsReal("-ksp_gltr_eigen_tol", "Eigenvalue tolerance", "", cg->eigen_tol, &cg->eigen_tol, PETSC_NULL);
1563: PetscOptionsReal("-ksp_gltr_newton_tol", "Newton tolerance", "", cg->newton_tol, &cg->newton_tol, PETSC_NULL);
1565: PetscOptionsInt("-ksp_gltr_max_lanczos_its", "Maximum Lanczos Iters", "", cg->max_lanczos_its, &cg->max_lanczos_its, PETSC_NULL);
1566: PetscOptionsInt("-ksp_gltr_max_newton_its", "Maximum Newton Iters", "", cg->max_newton_its, &cg->max_newton_its, PETSC_NULL);
1568: PetscOptionsEList("-ksp_gltr_dtype", "Norm used for direction", "", DType_Table, GLTR_DIRECTION_TYPES, DType_Table[cg->dtype], &cg->dtype, PETSC_NULL);
1570: PetscOptionsTail();
1571: return(0);
1572: }
1574: /*MC
1575: KSPGLTR - Code to run conjugate gradient method subject to a constraint
1576: on the solution norm. This is used in Trust Region methods for
1577: nonlinear equations, SNESTR
1579: Options Database Keys:
1580: . -ksp_gltr_radius <r> - Trust Region Radius
1582: Notes: This is rarely used directly
1584: Level: developer
1586: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPGLTRSetRadius(), KSPGLTRGetNormD(), KSPGLTRGetObjFcn(), KSPGLTRGetMinEig(), KSPGLTRGetLambda()
1587: M*/
1592: PetscErrorCode KSPCreate_GLTR(KSP ksp)
1593: {
1595: KSP_GLTR *cg;
1599: PetscNewLog(ksp, KSP_GLTR, &cg);
1601: cg->radius = 0.0;
1602: cg->dtype = GLTR_UNPRECONDITIONED_DIRECTION;
1604: cg->init_pert = 1.0e-8;
1605: cg->eigen_tol = 1.0e-10;
1606: cg->newton_tol = 1.0e-6;
1608: cg->alloced = 0;
1609: cg->init_alloc = 1024;
1611: cg->max_lanczos_its = 20;
1612: cg->max_newton_its = 10;
1614: ksp->data = (void *) cg;
1615: ksp->pc_side = PC_LEFT;
1616: ksp->normtype = KSP_NORM_UNPRECONDITIONED;
1618: /***************************************************************************/
1619: /* Sets the functions that are associated with this data structure */
1620: /* (in C++ this is the same as defining virtual functions). */
1621: /***************************************************************************/
1623: ksp->ops->setup = KSPSetUp_GLTR;
1624: ksp->ops->solve = KSPSolve_GLTR;
1625: ksp->ops->destroy = KSPDestroy_GLTR;
1626: ksp->ops->setfromoptions = KSPSetFromOptions_GLTR;
1627: ksp->ops->buildsolution = KSPDefaultBuildSolution;
1628: ksp->ops->buildresidual = KSPDefaultBuildResidual;
1629: ksp->ops->view = 0;
1631: PetscObjectComposeFunctionDynamic((PetscObject)ksp,
1632: "KSPGLTRSetRadius_C",
1633: "KSPGLTRSetRadius_GLTR",
1634: KSPGLTRSetRadius_GLTR);
1635: PetscObjectComposeFunctionDynamic((PetscObject)ksp,
1636: "KSPGLTRGetNormD_C",
1637: "KSPGLTRGetNormD_GLTR",
1638: KSPGLTRGetNormD_GLTR);
1639: PetscObjectComposeFunctionDynamic((PetscObject)ksp,
1640: "KSPGLTRGetObjFcn_C",
1641: "KSPGLTRGetObjFcn_GLTR",
1642: KSPGLTRGetObjFcn_GLTR);
1643: PetscObjectComposeFunctionDynamic((PetscObject)ksp,
1644: "KSPGLTRGetMinEig_C",
1645: "KSPGLTRGetMinEig_GLTR",
1646: KSPGLTRGetMinEig_GLTR);
1647: PetscObjectComposeFunctionDynamic((PetscObject)ksp,
1648: "KSPGLTRGetLambda_C",
1649: "KSPGLTRGetLambda_GLTR",
1650: KSPGLTRGetLambda_GLTR);
1651: return(0);
1652: }