Actual source code: qcg.c
2: #include <../src/ksp/ksp/impls/qcg/qcgimpl.h>
4: static PetscErrorCode KSPQCGQuadraticRoots(Vec, Vec, PetscReal, PetscReal *, PetscReal *);
6: /*@
7: KSPQCGSetTrustRegionRadius - Sets the radius of the trust region.
9: Logically Collective on ksp
11: Input Parameters:
12: + ksp - the iterative context
13: - delta - the trust region radius (Infinity is the default)
15: Options Database Key:
16: . -ksp_qcg_trustregionradius <delta> - trust region radius
18: Level: advanced
20: @*/
21: PetscErrorCode KSPQCGSetTrustRegionRadius(KSP ksp, PetscReal delta)
22: {
25: PetscTryMethod(ksp, "KSPQCGSetTrustRegionRadius_C", (KSP, PetscReal), (ksp, delta));
26: return 0;
27: }
29: /*@
30: KSPQCGGetTrialStepNorm - Gets the norm of a trial step vector. The WCG step may be
31: constrained, so this is not necessarily the length of the ultimate step taken in QCG.
33: Not Collective
35: Input Parameter:
36: . ksp - the iterative context
38: Output Parameter:
39: . tsnorm - the norm
41: Level: advanced
42: @*/
43: PetscErrorCode KSPQCGGetTrialStepNorm(KSP ksp, PetscReal *tsnorm)
44: {
46: PetscUseMethod(ksp, "KSPQCGGetTrialStepNorm_C", (KSP, PetscReal *), (ksp, tsnorm));
47: return 0;
48: }
50: /*@
51: KSPQCGGetQuadratic - Gets the value of the quadratic function, evaluated at the new iterate:
53: q(s) = g^T * s + 0.5 * s^T * H * s
55: which satisfies the Euclidian Norm trust region constraint
57: || D * s || <= delta,
59: where
61: delta is the trust region radius,
62: g is the gradient vector, and
63: H is Hessian matrix,
64: D is a scaling matrix.
66: Collective on ksp
68: Input Parameter:
69: . ksp - the iterative context
71: Output Parameter:
72: . quadratic - the quadratic function evaluated at the new iterate
74: Level: advanced
75: @*/
76: PetscErrorCode KSPQCGGetQuadratic(KSP ksp, PetscReal *quadratic)
77: {
79: PetscUseMethod(ksp, "KSPQCGGetQuadratic_C", (KSP, PetscReal *), (ksp, quadratic));
80: return 0;
81: }
83: PetscErrorCode KSPSolve_QCG(KSP ksp)
84: {
85: /*
86: Correpondence with documentation above:
87: B = g = gradient,
88: X = s = step
89: Note: This is not coded correctly for complex arithmetic!
90: */
92: KSP_QCG *pcgP = (KSP_QCG *)ksp->data;
93: Mat Amat, Pmat;
94: Vec W, WA, WA2, R, P, ASP, BS, X, B;
95: PetscScalar scal, beta, rntrn, step;
96: PetscReal q1, q2, xnorm, step1, step2, rnrm = 0.0, btx, xtax;
97: PetscReal ptasp, rtr, wtasp, bstp;
98: PetscReal dzero = 0.0, bsnrm = 0.0;
99: PetscInt i, maxit;
100: PC pc = ksp->pc;
101: PetscBool diagonalscale;
103: PCGetDiagonalScale(ksp->pc, &diagonalscale);
107: ksp->its = 0;
108: maxit = ksp->max_it;
109: WA = ksp->work[0];
110: R = ksp->work[1];
111: P = ksp->work[2];
112: ASP = ksp->work[3];
113: BS = ksp->work[4];
114: W = ksp->work[5];
115: WA2 = ksp->work[6];
116: X = ksp->vec_sol;
117: B = ksp->vec_rhs;
121: /* Initialize variables */
122: VecSet(W, 0.0); /* W = 0 */
123: VecSet(X, 0.0); /* X = 0 */
124: PCGetOperators(pc, &Amat, &Pmat);
126: /* Compute: BS = D^{-1} B */
127: PCApplySymmetricLeft(pc, B, BS);
129: if (ksp->normtype != KSP_NORM_NONE) {
130: VecNorm(BS, NORM_2, &bsnrm);
131: KSPCheckNorm(ksp, bsnrm);
132: }
133: PetscObjectSAWsTakeAccess((PetscObject)ksp);
134: ksp->its = 0;
135: ksp->rnorm = bsnrm;
136: PetscObjectSAWsGrantAccess((PetscObject)ksp);
137: KSPLogResidualHistory(ksp, bsnrm);
138: KSPMonitor(ksp, 0, bsnrm);
139: (*ksp->converged)(ksp, 0, bsnrm, &ksp->reason, ksp->cnvP);
140: if (ksp->reason) return 0;
142: /* Compute the initial scaled direction and scaled residual */
143: VecCopy(BS, R);
144: VecScale(R, -1.0);
145: VecCopy(R, P);
146: VecDotRealPart(R, R, &rtr);
148: for (i = 0; i <= maxit; i++) {
149: PetscObjectSAWsTakeAccess((PetscObject)ksp);
150: ksp->its++;
151: PetscObjectSAWsGrantAccess((PetscObject)ksp);
153: /* Compute: asp = D^{-T}*A*D^{-1}*p */
154: PCApplySymmetricRight(pc, P, WA);
155: KSP_MatMult(ksp, Amat, WA, WA2);
156: PCApplySymmetricLeft(pc, WA2, ASP);
158: /* Check for negative curvature */
159: VecDotRealPart(P, ASP, &ptasp);
160: if (ptasp <= dzero) {
161: /* Scaled negative curvature direction: Compute a step so that
162: ||w + step*p|| = delta and QS(w + step*p) is least */
164: if (!i) {
165: VecCopy(P, X);
166: VecNorm(X, NORM_2, &xnorm);
167: KSPCheckNorm(ksp, xnorm);
168: scal = pcgP->delta / xnorm;
169: VecScale(X, scal);
170: } else {
171: /* Compute roots of quadratic */
172: KSPQCGQuadraticRoots(W, P, pcgP->delta, &step1, &step2);
173: VecDotRealPart(W, ASP, &wtasp);
174: VecDotRealPart(BS, P, &bstp);
175: VecCopy(W, X);
176: q1 = step1 * (bstp + wtasp + .5 * step1 * ptasp);
177: q2 = step2 * (bstp + wtasp + .5 * step2 * ptasp);
178: if (q1 <= q2) {
179: VecAXPY(X, step1, P);
180: } else {
181: VecAXPY(X, step2, P);
182: }
183: }
184: pcgP->ltsnrm = pcgP->delta; /* convergence in direction of */
185: ksp->reason = KSP_CONVERGED_CG_NEG_CURVE; /* negative curvature */
186: if (!i) {
187: PetscInfo(ksp, "negative curvature: delta=%g\n", (double)pcgP->delta);
188: } else {
189: PetscInfo(ksp, "negative curvature: step1=%g, step2=%g, delta=%g\n", (double)step1, (double)step2, (double)pcgP->delta);
190: }
192: } else {
193: /* Compute step along p */
194: step = rtr / ptasp;
195: VecCopy(W, X); /* x = w */
196: VecAXPY(X, step, P); /* x <- step*p + x */
197: VecNorm(X, NORM_2, &pcgP->ltsnrm);
198: KSPCheckNorm(ksp, pcgP->ltsnrm);
200: if (pcgP->ltsnrm > pcgP->delta) {
201: /* Since the trial iterate is outside the trust region,
202: evaluate a constrained step along p so that
203: ||w + step*p|| = delta
204: The positive step is always better in this case. */
205: if (!i) {
206: scal = pcgP->delta / pcgP->ltsnrm;
207: VecScale(X, scal);
208: } else {
209: /* Compute roots of quadratic */
210: KSPQCGQuadraticRoots(W, P, pcgP->delta, &step1, &step2);
211: VecCopy(W, X);
212: VecAXPY(X, step1, P); /* x <- step1*p + x */
213: }
214: pcgP->ltsnrm = pcgP->delta;
215: ksp->reason = KSP_CONVERGED_CG_CONSTRAINED; /* convergence along constrained step */
216: if (!i) {
217: PetscInfo(ksp, "constrained step: delta=%g\n", (double)pcgP->delta);
218: } else {
219: PetscInfo(ksp, "constrained step: step1=%g, step2=%g, delta=%g\n", (double)step1, (double)step2, (double)pcgP->delta);
220: }
222: } else {
223: /* Evaluate the current step */
224: VecCopy(X, W); /* update interior iterate */
225: VecAXPY(R, -step, ASP); /* r <- -step*asp + r */
226: if (ksp->normtype != KSP_NORM_NONE) {
227: VecNorm(R, NORM_2, &rnrm);
228: KSPCheckNorm(ksp, rnrm);
229: }
230: PetscObjectSAWsTakeAccess((PetscObject)ksp);
231: ksp->rnorm = rnrm;
232: PetscObjectSAWsGrantAccess((PetscObject)ksp);
233: KSPLogResidualHistory(ksp, rnrm);
234: KSPMonitor(ksp, i + 1, rnrm);
235: (*ksp->converged)(ksp, i + 1, rnrm, &ksp->reason, ksp->cnvP);
236: if (ksp->reason) { /* convergence for */
237: PetscInfo(ksp, "truncated step: step=%g, rnrm=%g, delta=%g\n", (double)PetscRealPart(step), (double)rnrm, (double)pcgP->delta);
238: }
239: }
240: }
241: if (ksp->reason) break; /* Convergence has been attained */
242: else { /* Compute a new AS-orthogonal direction */ VecDot(R, R, &rntrn);
243: beta = rntrn / rtr;
244: VecAYPX(P, beta, R); /* p <- r + beta*p */
245: rtr = PetscRealPart(rntrn);
246: }
247: }
248: if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
250: /* Unscale x */
251: VecCopy(X, WA2);
252: PCApplySymmetricRight(pc, WA2, X);
254: KSP_MatMult(ksp, Amat, X, WA);
255: VecDotRealPart(B, X, &btx);
256: VecDotRealPart(X, WA, &xtax);
258: pcgP->quadratic = btx + .5 * xtax;
259: return 0;
260: }
262: PetscErrorCode KSPSetUp_QCG(KSP ksp)
263: {
264: /* Get work vectors from user code */
265: KSPSetWorkVecs(ksp, 7);
266: return 0;
267: }
269: PetscErrorCode KSPDestroy_QCG(KSP ksp)
270: {
271: PetscObjectComposeFunction((PetscObject)ksp, "KSPQCGGetQuadratic_C", NULL);
272: PetscObjectComposeFunction((PetscObject)ksp, "KSPQCGGetTrialStepNorm_C", NULL);
273: PetscObjectComposeFunction((PetscObject)ksp, "KSPQCGSetTrustRegionRadius_C", NULL);
274: KSPDestroyDefault(ksp);
275: return 0;
276: }
278: static PetscErrorCode KSPQCGSetTrustRegionRadius_QCG(KSP ksp, PetscReal delta)
279: {
280: KSP_QCG *cgP = (KSP_QCG *)ksp->data;
282: cgP->delta = delta;
283: return 0;
284: }
286: static PetscErrorCode KSPQCGGetTrialStepNorm_QCG(KSP ksp, PetscReal *ltsnrm)
287: {
288: KSP_QCG *cgP = (KSP_QCG *)ksp->data;
290: *ltsnrm = cgP->ltsnrm;
291: return 0;
292: }
294: static PetscErrorCode KSPQCGGetQuadratic_QCG(KSP ksp, PetscReal *quadratic)
295: {
296: KSP_QCG *cgP = (KSP_QCG *)ksp->data;
298: *quadratic = cgP->quadratic;
299: return 0;
300: }
302: PetscErrorCode KSPSetFromOptions_QCG(KSP ksp, PetscOptionItems *PetscOptionsObject)
303: {
304: PetscReal delta;
305: KSP_QCG *cgP = (KSP_QCG *)ksp->data;
306: PetscBool flg;
308: PetscOptionsHeadBegin(PetscOptionsObject, "KSP QCG Options");
309: PetscOptionsReal("-ksp_qcg_trustregionradius", "Trust Region Radius", "KSPQCGSetTrustRegionRadius", cgP->delta, &delta, &flg);
310: if (flg) KSPQCGSetTrustRegionRadius(ksp, delta);
311: PetscOptionsHeadEnd();
312: return 0;
313: }
315: /*MC
316: KSPQCG - Code to run conjugate gradient method subject to a constraint
317: on the solution norm. This is used in Trust Region methods for nonlinear equations, SNESNEWTONTR
319: Options Database Keys:
320: . -ksp_qcg_trustregionradius <r> - Trust Region Radius
322: Notes:
323: This is rarely used directly
325: Level: developer
327: Notes:
328: Use preconditioned conjugate gradient to compute
329: an approximate minimizer of the quadratic function
331: q(s) = g^T * s + .5 * s^T * H * s
333: subject to the Euclidean norm trust region constraint
335: || D * s || <= delta,
337: where
339: delta is the trust region radius,
340: g is the gradient vector, and
341: H is Hessian matrix,
342: D is a scaling matrix.
344: KSPConvergedReason may be
345: $ KSP_CONVERGED_CG_NEG_CURVE if convergence is reached along a negative curvature direction,
346: $ KSP_CONVERGED_CG_CONSTRAINED if convergence is reached along a constrained step,
347: $ other KSP converged/diverged reasons
349: Notes:
350: Currently we allow symmetric preconditioning with the following scaling matrices:
351: PCNONE: D = Identity matrix
352: PCJACOBI: D = diag [d_1, d_2, ...., d_n], where d_i = sqrt(H[i,i])
353: PCICC: D = L^T, implemented with forward and backward solves.
354: Here L is an incomplete Cholesky factor of H.
356: References:
357: . * - Trond Steihaug, The Conjugate Gradient Method and Trust Regions in Large Scale Optimization,
358: SIAM Journal on Numerical Analysis, Vol. 20, No. 3 (Jun., 1983).
360: .seealso: `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPQCGSetTrustRegionRadius()`
361: `KSPQCGGetTrialStepNorm()`, `KSPQCGGetQuadratic()`
362: M*/
364: PETSC_EXTERN PetscErrorCode KSPCreate_QCG(KSP ksp)
365: {
366: KSP_QCG *cgP;
368: KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_SYMMETRIC, 3);
369: KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_SYMMETRIC, 1);
370: PetscNew(&cgP);
372: ksp->data = (void *)cgP;
373: ksp->ops->setup = KSPSetUp_QCG;
374: ksp->ops->setfromoptions = KSPSetFromOptions_QCG;
375: ksp->ops->solve = KSPSolve_QCG;
376: ksp->ops->destroy = KSPDestroy_QCG;
377: ksp->ops->buildsolution = KSPBuildSolutionDefault;
378: ksp->ops->buildresidual = KSPBuildResidualDefault;
379: ksp->ops->view = NULL;
381: PetscObjectComposeFunction((PetscObject)ksp, "KSPQCGGetQuadratic_C", KSPQCGGetQuadratic_QCG);
382: PetscObjectComposeFunction((PetscObject)ksp, "KSPQCGGetTrialStepNorm_C", KSPQCGGetTrialStepNorm_QCG);
383: PetscObjectComposeFunction((PetscObject)ksp, "KSPQCGSetTrustRegionRadius_C", KSPQCGSetTrustRegionRadius_QCG);
384: cgP->delta = PETSC_MAX_REAL; /* default trust region radius is infinite */
385: return 0;
386: }
388: /* ---------------------------------------------------------- */
389: /*
390: KSPQCGQuadraticRoots - Computes the roots of the quadratic,
391: ||s + step*p|| - delta = 0
392: such that step1 >= 0 >= step2.
393: where
394: delta:
395: On entry delta must contain scalar delta.
396: On exit delta is unchanged.
397: step1:
398: On entry step1 need not be specified.
399: On exit step1 contains the non-negative root.
400: step2:
401: On entry step2 need not be specified.
402: On exit step2 contains the non-positive root.
403: C code is translated from the Fortran version of the MINPACK-2 Project,
404: Argonne National Laboratory, Brett M. Averick and Richard G. Carter.
405: */
406: static PetscErrorCode KSPQCGQuadraticRoots(Vec s, Vec p, PetscReal delta, PetscReal *step1, PetscReal *step2)
407: {
408: PetscReal dsq, ptp, pts, rad, sts;
410: VecDotRealPart(p, s, &pts);
411: VecDotRealPart(p, p, &ptp);
412: VecDotRealPart(s, s, &sts);
413: dsq = delta * delta;
414: rad = PetscSqrtReal((pts * pts) - ptp * (sts - dsq));
415: if (pts > 0.0) {
416: *step2 = -(pts + rad) / ptp;
417: *step1 = (sts - dsq) / (ptp * *step2);
418: } else {
419: *step1 = -(pts - rad) / ptp;
420: *step2 = (sts - dsq) / (ptp * *step1);
421: }
422: return 0;
423: }