Actual source code: ntrdc.c
2: #include <../src/snes/impls/ntrdc/ntrdcimpl.h>
4: typedef struct {
5: SNES snes;
6: /* Information on the regular SNES convergence test; which may have been user provided
7: Copied from tr.c (maybe able to disposed, but this is a private function) - Heeho
8: Same with SNESTR_KSPConverged_Private, SNESTR_KSPConverged_Destroy, and SNESTR_Converged_Private
9: */
11: PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
12: PetscErrorCode (*convdestroy)(void *);
13: void *convctx;
14: } SNES_TRDC_KSPConverged_Ctx;
16: static PetscErrorCode SNESTRDC_KSPConverged_Private(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *cctx)
17: {
18: SNES_TRDC_KSPConverged_Ctx *ctx = (SNES_TRDC_KSPConverged_Ctx *)cctx;
19: SNES snes = ctx->snes;
20: SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC *)snes->data;
21: Vec x;
22: PetscReal nrm;
24: PetscFunctionBegin;
25: PetscCall((*ctx->convtest)(ksp, n, rnorm, reason, ctx->convctx));
26: if (*reason) PetscCall(PetscInfo(snes, "Default or user provided convergence test KSP iterations=%" PetscInt_FMT ", rnorm=%g\n", n, (double)rnorm));
27: /* Determine norm of solution */
28: PetscCall(KSPBuildSolution(ksp, NULL, &x));
29: PetscCall(VecNorm(x, NORM_2, &nrm));
30: if (nrm >= neP->delta) {
31: PetscCall(PetscInfo(snes, "Ending linear iteration early, delta=%g, length=%g\n", (double)neP->delta, (double)nrm));
32: *reason = KSP_CONVERGED_STEP_LENGTH;
33: }
34: PetscFunctionReturn(PETSC_SUCCESS);
35: }
37: static PetscErrorCode SNESTRDC_KSPConverged_Destroy(void *cctx)
38: {
39: SNES_TRDC_KSPConverged_Ctx *ctx = (SNES_TRDC_KSPConverged_Ctx *)cctx;
41: PetscFunctionBegin;
42: PetscCall((*ctx->convdestroy)(ctx->convctx));
43: PetscCall(PetscFree(ctx));
45: PetscFunctionReturn(PETSC_SUCCESS);
46: }
48: /*
49: SNESTRDC_Converged_Private -test convergence JUST for
50: the trust region tolerance.
52: */
53: static PetscErrorCode SNESTRDC_Converged_Private(SNES snes, PetscInt it, PetscReal xnorm, PetscReal pnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy)
54: {
55: SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC *)snes->data;
57: PetscFunctionBegin;
58: *reason = SNES_CONVERGED_ITERATING;
59: if (neP->delta < xnorm * snes->deltatol) {
60: PetscCall(PetscInfo(snes, "Diverged due to too small a trust region %g<%g*%g\n", (double)neP->delta, (double)xnorm, (double)snes->deltatol));
61: *reason = SNES_DIVERGED_TR_DELTA;
62: } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
63: PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n", snes->max_funcs));
64: *reason = SNES_DIVERGED_FUNCTION_COUNT;
65: }
66: PetscFunctionReturn(PETSC_SUCCESS);
67: }
69: /*@
70: SNESNewtonTRDCGetRhoFlag - Get whether the current solution update is within the trust-region.
72: Input Parameter:
73: . snes - the nonlinear solver object
75: Output Parameter:
76: . rho_flag: `PETSC_TRUE` if the solution update is in the trust-region; otherwise, `PETSC_FALSE`
78: Level: developer
80: .seealso: `SNESNEWTONTRDC`, `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, , `SNESNewtonTRDCSetPreCheck()`,
81: `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`
82: @*/
83: PetscErrorCode SNESNewtonTRDCGetRhoFlag(SNES snes, PetscBool *rho_flag)
84: {
85: SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
87: PetscFunctionBegin;
90: *rho_flag = tr->rho_satisfied;
91: PetscFunctionReturn(PETSC_SUCCESS);
92: }
94: /*@C
95: SNESNewtonTRDCSetPreCheck - Sets a user function that is called before the search step has been determined.
96: Allows the user a chance to change or override the trust region decision.
98: Logically Collective
100: Input Parameters:
101: + snes - the nonlinear solver object
102: . func - [optional] function evaluation routine, see `SNESNewtonTRDCPreCheck()` for the calling sequence
103: - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
105: Level: intermediate
107: Note:
108: This function is called BEFORE the function evaluation within the `SNESNEWTONTRDC` solver.
110: .seealso: `SNESNEWTONTRDC`, `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`,
111: `SNESNewtonTRDCGetRhoFlag()`
112: @*/
113: PetscErrorCode SNESNewtonTRDCSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, PetscBool *, void *), void *ctx)
114: {
115: SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
117: PetscFunctionBegin;
119: if (func) tr->precheck = func;
120: if (ctx) tr->precheckctx = ctx;
121: PetscFunctionReturn(PETSC_SUCCESS);
122: }
124: /*@C
125: SNESNewtonTRDCGetPreCheck - Gets the pre-check function
127: Not collective
129: Input Parameter:
130: . snes - the nonlinear solver context
132: Output Parameters:
133: + func - [optional] function evaluation routine, see for the calling sequence `SNESNewtonTRDCPreCheck()`
134: - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
136: Level: intermediate
138: .seealso: `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCPreCheck()`
139: @*/
140: PetscErrorCode SNESNewtonTRDCGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, PetscBool *, void *), void **ctx)
141: {
142: SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
144: PetscFunctionBegin;
146: if (func) *func = tr->precheck;
147: if (ctx) *ctx = tr->precheckctx;
148: PetscFunctionReturn(PETSC_SUCCESS);
149: }
151: /*@C
152: SNESNewtonTRDCSetPostCheck - Sets a user function that is called after the search step has been determined but before the next
153: function evaluation. Allows the user a chance to change or override the decision of the line search routine
155: Logically Collective
157: Input Parameters:
158: + snes - the nonlinear solver object
159: . func - [optional] function evaluation routine, see `SNESNewtonTRDCPostCheck()` for the calling sequence
160: - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
162: Level: intermediate
164: Note:
165: This function is called BEFORE the function evaluation within the `SNESNEWTONTRDC` solver while the function set in
166: `SNESLineSearchSetPostCheck()` is called AFTER the function evaluation.
168: .seealso: `SNESNEWTONTRDC`, `SNESNewtonTRDCPostCheck()`, `SNESNewtonTRDCGetPostCheck()`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCGetPreCheck()`
169: @*/
170: PetscErrorCode SNESNewtonTRDCSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx)
171: {
172: SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
174: PetscFunctionBegin;
176: if (func) tr->postcheck = func;
177: if (ctx) tr->postcheckctx = ctx;
178: PetscFunctionReturn(PETSC_SUCCESS);
179: }
181: /*@C
182: SNESNewtonTRDCGetPostCheck - Gets the post-check function
184: Not collective
186: Input Parameter:
187: . snes - the nonlinear solver context
189: Output Parameters:
190: + func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRDCPostCheck()
191: - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
193: Level: intermediate
195: .seealso: `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCPostCheck()`
196: @*/
197: PetscErrorCode SNESNewtonTRDCGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx)
198: {
199: SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
201: PetscFunctionBegin;
203: if (func) *func = tr->postcheck;
204: if (ctx) *ctx = tr->postcheckctx;
205: PetscFunctionReturn(PETSC_SUCCESS);
206: }
208: /*@C
209: SNESNewtonTRDCPreCheck - Called before the step has been determined in `SNESNEWTONTRDC`
211: Logically Collective
213: Input Parameters:
214: + snes - the solver
215: . X - The last solution
216: - Y - The step direction
218: Output Parameters:
219: . changed_Y - Indicator that the step direction Y has been changed.
221: Level: developer
223: .seealso: `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCPostCheck()`
224: @*/
225: static PetscErrorCode SNESNewtonTRDCPreCheck(SNES snes, Vec X, Vec Y, PetscBool *changed_Y)
226: {
227: SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
229: PetscFunctionBegin;
230: *changed_Y = PETSC_FALSE;
231: if (tr->precheck) {
232: PetscCall((*tr->precheck)(snes, X, Y, changed_Y, tr->precheckctx));
234: }
235: PetscFunctionReturn(PETSC_SUCCESS);
236: }
238: /*@C
239: SNESNewtonTRDCPostCheck - Called after the step has been determined in `SNESNEWTONTRDC` but before the function evaluation at that step
241: Logically Collective
243: Input Parameters:
244: + snes - the solver
245: . X - The last solution
246: . Y - The full step direction
247: - W - The updated solution, W = X - Y
249: Output Parameters:
250: + changed_Y - indicator if step has been changed
251: - changed_W - Indicator if the new candidate solution W has been changed.
253: Note:
254: If Y is changed then W is recomputed as X - Y
256: Level: developer
258: .seealso: `SNESNEWTONTRDC`, `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`, `SNESNewtonTRDCPreCheck()
259: @*/
260: static PetscErrorCode SNESNewtonTRDCPostCheck(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W)
261: {
262: SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
264: PetscFunctionBegin;
265: *changed_Y = PETSC_FALSE;
266: *changed_W = PETSC_FALSE;
267: if (tr->postcheck) {
268: PetscCall((*tr->postcheck)(snes, X, Y, W, changed_Y, changed_W, tr->postcheckctx));
271: }
272: PetscFunctionReturn(PETSC_SUCCESS);
273: }
275: /*
276: SNESSolve_NEWTONTRDC - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy
277: (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of
278: nonlinear equations
280: */
281: static PetscErrorCode SNESSolve_NEWTONTRDC(SNES snes)
282: {
283: SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC *)snes->data;
284: Vec X, F, Y, G, W, GradF, YNtmp;
285: Vec YCtmp;
286: Mat jac;
287: PetscInt maxits, i, j, lits, inner_count, bs;
288: PetscReal rho, fnorm, gnorm, xnorm = 0, delta, ynorm, temp_xnorm, temp_ynorm; /* TRDC inner iteration */
289: PetscReal inorms[99]; /* need to make it dynamic eventually, fixed max block size of 99 for now */
290: PetscReal deltaM, ynnorm, f0, mp, gTy, g, yTHy; /* rho calculation */
291: PetscReal auk, gfnorm, ycnorm, c0, c1, c2, tau, tau_pos, tau_neg, gTBg; /* Cauchy Point */
292: KSP ksp;
293: SNESConvergedReason reason = SNES_CONVERGED_ITERATING;
294: PetscBool breakout = PETSC_FALSE;
295: SNES_TRDC_KSPConverged_Ctx *ctx;
296: PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), (*convdestroy)(void *);
297: void *convctx;
299: PetscFunctionBegin;
300: maxits = snes->max_its; /* maximum number of iterations */
301: X = snes->vec_sol; /* solution vector */
302: F = snes->vec_func; /* residual vector */
303: Y = snes->work[0]; /* update vector */
304: G = snes->work[1]; /* updated residual */
305: W = snes->work[2]; /* temporary vector */
306: GradF = snes->work[3]; /* grad f = J^T F */
307: YNtmp = snes->work[4]; /* Newton solution */
308: YCtmp = snes->work[5]; /* Cauchy solution */
310: PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
312: PetscCall(VecGetBlockSize(YNtmp, &bs));
314: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
315: snes->iter = 0;
316: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
318: /* Set the linear stopping criteria to use the More' trick. From tr.c */
319: PetscCall(SNESGetKSP(snes, &ksp));
320: PetscCall(KSPGetConvergenceTest(ksp, &convtest, &convctx, &convdestroy));
321: if (convtest != SNESTRDC_KSPConverged_Private) {
322: PetscCall(PetscNew(&ctx));
323: ctx->snes = snes;
324: PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy));
325: PetscCall(KSPSetConvergenceTest(ksp, SNESTRDC_KSPConverged_Private, ctx, SNESTRDC_KSPConverged_Destroy));
326: PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTRDC_KSPConverged_Private\n"));
327: }
329: if (!snes->vec_func_init_set) {
330: PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */
331: } else snes->vec_func_init_set = PETSC_FALSE;
333: PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */
334: SNESCheckFunctionNorm(snes, fnorm);
335: PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */
337: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
338: snes->norm = fnorm;
339: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
340: delta = xnorm ? neP->delta0 * xnorm : neP->delta0; /* initial trust region size scaled by xnorm */
341: deltaM = xnorm ? neP->deltaM * xnorm : neP->deltaM; /* maximum trust region size scaled by xnorm */
342: neP->delta = delta;
343: PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
344: PetscCall(SNESMonitor(snes, 0, fnorm));
346: neP->rho_satisfied = PETSC_FALSE;
348: /* test convergence */
349: PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
350: if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
352: for (i = 0; i < maxits; i++) {
353: PetscBool changed_y;
354: PetscBool changed_w;
356: /* dogleg method */
357: PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
358: SNESCheckJacobianDomainerror(snes);
359: PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian));
360: PetscCall(KSPSolve(snes->ksp, F, YNtmp)); /* Quasi Newton Solution */
361: SNESCheckKSPSolve(snes); /* this is necessary but old tr.c did not have it*/
362: PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
363: PetscCall(SNESGetJacobian(snes, &jac, NULL, NULL, NULL));
365: /* rescale Jacobian, Newton solution update, and re-calculate delta for multiphase (multivariable)
366: for inner iteration and Cauchy direction calculation
367: */
368: if (bs > 1 && neP->auto_scale_multiphase) {
369: PetscCall(VecStrideNormAll(YNtmp, NORM_INFINITY, inorms));
370: for (j = 0; j < bs; j++) {
371: if (neP->auto_scale_max > 1.0) {
372: if (inorms[j] < 1.0 / neP->auto_scale_max) inorms[j] = 1.0 / neP->auto_scale_max;
373: }
374: PetscCall(VecStrideSet(W, j, inorms[j]));
375: PetscCall(VecStrideScale(YNtmp, j, 1.0 / inorms[j]));
376: PetscCall(VecStrideScale(X, j, 1.0 / inorms[j]));
377: }
378: PetscCall(VecNorm(X, NORM_2, &xnorm));
379: if (i == 0) {
380: delta = neP->delta0 * xnorm;
381: } else {
382: delta = neP->delta * xnorm;
383: }
384: deltaM = neP->deltaM * xnorm;
385: PetscCall(MatDiagonalScale(jac, NULL, W));
386: }
388: /* calculating GradF of minimization function */
389: PetscCall(MatMultTranspose(jac, F, GradF)); /* grad f = J^T F */
390: PetscCall(VecNorm(YNtmp, NORM_2, &ynnorm)); /* ynnorm <- || Y_newton || */
392: inner_count = 0;
393: neP->rho_satisfied = PETSC_FALSE;
394: while (1) {
395: if (ynnorm <= delta) { /* see if the Newton solution is within the trust region */
396: PetscCall(VecCopy(YNtmp, Y));
397: } else if (neP->use_cauchy) { /* use Cauchy direction if enabled */
398: PetscCall(MatMult(jac, GradF, W));
399: PetscCall(VecDotRealPart(W, W, &gTBg)); /* completes GradF^T J^T J GradF */
400: PetscCall(VecNorm(GradF, NORM_2, &gfnorm)); /* grad f norm <- || grad f || */
401: if (gTBg <= 0.0) {
402: auk = PETSC_MAX_REAL;
403: } else {
404: auk = PetscSqr(gfnorm) / gTBg;
405: }
406: auk = PetscMin(delta / gfnorm, auk);
407: PetscCall(VecCopy(GradF, YCtmp)); /* this could be improved */
408: PetscCall(VecScale(YCtmp, auk)); /* YCtmp, Cauchy solution*/
409: PetscCall(VecNorm(YCtmp, NORM_2, &ycnorm)); /* ycnorm <- || Y_cauchy || */
410: if (ycnorm >= delta) { /* see if the Cauchy solution meets the criteria */
411: PetscCall(VecCopy(YCtmp, Y));
412: PetscCall(PetscInfo(snes, "DL evaluated. delta: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n", (double)delta, (double)ynnorm, (double)ycnorm));
413: } else { /* take ratio, tau, of Cauchy and Newton direction and step */
414: PetscCall(VecAXPY(YNtmp, -1.0, YCtmp)); /* YCtmp = A, YNtmp = B */
415: PetscCall(VecNorm(YNtmp, NORM_2, &c0)); /* this could be improved */
416: c0 = PetscSqr(c0);
417: PetscCall(VecDotRealPart(YCtmp, YNtmp, &c1));
418: c1 = 2.0 * c1;
419: PetscCall(VecNorm(YCtmp, NORM_2, &c2)); /* this could be improved */
420: c2 = PetscSqr(c2) - PetscSqr(delta);
421: tau_pos = (c1 + PetscSqrtReal(PetscSqr(c1) - 4. * c0 * c2)) / (2. * c0); /* quadratic formula */
422: tau_neg = (c1 - PetscSqrtReal(PetscSqr(c1) - 4. * c0 * c2)) / (2. * c0);
423: tau = PetscMax(tau_pos, tau_neg); /* can tau_neg > tau_pos? I don't think so, but just in case. */
424: PetscCall(PetscInfo(snes, "DL evaluated. tau: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n", (double)tau, (double)ynnorm, (double)ycnorm));
425: PetscCall(VecWAXPY(W, tau, YNtmp, YCtmp));
426: PetscCall(VecAXPY(W, -tau, YCtmp));
427: PetscCall(VecCopy(W, Y)); /* this could be improved */
428: }
429: } else {
430: /* if Cauchy is disabled, only use Newton direction */
431: auk = delta / ynnorm;
432: PetscCall(VecScale(YNtmp, auk));
433: PetscCall(VecCopy(YNtmp, Y)); /* this could be improved (many VecCopy, VecNorm)*/
434: }
436: PetscCall(VecNorm(Y, NORM_2, &ynorm)); /* compute the final ynorm */
437: f0 = 0.5 * PetscSqr(fnorm); /* minimizing function f(X) */
438: PetscCall(MatMult(jac, Y, W));
439: PetscCall(VecDotRealPart(W, W, &yTHy)); /* completes GradY^T J^T J GradY */
440: PetscCall(VecDotRealPart(GradF, Y, &gTy));
441: mp = f0 - gTy + 0.5 * yTHy; /* quadratic model to satisfy, -gTy because our update is X-Y*/
443: /* scale back solution update */
444: if (bs > 1 && neP->auto_scale_multiphase) {
445: for (j = 0; j < bs; j++) {
446: PetscCall(VecStrideScale(Y, j, inorms[j]));
447: if (inner_count == 0) {
448: /* TRDC inner algorithm does not need scaled X after calculating delta in the outer iteration */
449: /* need to scale back X to match Y and provide proper update to the external code */
450: PetscCall(VecStrideScale(X, j, inorms[j]));
451: }
452: }
453: if (inner_count == 0) PetscCall(VecNorm(X, NORM_2, &temp_xnorm)); /* only in the first iteration */
454: PetscCall(VecNorm(Y, NORM_2, &temp_ynorm));
455: } else {
456: temp_xnorm = xnorm;
457: temp_ynorm = ynorm;
458: }
459: inner_count++;
461: /* Evaluate the solution to meet the improvement ratio criteria */
462: PetscCall(SNESNewtonTRDCPreCheck(snes, X, Y, &changed_y));
463: PetscCall(VecWAXPY(W, -1.0, Y, X));
464: PetscCall(SNESNewtonTRDCPostCheck(snes, X, Y, W, &changed_y, &changed_w));
465: if (changed_y) PetscCall(VecWAXPY(W, -1.0, Y, X));
466: PetscCall(VecCopy(Y, snes->vec_sol_update));
467: PetscCall(SNESComputeFunction(snes, W, G)); /* F(X-Y) = G */
468: PetscCall(VecNorm(G, NORM_2, &gnorm)); /* gnorm <- || g || */
469: SNESCheckFunctionNorm(snes, gnorm);
470: g = 0.5 * PetscSqr(gnorm); /* minimizing function g(W) */
471: if (f0 == mp) rho = 0.0;
472: else rho = (f0 - g) / (f0 - mp); /* actual improvement over predicted improvement */
474: if (rho < neP->eta2) {
475: delta *= neP->t1; /* shrink the region */
476: } else if (rho > neP->eta3) {
477: delta = PetscMin(neP->t2 * delta, deltaM); /* expand the region, but not greater than deltaM */
478: }
480: neP->delta = delta;
481: if (rho >= neP->eta1) {
482: /* unscale delta and xnorm before going to the next outer iteration */
483: if (bs > 1 && neP->auto_scale_multiphase) {
484: neP->delta = delta / xnorm;
485: xnorm = temp_xnorm;
486: ynorm = temp_ynorm;
487: }
488: neP->rho_satisfied = PETSC_TRUE;
489: break; /* the improvement ratio is satisfactory */
490: }
491: PetscCall(PetscInfo(snes, "Trying again in smaller region\n"));
493: /* check to see if progress is hopeless */
494: neP->itflag = PETSC_FALSE;
495: /* both delta, ynorm, and xnorm are either scaled or unscaled */
496: PetscCall(SNESTRDC_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP));
497: /* if multiphase state changes, break out inner iteration */
498: if (reason == SNES_BREAKOUT_INNER_ITER) {
499: if (bs > 1 && neP->auto_scale_multiphase) {
500: /* unscale delta and xnorm before going to the next outer iteration */
501: neP->delta = delta / xnorm;
502: xnorm = temp_xnorm;
503: ynorm = temp_ynorm;
504: }
505: reason = SNES_CONVERGED_ITERATING;
506: break;
507: }
508: if (reason == SNES_CONVERGED_SNORM_RELATIVE) reason = SNES_DIVERGED_INNER;
509: if (reason) {
510: if (reason < 0) {
511: /* We're not progressing, so return with the current iterate */
512: PetscCall(SNESMonitor(snes, i + 1, fnorm));
513: breakout = PETSC_TRUE;
514: break;
515: } else if (reason > 0) {
516: /* We're converged, so return with the current iterate and update solution */
517: PetscCall(SNESMonitor(snes, i + 1, fnorm));
518: breakout = PETSC_FALSE;
519: break;
520: }
521: }
522: snes->numFailures++;
523: }
524: if (!breakout) {
525: /* Update function and solution vectors */
526: fnorm = gnorm;
527: PetscCall(VecCopy(G, F));
528: PetscCall(VecCopy(W, X));
529: /* Monitor convergence */
530: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
531: snes->iter = i + 1;
532: snes->norm = fnorm;
533: snes->xnorm = xnorm;
534: snes->ynorm = ynorm;
535: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
536: PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
537: PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
538: /* Test for convergence, xnorm = || X || */
539: neP->itflag = PETSC_TRUE;
540: if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm));
541: PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP);
542: if (reason) break;
543: } else break;
544: }
546: /* PetscCall(PetscFree(inorms)); */
547: if (i == maxits) {
548: PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
549: if (!reason) reason = SNES_DIVERGED_MAX_IT;
550: }
551: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
552: snes->reason = reason;
553: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
554: if (convtest != SNESTRDC_KSPConverged_Private) {
555: PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy));
556: PetscCall(PetscFree(ctx));
557: PetscCall(KSPSetConvergenceTest(ksp, convtest, convctx, convdestroy));
558: }
559: PetscFunctionReturn(PETSC_SUCCESS);
560: }
562: static PetscErrorCode SNESSetUp_NEWTONTRDC(SNES snes)
563: {
564: PetscFunctionBegin;
565: PetscCall(SNESSetWorkVecs(snes, 6));
566: PetscCall(SNESSetUpMatrices(snes));
567: PetscFunctionReturn(PETSC_SUCCESS);
568: }
570: PetscErrorCode SNESReset_NEWTONTRDC(SNES snes)
571: {
572: PetscFunctionBegin;
573: PetscFunctionReturn(PETSC_SUCCESS);
574: }
576: static PetscErrorCode SNESDestroy_NEWTONTRDC(SNES snes)
577: {
578: PetscFunctionBegin;
579: PetscCall(SNESReset_NEWTONTRDC(snes));
580: PetscCall(PetscFree(snes->data));
581: PetscFunctionReturn(PETSC_SUCCESS);
582: }
584: static PetscErrorCode SNESSetFromOptions_NEWTONTRDC(SNES snes, PetscOptionItems *PetscOptionsObject)
585: {
586: SNES_NEWTONTRDC *ctx = (SNES_NEWTONTRDC *)snes->data;
588: PetscFunctionBegin;
589: PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations");
590: PetscCall(PetscOptionsReal("-snes_trdc_tol", "Trust region tolerance", "SNESSetTrustRegionTolerance", snes->deltatol, &snes->deltatol, NULL));
591: PetscCall(PetscOptionsReal("-snes_trdc_eta1", "eta1", "None", ctx->eta1, &ctx->eta1, NULL));
592: PetscCall(PetscOptionsReal("-snes_trdc_eta2", "eta2", "None", ctx->eta2, &ctx->eta2, NULL));
593: PetscCall(PetscOptionsReal("-snes_trdc_eta3", "eta3", "None", ctx->eta3, &ctx->eta3, NULL));
594: PetscCall(PetscOptionsReal("-snes_trdc_t1", "t1", "None", ctx->t1, &ctx->t1, NULL));
595: PetscCall(PetscOptionsReal("-snes_trdc_t2", "t2", "None", ctx->t2, &ctx->t2, NULL));
596: PetscCall(PetscOptionsReal("-snes_trdc_deltaM", "deltaM", "None", ctx->deltaM, &ctx->deltaM, NULL));
597: PetscCall(PetscOptionsReal("-snes_trdc_delta0", "delta0", "None", ctx->delta0, &ctx->delta0, NULL));
598: PetscCall(PetscOptionsReal("-snes_trdc_auto_scale_max", "auto_scale_max", "None", ctx->auto_scale_max, &ctx->auto_scale_max, NULL));
599: PetscCall(PetscOptionsBool("-snes_trdc_use_cauchy", "use_cauchy", "use Cauchy step and direction", ctx->use_cauchy, &ctx->use_cauchy, NULL));
600: PetscCall(PetscOptionsBool("-snes_trdc_auto_scale_multiphase", "auto_scale_multiphase", "Auto scaling for proper cauchy direction", ctx->auto_scale_multiphase, &ctx->auto_scale_multiphase, NULL));
601: PetscOptionsHeadEnd();
602: PetscFunctionReturn(PETSC_SUCCESS);
603: }
605: static PetscErrorCode SNESView_NEWTONTRDC(SNES snes, PetscViewer viewer)
606: {
607: SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
608: PetscBool iascii;
610: PetscFunctionBegin;
611: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
612: if (iascii) {
613: PetscCall(PetscViewerASCIIPrintf(viewer, " Trust region tolerance %g (-snes_trtol)\n", (double)snes->deltatol));
614: PetscCall(PetscViewerASCIIPrintf(viewer, " eta1=%g, eta2=%g, eta3=%g\n", (double)tr->eta1, (double)tr->eta2, (double)tr->eta3));
615: PetscCall(PetscViewerASCIIPrintf(viewer, " delta0=%g, t1=%g, t2=%g, deltaM=%g\n", (double)tr->delta0, (double)tr->t1, (double)tr->t2, (double)tr->deltaM));
616: }
617: PetscFunctionReturn(PETSC_SUCCESS);
618: }
620: /*MC
621: SNESNEWTONTRDC - Newton based nonlinear solver that uses trust-region dogleg method with Cauchy direction
623: Options Database Keys:
624: + -snes_trdc_tol <tol> - trust region tolerance
625: . -snes_trdc_eta1 <eta1> - trust region parameter 0.0 <= eta1 <= eta2, rho >= eta1 breaks out of the inner iteration (default: eta1=0.001)
626: . -snes_trdc_eta2 <eta2> - trust region parameter 0.0 <= eta1 <= eta2, rho <= eta2 shrinks the trust region (default: eta2=0.25)
627: . -snes_trdc_eta3 <eta3> - trust region parameter eta3 > eta2, rho >= eta3 expands the trust region (default: eta3=0.75)
628: . -snes_trdc_t1 <t1> - trust region parameter, shrinking factor of trust region (default: 0.25)
629: . -snes_trdc_t2 <t2> - trust region parameter, expanding factor of trust region (default: 2.0)
630: . -snes_trdc_deltaM <deltaM> - trust region parameter, max size of trust region, deltaM*norm2(x) (default: 0.5)
631: . -snes_trdc_delta0 <delta0> - trust region parameter, initial size of trust region, delta0*norm2(x) (default: 0.1)
632: . -snes_trdc_auto_scale_max <auto_scale_max> - used with auto_scale_multiphase, caps the maximum auto-scaling factor
633: . -snes_trdc_use_cauchy <use_cauchy> - True uses dogleg Cauchy (Steepest Descent direction) step & direction in the trust region algorithm
634: - -snes_trdc_auto_scale_multiphase <auto_scale_multiphase> - True turns on auto-scaling for multivariable block matrix for Cauchy and trust region
636: Reference:
637: . - * "Linear and Nonlinear Solvers for Simulating Multiphase Flow
638: within Large-Scale Engineered Subsurface Systems" by Heeho D. Park, Glenn E. Hammond,
639: Albert J. Valocchi, Tara LaForce.
641: Level: intermediate
643: .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESSetTrustRegionTolerance()`,
644: `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`,
645: `SNESNewtonTRDCGetRhoFlag()`, `SNESNewtonTRDCSetPreCheck()`
646: M*/
647: PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTRDC(SNES snes)
648: {
649: SNES_NEWTONTRDC *neP;
651: PetscFunctionBegin;
652: snes->ops->setup = SNESSetUp_NEWTONTRDC;
653: snes->ops->solve = SNESSolve_NEWTONTRDC;
654: snes->ops->destroy = SNESDestroy_NEWTONTRDC;
655: snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTRDC;
656: snes->ops->view = SNESView_NEWTONTRDC;
657: snes->ops->reset = SNESReset_NEWTONTRDC;
659: snes->usesksp = PETSC_TRUE;
660: snes->usesnpc = PETSC_FALSE;
662: snes->alwayscomputesfinalresidual = PETSC_TRUE;
664: PetscCall(PetscNew(&neP));
665: snes->data = (void *)neP;
666: neP->delta = 0.0;
667: neP->delta0 = 0.1;
668: neP->eta1 = 0.001;
669: neP->eta2 = 0.25;
670: neP->eta3 = 0.75;
671: neP->t1 = 0.25;
672: neP->t2 = 2.0;
673: neP->deltaM = 0.5;
674: neP->sigma = 0.0001;
675: neP->itflag = PETSC_FALSE;
676: neP->rnorm0 = 0.0;
677: neP->ttol = 0.0;
678: neP->use_cauchy = PETSC_TRUE;
679: neP->auto_scale_multiphase = PETSC_FALSE;
680: neP->auto_scale_max = -1.0;
681: neP->rho_satisfied = PETSC_FALSE;
682: snes->deltatol = 1.e-12;
684: /* for multiphase (multivariable) scaling */
685: /* may be used for dynamic allocation of inorms, but it fails snes_tutorials-ex3_13
686: on test forced DIVERGED_JACOBIAN_DOMAIN test. I will use static array for now.
687: PetscCall(VecGetBlockSize(snes->work[0],&neP->bs));
688: PetscCall(PetscCalloc1(neP->bs,&neP->inorms));
689: */
691: PetscFunctionReturn(PETSC_SUCCESS);
692: }