Actual source code: tr.c
1: #include <../src/snes/impls/tr/trimpl.h>
3: typedef struct {
4: SNES snes;
5: KSPConvergenceTestFn *convtest;
6: PetscCtxDestroyFn *convdestroy;
7: void *convctx;
8: } SNES_TR_KSPConverged_Ctx;
10: const char *const SNESNewtonTRFallbackTypes[] = {"NEWTON", "CAUCHY", "DOGLEG", "SNESNewtonTRFallbackType", "SNES_TR_FALLBACK_", NULL};
11: const char *const SNESNewtonTRQNTypes[] = {"NONE", "SAME", "DIFFERENT", "SNESNewtonTRQNType", "SNES_TR_QN_", NULL};
13: static PetscErrorCode SNESNewtonTRSetTolerances_TR(SNES snes, PetscReal delta_min, PetscReal delta_max, PetscReal delta_0)
14: {
15: SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
17: PetscFunctionBegin;
18: if (delta_min == PETSC_DETERMINE) delta_min = tr->default_deltam;
19: if (delta_max == PETSC_DETERMINE) delta_max = tr->default_deltaM;
20: if (delta_0 == PETSC_DETERMINE) delta_0 = tr->default_delta0;
21: if (delta_min != PETSC_CURRENT) tr->deltam = delta_min;
22: if (delta_max != PETSC_CURRENT) tr->deltaM = delta_max;
23: if (delta_0 != PETSC_CURRENT) tr->delta0 = delta_0;
24: PetscFunctionReturn(PETSC_SUCCESS);
25: }
27: static PetscErrorCode SNESNewtonTRGetTolerances_TR(SNES snes, PetscReal *delta_min, PetscReal *delta_max, PetscReal *delta_0)
28: {
29: SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
31: PetscFunctionBegin;
32: if (delta_min) *delta_min = tr->deltam;
33: if (delta_max) *delta_max = tr->deltaM;
34: if (delta_0) *delta_0 = tr->delta0;
35: PetscFunctionReturn(PETSC_SUCCESS);
36: }
38: static PetscErrorCode SNESComputeJacobian_MATLMVM(SNES snes, Vec X, Mat J, Mat B, void *dummy)
39: {
40: PetscFunctionBegin;
41: // PetscCall(MatLMVMSymBroydenSetDelta(B, _some_delta));
42: PetscCall(MatLMVMUpdate(B, X, snes->vec_func));
43: PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
44: PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
45: if (J != B) {
46: // PetscCall(MatLMVMSymBroydenSetDelta(J, _some_delta));
47: PetscCall(MatLMVMUpdate(J, X, snes->vec_func));
48: PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
49: PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
50: }
51: PetscFunctionReturn(PETSC_SUCCESS);
52: }
54: static PetscErrorCode SNESTR_KSPConverged_Private(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *cctx)
55: {
56: SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx *)cctx;
57: SNES snes = ctx->snes;
58: SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data;
59: Vec x;
60: PetscReal nrm;
62: PetscFunctionBegin;
63: /* Determine norm of solution */
64: PetscCall(KSPBuildSolution(ksp, NULL, &x));
65: PetscCall(VecNorm(x, neP->norm, &nrm));
66: if (nrm >= neP->delta) {
67: PetscCall(PetscInfo(snes, "Ending linear iteration early due to exiting trust region, delta=%g, length=%g\n", (double)neP->delta, (double)nrm));
68: *reason = KSP_CONVERGED_STEP_LENGTH;
69: PetscFunctionReturn(PETSC_SUCCESS);
70: }
71: PetscCall((*ctx->convtest)(ksp, n, rnorm, reason, ctx->convctx));
72: if (*reason) PetscCall(PetscInfo(snes, "Default or user provided convergence test KSP iterations=%" PetscInt_FMT ", rnorm=%g\n", n, (double)rnorm));
73: PetscFunctionReturn(PETSC_SUCCESS);
74: }
76: static PetscErrorCode SNESTR_KSPConverged_Destroy(PetscCtxRt cctx)
77: {
78: SNES_TR_KSPConverged_Ctx *ctx = *(SNES_TR_KSPConverged_Ctx **)cctx;
80: PetscFunctionBegin;
81: PetscCall((*ctx->convdestroy)(&ctx->convctx));
82: PetscCall(PetscFree(ctx));
83: PetscFunctionReturn(PETSC_SUCCESS);
84: }
86: static PetscErrorCode SNESTR_Converged_Private(SNES snes, PetscInt it, PetscReal xnorm, PetscReal pnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy)
87: {
88: SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data;
90: PetscFunctionBegin;
91: *reason = SNES_CONVERGED_ITERATING;
92: if (neP->delta < neP->deltam) {
93: PetscCall(PetscInfo(snes, "Diverged due to too small a trust region %g<%g\n", (double)neP->delta, (double)neP->deltam));
94: *reason = SNES_DIVERGED_TR_DELTA;
95: } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
96: PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n", snes->max_funcs));
97: *reason = SNES_DIVERGED_FUNCTION_COUNT;
98: }
99: PetscFunctionReturn(PETSC_SUCCESS);
100: }
102: /*@
103: SNESNewtonTRSetNormType - Specify the type of norm to use for the computation of the trust region.
105: Input Parameters:
106: + snes - the nonlinear solver object
107: - norm - the norm type
109: Level: intermediate
111: .seealso: `SNESNEWTONTR`, `NormType`
112: @*/
113: PetscErrorCode SNESNewtonTRSetNormType(SNES snes, NormType norm)
114: {
115: PetscBool flg;
117: PetscFunctionBegin;
120: PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
121: if (flg) {
122: SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
124: tr->norm = norm;
125: }
126: PetscFunctionReturn(PETSC_SUCCESS);
127: }
129: /*@
130: SNESNewtonTRSetQNType - Specify to use a quasi-Newton model.
132: Input Parameters:
133: + snes - the nonlinear solver object
134: - use - the type of approximations to be used
136: Level: intermediate
138: Notes:
139: Options for the approximations can be set with the `snes_tr_qn_` and `snes_tr_qn_pre_` prefixes.
141: .seealso: `SNESNEWTONTR`, `SNESNewtonTRQNType`, `MATLMVM`
142: @*/
143: PetscErrorCode SNESNewtonTRSetQNType(SNES snes, SNESNewtonTRQNType use)
144: {
145: PetscBool flg;
147: PetscFunctionBegin;
150: PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
151: if (flg) {
152: SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
154: tr->qn = use;
155: }
156: PetscFunctionReturn(PETSC_SUCCESS);
157: }
159: /*@
160: SNESNewtonTRSetFallbackType - Set the type of fallback to use if the solution of the trust region subproblem is outside the radius
162: Input Parameters:
163: + snes - the nonlinear solver object
164: - ftype - the fallback type, see `SNESNewtonTRFallbackType`
166: Level: intermediate
168: .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPreCheck()`,
169: `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`
170: @*/
171: PetscErrorCode SNESNewtonTRSetFallbackType(SNES snes, SNESNewtonTRFallbackType ftype)
172: {
173: SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
174: PetscBool flg;
176: PetscFunctionBegin;
179: PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
180: if (flg) tr->fallback = ftype;
181: PetscFunctionReturn(PETSC_SUCCESS);
182: }
184: /*@C
185: SNESNewtonTRSetPreCheck - Sets a user function that is called before the search step has been determined.
186: Allows the user a chance to change or override the trust region decision.
188: Logically Collective
190: Input Parameters:
191: + snes - the nonlinear solver object
192: . func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPreCheck()`
193: - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
195: Calling sequence of `func`:
196: + snes - the nonlinear solver object
197: . X - the current solution value
198: . Y - the tentative update step
199: . changed - output, flag indicating `Y` has been changed by the pre-check
200: - ctx - the optional application context
202: Level: intermediate
204: Note:
205: This function is called BEFORE the function evaluation within the solver.
207: .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`,
208: @*/
209: PetscErrorCode SNESNewtonTRSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES snes, Vec X, Vec Y, PetscBool *changed, PetscCtx ctx), PetscCtx ctx)
210: {
211: SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
212: PetscBool flg;
214: PetscFunctionBegin;
216: PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
217: if (flg) {
218: if (func) tr->precheck = func;
219: if (ctx) tr->precheckctx = ctx;
220: }
221: PetscFunctionReturn(PETSC_SUCCESS);
222: }
224: /*@C
225: SNESNewtonTRGetPreCheck - Gets the pre-check function
227: Not Collective
229: Input Parameter:
230: . snes - the nonlinear solver context
232: Output Parameters:
233: + func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPreCheck()`
234: - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
236: Calling sequence of `func`:
237: + snes - the nonlinear solver object
238: . X - the current solution value
239: . Y - the tentative update step
240: . changed - output, flag indicating `Y` has been changed by the pre-check
241: - ctx - the optional application context
243: Level: intermediate
245: .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRPreCheck()`
246: @*/
247: PetscErrorCode SNESNewtonTRGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES snes, Vec X, Vec Y, PetscBool *changed, PetscCtx ctx), PetscCtxRt ctx)
248: {
249: SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
250: PetscBool flg;
252: PetscFunctionBegin;
254: PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
255: PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name);
256: if (func) *func = tr->precheck;
257: if (ctx) *(void **)ctx = tr->precheckctx;
258: PetscFunctionReturn(PETSC_SUCCESS);
259: }
261: /*@C
262: SNESNewtonTRSetPostCheck - Sets a user function that is called after the search step has been determined but before the next
263: function evaluation. Allows the user a chance to change or override the internal decision of the solver
265: Logically Collective
267: Input Parameters:
268: + snes - the nonlinear solver object
269: . func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPostCheck()`
270: - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
272: Calling sequence of `func`:
273: + snes - the nonlinear solver object
274: . X - the current solution value
275: . Y - the tentative update step
276: . W - the tentative new solution value
277: . changed_Y - output, flag indicated `Y` has been changed by the post-check
278: . changed_W - output, flag indicated `W` has been changed by the post-check
279: - ctx - the optional application context
281: Level: intermediate
283: Note:
284: This function is called BEFORE the function evaluation within the solver while the function set in
285: `SNESLineSearchSetPostCheck()` is called AFTER the function evaluation.
287: .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRPostCheck()`, `SNESNewtonTRGetPostCheck()`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()`
288: @*/
289: PetscErrorCode SNESNewtonTRSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W, PetscCtx ctx), PetscCtx ctx)
290: {
291: SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
292: PetscBool flg;
294: PetscFunctionBegin;
296: PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
297: if (flg) {
298: if (func) tr->postcheck = func;
299: if (ctx) tr->postcheckctx = ctx;
300: }
301: PetscFunctionReturn(PETSC_SUCCESS);
302: }
304: /*@C
305: SNESNewtonTRGetPostCheck - Gets the post-check function
307: Not Collective
309: Input Parameter:
310: . snes - the nonlinear solver context
312: Output Parameters:
313: + func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPostCheck()`
314: - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
316: Calling sequence of `func`:
317: + snes - the nonlinear solver object
318: . X - the current solution value
319: . Y - the tentative update step
320: . W - the tentative new solution value
321: . changed_Y - output, flag indicated `Y` has been changed by the post-check
322: . changed_W - output, flag indicated `W` has been changed by the post-check
323: - ctx - the optional application context
325: Level: intermediate
327: .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRPostCheck()`
328: @*/
329: PetscErrorCode SNESNewtonTRGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W, PetscCtx ctx), PetscCtxRt ctx)
330: {
331: SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
332: PetscBool flg;
334: PetscFunctionBegin;
336: PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
337: PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name);
338: if (func) *func = tr->postcheck;
339: if (ctx) *(void **)ctx = tr->postcheckctx;
340: PetscFunctionReturn(PETSC_SUCCESS);
341: }
343: /*@C
344: SNESNewtonTRPreCheck - Runs the precheck routine
346: Logically Collective
348: Input Parameters:
349: + snes - the solver
350: . X - The last solution
351: - Y - The step direction
353: Output Parameter:
354: . changed_Y - Indicator that the step direction `Y` has been changed.
356: Level: intermediate
358: .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRPostCheck()`
359: @*/
360: PetscErrorCode SNESNewtonTRPreCheck(SNES snes, Vec X, Vec Y, PetscBool *changed_Y)
361: {
362: SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
363: PetscBool flg;
365: PetscFunctionBegin;
367: PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
368: PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name);
369: *changed_Y = PETSC_FALSE;
370: if (tr->precheck) {
371: PetscCall((*tr->precheck)(snes, X, Y, changed_Y, tr->precheckctx));
373: }
374: PetscFunctionReturn(PETSC_SUCCESS);
375: }
377: /*@C
378: SNESNewtonTRPostCheck - Runs the postcheck routine
380: Logically Collective
382: Input Parameters:
383: + snes - the solver
384: . X - The last solution
385: . Y - The full step direction
386: - W - The updated solution, W = X - Y
388: Output Parameters:
389: + changed_Y - indicator if step has been changed
390: - changed_W - Indicator if the new candidate solution W has been changed.
392: Note:
393: If Y is changed then W is recomputed as X - Y
395: Level: intermediate
397: .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`, `SNESNewtonTRPreCheck()`
398: @*/
399: PetscErrorCode SNESNewtonTRPostCheck(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W)
400: {
401: SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
402: PetscBool flg;
404: PetscFunctionBegin;
406: PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
407: PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name);
408: *changed_Y = PETSC_FALSE;
409: *changed_W = PETSC_FALSE;
410: if (tr->postcheck) {
411: PetscCall((*tr->postcheck)(snes, X, Y, W, changed_Y, changed_W, tr->postcheckctx));
414: }
415: PetscFunctionReturn(PETSC_SUCCESS);
416: }
418: /* stable implementation of roots of a*x^2 + b*x + c = 0 */
419: static inline void PetscQuadraticRoots(PetscReal a, PetscReal b, PetscReal c, PetscReal *xm, PetscReal *xp)
420: {
421: PetscReal temp = -0.5 * (b + PetscCopysignReal(1.0, b) * PetscSqrtReal(b * b - 4 * a * c));
422: PetscReal x1 = temp / a;
423: PetscReal x2 = c / temp;
424: *xm = PetscMin(x1, x2);
425: *xp = PetscMax(x1, x2);
426: }
428: /* Computes the quadratic model difference */
429: static PetscErrorCode SNESNewtonTRQuadraticDelta(SNES snes, Mat J, PetscBool has_objective, Vec Y, Vec GradF, Vec W, PetscReal *yTHy_, PetscReal *gTy_, PetscReal *deltaqm)
430: {
431: PetscReal yTHy, gTy;
433: PetscFunctionBegin;
434: PetscCall(MatMult(J, Y, W));
435: if (has_objective) PetscCall(VecDotRealPart(Y, W, &yTHy));
436: else PetscCall(VecDotRealPart(W, W, &yTHy)); /* Gauss-Newton approximation J^t * J */
437: PetscCall(VecDotRealPart(GradF, Y, &gTy));
438: *deltaqm = -(-(gTy) + 0.5 * (yTHy)); /* difference in quadratic model, -gTy because SNES solves it this way */
439: if (yTHy_) *yTHy_ = yTHy;
440: if (gTy_) *gTy_ = gTy;
441: PetscFunctionReturn(PETSC_SUCCESS);
442: }
444: /* Computes the new objective given X = Xk, Y = direction
445: W work vector, on output W = X - Y
446: G work vector, on output G = SNESFunction(W) */
447: static PetscErrorCode SNESNewtonTRObjective(SNES snes, PetscBool has_objective, Vec X, Vec Y, Vec W, Vec G, PetscReal *gnorm, PetscReal *fkp1)
448: {
449: PetscBool changed_Y, changed_W;
451: PetscFunctionBegin;
452: /* TODO: we can add a linesearch here */
453: PetscCall(SNESNewtonTRPreCheck(snes, X, Y, &changed_Y));
454: PetscCall(VecWAXPY(W, -1.0, Y, X)); /* Xkp1 */
455: PetscCall(SNESNewtonTRPostCheck(snes, X, Y, W, &changed_Y, &changed_W));
456: if (changed_Y && !changed_W) PetscCall(VecWAXPY(W, -1.0, Y, X));
458: PetscCall(SNESComputeFunction(snes, W, G)); /* F(Xkp1) = G */
459: PetscCall(VecNorm(G, NORM_2, gnorm));
460: SNESCheckFunctionDomainError(snes, *gnorm);
461: if (has_objective) {
462: PetscCall(SNESComputeObjective(snes, W, fkp1));
463: SNESCheckObjectiveDomainError(snes, *fkp1);
464: } else *fkp1 = 0.5 * PetscSqr(*gnorm);
465: PetscFunctionReturn(PETSC_SUCCESS);
466: }
468: static PetscErrorCode SNESSetUpQN_NEWTONTR(SNES snes)
469: {
470: SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
472: PetscFunctionBegin;
473: PetscCall(MatDestroy(&tr->qnB));
474: PetscCall(MatDestroy(&tr->qnB_pre));
475: if (tr->qn) {
476: PetscInt n, N;
477: const char *optionsprefix;
478: Mat B;
480: PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &B));
481: PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
482: PetscCall(MatSetOptionsPrefix(B, "snes_tr_qn_"));
483: PetscCall(MatAppendOptionsPrefix(B, optionsprefix));
484: PetscCall(MatSetType(B, MATLMVMBFGS));
485: PetscCall(VecGetLocalSize(snes->vec_sol, &n));
486: PetscCall(VecGetSize(snes->vec_sol, &N));
487: PetscCall(MatSetSizes(B, n, n, N, N));
488: PetscCall(MatSetUp(B));
489: PetscCall(MatSetFromOptions(B));
490: PetscCall(MatLMVMAllocate(B, snes->vec_sol, snes->vec_func));
491: tr->qnB = B;
492: if (tr->qn == SNES_TR_QN_DIFFERENT) {
493: PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &B));
494: PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
495: PetscCall(MatSetOptionsPrefix(B, "snes_tr_qn_pre_"));
496: PetscCall(MatAppendOptionsPrefix(B, optionsprefix));
497: PetscCall(MatSetType(B, MATLMVMBFGS));
498: PetscCall(MatSetSizes(B, n, n, N, N));
499: PetscCall(MatSetUp(B));
500: PetscCall(MatSetFromOptions(B));
501: PetscCall(MatLMVMAllocate(B, snes->vec_sol, snes->vec_func));
502: tr->qnB_pre = B;
503: } else {
504: PetscCall(PetscObjectReference((PetscObject)tr->qnB));
505: tr->qnB_pre = tr->qnB;
506: }
507: }
508: PetscFunctionReturn(PETSC_SUCCESS);
509: }
511: /*
512: SNESSolve_NEWTONTR - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy
513: (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of
514: nonlinear equations
516: */
517: static PetscErrorCode SNESSolve_NEWTONTR(SNES snes)
518: {
519: SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data;
520: Vec X, F, Y, G, W, GradF, YU, Yc;
521: PetscInt maxits, lits;
522: PetscReal rho, fnorm, gnorm = 0.0, xnorm = 0.0, delta, ynorm;
523: PetscReal fk, fkp1, deltaqm = 0.0, gTy = 0.0, yTHy = 0.0;
524: PetscReal auk, tauk, gfnorm, gfnorm_k, ycnorm, gTBg, objmin = 0.0, beta_k = 1.0;
525: PC pc;
526: Mat J, Jp;
527: PetscBool already_done = PETSC_FALSE, on_boundary, use_cauchy;
528: PetscBool clear_converged_test, rho_satisfied, has_objective;
529: SNES_TR_KSPConverged_Ctx *ctx;
530: void *convctx;
531: SNESObjectiveFn *objective;
532: KSPConvergenceTestFn *convtest;
533: PetscCtxDestroyFn *convdestroy;
535: PetscFunctionBegin;
536: PetscCall(SNESGetObjective(snes, &objective, NULL));
537: has_objective = objective ? PETSC_TRUE : PETSC_FALSE;
539: maxits = snes->max_its; /* maximum number of iterations */
540: X = snes->vec_sol; /* solution vector */
541: F = snes->vec_func; /* residual vector */
542: Y = snes->vec_sol_update; /* update vector */
543: G = snes->work[0]; /* updated residual */
544: W = snes->work[1]; /* temporary vector */
545: GradF = !has_objective ? snes->work[2] : snes->vec_func; /* grad f = J^T F */
546: YU = snes->work[3]; /* work vector for dogleg method */
547: Yc = snes->work[4]; /* Cauchy point */
549: 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);
551: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
552: snes->iter = 0;
553: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
555: /* setup QN matrices if needed */
556: PetscCall(SNESSetUpQN_NEWTONTR(snes));
558: /* Set the linear stopping criteria to use the More' trick if needed */
559: clear_converged_test = PETSC_FALSE;
560: PetscCall(SNESGetKSP(snes, &snes->ksp));
561: PetscCall(KSPGetConvergenceTest(snes->ksp, &convtest, &convctx, &convdestroy));
562: if (convtest != SNESTR_KSPConverged_Private) {
563: clear_converged_test = PETSC_TRUE;
564: PetscCall(PetscNew(&ctx));
565: ctx->snes = snes;
566: PetscCall(KSPGetAndClearConvergenceTest(snes->ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy));
567: PetscCall(KSPSetConvergenceTest(snes->ksp, SNESTR_KSPConverged_Private, ctx, SNESTR_KSPConverged_Destroy));
568: PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTR_KSPConverged_Private\n"));
569: }
571: if (!snes->vec_func_init_set) {
572: PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */
573: } else snes->vec_func_init_set = PETSC_FALSE;
575: PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */
576: SNESCheckFunctionDomainError(snes, fnorm);
577: PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */
579: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
580: snes->norm = fnorm;
581: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
582: delta = neP->delta0;
583: neP->delta = delta;
584: PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
586: /* test convergence */
587: rho_satisfied = PETSC_FALSE;
588: PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
589: PetscCall(SNESMonitor(snes, 0, fnorm));
590: if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
592: if (has_objective) {
593: PetscCall(SNESComputeObjective(snes, X, &fk));
594: SNESCheckObjectiveDomainError(snes, fk);
595: } else fk = 0.5 * PetscSqr(fnorm); /* obj(x) = 0.5 * ||F(x)||^2 */
597: /* hook state vector to BFGS preconditioner */
598: PetscCall(KSPGetPC(snes->ksp, &pc));
599: PetscCall(PCLMVMSetUpdateVec(pc, X));
601: if (neP->kmdc) PetscCall(KSPSetComputeEigenvalues(snes->ksp, PETSC_TRUE));
603: while (snes->iter < maxits) {
604: /* calculating Jacobian and GradF of minimization function only once */
605: if (!already_done) {
606: /* Call general purpose update function */
607: PetscTryTypeMethod(snes, update, snes->iter);
609: /* apply the nonlinear preconditioner */
610: if (snes->npc && snes->npcside == PC_RIGHT) {
611: SNESConvergedReason reason;
613: PetscCall(SNESSetInitialFunction(snes->npc, F));
614: PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0));
615: PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X));
616: PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0));
617: PetscCall(SNESGetConvergedReason(snes->npc, &reason));
618: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) {
619: snes->reason = SNES_DIVERGED_INNER;
620: PetscFunctionReturn(PETSC_SUCCESS);
621: }
622: // XXX
623: PetscCall(SNESGetNPCFunction(snes, F, &fnorm));
624: }
626: /* Jacobian */
627: J = NULL;
628: Jp = NULL;
629: if (!neP->qnB) {
630: PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
631: SNESCheckJacobianDomainError(snes);
632: J = snes->jacobian;
633: Jp = snes->jacobian_pre;
634: } else { /* QN model */
635: PetscCall(SNESComputeJacobian_MATLMVM(snes, X, neP->qnB, neP->qnB_pre, NULL));
636: J = neP->qnB;
637: Jp = neP->qnB_pre;
638: }
639: SNESCheckJacobianDomainError(snes);
641: /* objective function */
642: PetscCall(VecNorm(F, NORM_2, &fnorm));
643: SNESCheckFunctionDomainError(snes, fnorm);
644: if (has_objective) {
645: PetscCall(SNESComputeObjective(snes, X, &fk));
646: SNESCheckObjectiveDomainError(snes, fk);
647: } else fk = 0.5 * PetscSqr(fnorm); /* obj(x) = 0.5 * ||F(x)||^2 */
649: /* GradF */
650: if (has_objective) gfnorm = fnorm;
651: else {
652: PetscCall(MatMultTranspose(J, F, GradF)); /* grad f = J^T F */
653: PetscCall(VecNorm(GradF, NORM_2, &gfnorm));
654: }
655: PetscCall(VecNorm(GradF, neP->norm, &gfnorm_k));
656: }
657: already_done = PETSC_TRUE;
659: /* solve trust-region subproblem */
661: /* first compute Cauchy Point */
662: PetscCall(MatMult(J, GradF, W));
663: if (has_objective) PetscCall(VecDotRealPart(GradF, W, &gTBg));
664: else PetscCall(VecDotRealPart(W, W, &gTBg)); /* B = J^t * J */
665: /* Eqs 4.11 and 4.12 in Nocedal and Wright 2nd Edition (4.7 and 4.8 in 1st Edition) */
666: auk = delta / gfnorm_k;
667: if (gTBg < 0.0) tauk = 1.0;
668: else tauk = PetscMin(gfnorm * gfnorm * gfnorm_k / (delta * gTBg), 1);
669: auk *= tauk;
670: ycnorm = auk * gfnorm;
671: PetscCall(VecAXPBY(Yc, auk, 0.0, GradF));
673: on_boundary = PETSC_FALSE;
674: use_cauchy = (PetscBool)(tauk == 1.0 && has_objective);
675: if (!use_cauchy) {
676: KSPConvergedReason reason;
678: /* sufficient decrease (see 6.3.27 in Conn, Gould, Toint "Trust Region Methods")
679: beta_k the largest eigenvalue of the Hessian. Here we use the previous estimated value */
680: objmin = -neP->kmdc * gnorm * PetscMin(gnorm / beta_k, delta);
681: PetscCall(KSPCGSetObjectiveTarget(snes->ksp, objmin));
683: /* specify radius if looking for Newton step and trust region norm is the l2 norm */
684: PetscCall(KSPCGSetRadius(snes->ksp, neP->fallback == SNES_TR_FALLBACK_NEWTON && neP->norm == NORM_2 ? delta : 0.0));
685: PetscCall(KSPSetOperators(snes->ksp, J, Jp));
686: PetscCall(KSPSolve(snes->ksp, F, Y));
687: SNESCheckKSPSolve(snes);
688: PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
689: PetscCall(KSPGetConvergedReason(snes->ksp, &reason));
690: on_boundary = (PetscBool)(reason == KSP_CONVERGED_STEP_LENGTH);
691: PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits));
692: if (neP->kmdc) { /* update estimated Hessian largest eigenvalue */
693: PetscReal emax, emin;
694: PetscCall(KSPComputeExtremeSingularValues(snes->ksp, &emax, &emin));
695: if (emax > 0.0) beta_k = emax + 1;
696: }
697: } else { /* Cauchy point is on the boundary, accept it */
698: on_boundary = PETSC_TRUE;
699: PetscCall(VecCopy(Yc, Y));
700: PetscCall(PetscInfo(snes, "CP evaluated on boundary. delta: %g, ycnorm: %g, gTBg: %g\n", (double)delta, (double)ycnorm, (double)gTBg));
701: }
702: PetscCall(VecNorm(Y, neP->norm, &ynorm));
704: /* decide what to do when the update is outside of trust region */
705: if (!use_cauchy && (ynorm > delta || ynorm == 0.0)) {
706: SNESNewtonTRFallbackType fallback = ynorm > 0.0 ? neP->fallback : SNES_TR_FALLBACK_CAUCHY;
708: PetscCheck(neP->norm == NORM_2 || fallback != SNES_TR_FALLBACK_DOGLEG, PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "DOGLEG without l2 norm not implemented");
709: switch (fallback) {
710: case SNES_TR_FALLBACK_NEWTON:
711: auk = delta / ynorm;
712: PetscCall(VecScale(Y, auk));
713: PetscCall(PetscInfo(snes, "SN evaluated. delta: %g, ynorm: %g\n", (double)delta, (double)ynorm));
714: break;
715: case SNES_TR_FALLBACK_CAUCHY:
716: case SNES_TR_FALLBACK_DOGLEG:
717: if (fallback == SNES_TR_FALLBACK_CAUCHY || gTBg <= 0.0) {
718: PetscCall(VecCopy(Yc, Y));
719: PetscCall(PetscInfo(snes, "CP evaluated. delta: %g, ynorm: %g, ycnorm: %g, gTBg: %g\n", (double)delta, (double)ynorm, (double)ycnorm, (double)gTBg));
720: } else { /* take linear combination of Cauchy and Newton direction and step */
721: auk = gfnorm * gfnorm / gTBg;
722: if (gfnorm_k * auk >= delta) { /* first leg: Cauchy point outside of trust region */
723: PetscCall(VecAXPBY(Y, delta / gfnorm_k, 0.0, GradF));
724: PetscCall(PetscInfo(snes, "CP evaluated (outside region). delta: %g, ynorm: %g, ycnorm: %g\n", (double)delta, (double)ynorm, (double)ycnorm));
725: } else { /* second leg */
726: PetscReal c0, c1, c2, tau = 0.0, tpos, tneg;
727: PetscBool noroots;
729: /* Find solutions of (Eq. 4.16 in Nocedal and Wright)
730: ||p_U + lambda * (p_B - p_U)||^2 - delta^2 = 0,
731: where p_U the Cauchy direction, p_B the Newton direction */
732: PetscCall(VecAXPBY(YU, auk, 0.0, GradF));
733: PetscCall(VecAXPY(Y, -1.0, YU));
734: PetscCall(VecNorm(Y, NORM_2, &c0));
735: PetscCall(VecDotRealPart(YU, Y, &c1));
736: c0 = PetscSqr(c0);
737: c2 = PetscSqr(ycnorm) - PetscSqr(delta);
738: PetscQuadraticRoots(c0, 2 * c1, c2, &tneg, &tpos);
740: /* In principle the DL strategy as a unique solution in [0,1]
741: here we check that for some reason we numerically failed
742: to compute it. In that case, we use the Cauchy point */
743: noroots = PetscIsInfOrNanReal(tneg);
744: if (!noroots) {
745: if (tpos > 1) {
746: if (tneg >= 0 && tneg <= 1) {
747: tau = tneg;
748: } else noroots = PETSC_TRUE;
749: } else if (tpos >= 0) {
750: tau = tpos;
751: } else noroots = PETSC_TRUE;
752: }
753: if (noroots) { /* No roots, select Cauchy point */
754: PetscCall(VecCopy(Yc, Y));
755: } else {
756: PetscCall(VecAXPBY(Y, 1.0, tau, YU));
757: }
758: PetscCall(PetscInfo(snes, "%s evaluated. roots: (%g, %g), tau %g, ynorm: %g, ycnorm: %g, gTBg: %g\n", noroots ? "CP" : "DL", (double)tneg, (double)tpos, (double)tau, (double)ynorm, (double)ycnorm, (double)gTBg));
759: }
760: }
761: break;
762: default:
763: SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "Unknown fallback mode");
764: break;
765: }
766: }
768: /* compute the quadratic model difference */
769: PetscCall(SNESNewtonTRQuadraticDelta(snes, J, has_objective, Y, GradF, W, &yTHy, &gTy, &deltaqm));
771: /* Compute new objective function */
772: PetscCall(SNESNewtonTRObjective(snes, has_objective, X, Y, W, G, &gnorm, &fkp1));
773: if (PetscIsInfOrNanReal(fkp1)) rho = neP->eta1;
774: else {
775: if (deltaqm > 0.0) rho = (fk - fkp1) / deltaqm; /* actual improvement over predicted improvement */
776: else rho = neP->eta1; /* no reduction in quadratic model, step must be rejected */
777: }
779: PetscCall(VecNorm(Y, neP->norm, &ynorm));
780: PetscCall(PetscInfo(snes, "rho=%g, delta=%g, fk=%g, fkp1=%g, deltaqm=%g, gTy=%g, yTHy=%g, ynormk=%g\n", (double)rho, (double)delta, (double)fk, (double)fkp1, (double)deltaqm, (double)gTy, (double)yTHy, (double)ynorm));
782: /* update the size of the trust region */
783: if (rho < neP->eta2) delta *= neP->t1; /* shrink the region */
784: else if (rho > neP->eta3 && on_boundary) delta *= neP->t2; /* expand the region */
785: delta = PetscMin(delta, neP->deltaM); /* but not greater than deltaM */
787: /* log 2-norm of update for moniroting routines */
788: PetscCall(VecNorm(Y, NORM_2, &ynorm));
790: /* decide on new step */
791: neP->delta = delta;
792: if (rho > neP->eta1) {
793: rho_satisfied = PETSC_TRUE;
794: } else {
795: rho_satisfied = PETSC_FALSE;
796: PetscCall(PetscInfo(snes, "Trying again in smaller region\n"));
797: /* check to see if progress is hopeless */
798: PetscCall(SNESTR_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP));
799: if (!snes->reason) PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
800: if (snes->reason == SNES_CONVERGED_SNORM_RELATIVE) snes->reason = SNES_DIVERGED_TR_DELTA;
801: snes->numFailures++;
802: /* We're not progressing, so return with the current iterate */
803: if (snes->reason) break;
804: }
805: if (rho_satisfied) {
806: /* Update function values */
807: already_done = PETSC_FALSE;
808: fnorm = gnorm;
809: fk = fkp1;
811: /* New residual and linearization point */
812: PetscCall(VecCopy(G, F));
813: PetscCall(VecCopy(W, X));
815: /* Monitor convergence */
816: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
817: snes->iter++;
818: snes->norm = fnorm;
819: snes->xnorm = xnorm;
820: snes->ynorm = ynorm;
821: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
822: PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
824: /* Test for convergence, xnorm = || X || */
825: PetscCall(VecNorm(X, NORM_2, &xnorm));
826: PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
827: PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
828: if (snes->reason) break;
829: }
830: }
832: if (clear_converged_test) {
833: PetscCall(KSPGetAndClearConvergenceTest(snes->ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy));
834: PetscCall(PetscFree(ctx));
835: PetscCall(KSPSetConvergenceTest(snes->ksp, convtest, convctx, convdestroy));
836: }
837: PetscFunctionReturn(PETSC_SUCCESS);
838: }
840: static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes)
841: {
842: PetscFunctionBegin;
843: PetscCall(SNESSetWorkVecs(snes, 5));
844: PetscCall(SNESSetUpMatrices(snes));
845: PetscFunctionReturn(PETSC_SUCCESS);
846: }
848: static PetscErrorCode SNESReset_NEWTONTR(SNES snes)
849: {
850: SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
852: PetscFunctionBegin;
853: PetscCall(MatDestroy(&tr->qnB));
854: PetscCall(MatDestroy(&tr->qnB_pre));
855: PetscFunctionReturn(PETSC_SUCCESS);
856: }
858: static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes)
859: {
860: PetscFunctionBegin;
861: PetscCall(SNESReset_NEWTONTR(snes));
862: PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonTRSetTolerances_C", NULL));
863: PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonTRGetTolerances_C", NULL));
864: PetscCall(PetscFree(snes->data));
865: PetscFunctionReturn(PETSC_SUCCESS);
866: }
868: static PetscErrorCode SNESSetFromOptions_NEWTONTR(SNES snes, PetscOptionItems PetscOptionsObject)
869: {
870: SNES_NEWTONTR *ctx = (SNES_NEWTONTR *)snes->data;
871: SNESNewtonTRQNType qn;
872: SNESNewtonTRFallbackType fallback;
873: NormType norm;
874: PetscBool flg;
876: PetscFunctionBegin;
877: PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations");
878: PetscCall(PetscOptionsDeprecated("-snes_tr_deltaM", "-snes_tr_deltamax", "3.22", NULL));
879: PetscCall(PetscOptionsReal("-snes_tr_eta1", "eta1", "SNESNewtonTRSetUpdateParameters", ctx->eta1, &ctx->eta1, NULL));
880: PetscCall(PetscOptionsReal("-snes_tr_eta2", "eta2", "SNESNewtonTRSetUpdateParameters", ctx->eta2, &ctx->eta2, NULL));
881: PetscCall(PetscOptionsReal("-snes_tr_eta3", "eta3", "SNESNewtonTRSetUpdateParameters", ctx->eta3, &ctx->eta3, NULL));
882: PetscCall(PetscOptionsReal("-snes_tr_t1", "t1", "SNESNewtonTRSetUpdateParameters", ctx->t1, &ctx->t1, NULL));
883: PetscCall(PetscOptionsReal("-snes_tr_t2", "t2", "SNESNewtonTRSetUpdateParameters", ctx->t2, &ctx->t2, NULL));
884: PetscCall(PetscOptionsReal("-snes_tr_delta0", "Initial trust region size", "SNESNewtonTRSetTolerances", ctx->delta0, &ctx->delta0, NULL));
885: PetscCall(PetscOptionsReal("-snes_tr_deltamin", "Minimum allowed trust region size", "SNESNewtonTRSetTolerances", ctx->deltam, &ctx->deltam, NULL));
886: PetscCall(PetscOptionsReal("-snes_tr_deltamax", "Maximum allowed trust region size", "SNESNewtonTRSetTolerances", ctx->deltaM, &ctx->deltaM, NULL));
887: PetscCall(PetscOptionsReal("-snes_tr_kmdc", "sufficient decrease parameter", "None", ctx->kmdc, &ctx->kmdc, NULL));
889: fallback = ctx->fallback;
890: PetscCall(PetscOptionsEnum("-snes_tr_fallback_type", "Type of fallback if subproblem solution is outside of the trust region", "SNESNewtonTRSetFallbackType", SNESNewtonTRFallbackTypes, (PetscEnum)fallback, (PetscEnum *)&fallback, &flg));
891: if (flg) PetscCall(SNESNewtonTRSetFallbackType(snes, fallback));
893: qn = ctx->qn;
894: PetscCall(PetscOptionsEnum("-snes_tr_qn", "Use Quasi-Newton approximations for the model", "SNESNewtonTRSetQNType", SNESNewtonTRQNTypes, (PetscEnum)qn, (PetscEnum *)&qn, &flg));
895: if (flg) PetscCall(SNESNewtonTRSetQNType(snes, qn));
897: norm = ctx->norm;
898: PetscCall(PetscOptionsEnum("-snes_tr_norm_type", "Type of norm for trust region bounds", "SNESNewtonTRSetNormType", NormTypes, (PetscEnum)norm, (PetscEnum *)&norm, &flg));
899: if (flg) PetscCall(SNESNewtonTRSetNormType(snes, norm));
901: PetscOptionsHeadEnd();
902: PetscFunctionReturn(PETSC_SUCCESS);
903: }
905: static PetscErrorCode SNESView_NEWTONTR(SNES snes, PetscViewer viewer)
906: {
907: SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
908: PetscBool isascii;
910: PetscFunctionBegin;
911: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
912: if (isascii) {
913: PetscCall(PetscViewerASCIIPrintf(viewer, " Trust region parameters:\n"));
914: PetscCall(PetscViewerASCIIPrintf(viewer, " eta1=%g, eta2=%g, eta3=%g\n", (double)tr->eta1, (double)tr->eta2, (double)tr->eta3));
915: PetscCall(PetscViewerASCIIPrintf(viewer, " t1=%g, t2=%g\n", (double)tr->t1, (double)tr->t2));
916: PetscCall(PetscViewerASCIIPrintf(viewer, " delta_min=%g, delta_0=%g, delta_max=%g\n", (double)tr->deltam, (double)tr->delta0, (double)tr->deltaM));
917: PetscCall(PetscViewerASCIIPrintf(viewer, " kmdc=%g\n", (double)tr->kmdc));
918: PetscCall(PetscViewerASCIIPrintf(viewer, " fallback=%s\n", SNESNewtonTRFallbackTypes[tr->fallback]));
919: if (tr->qn) PetscCall(PetscViewerASCIIPrintf(viewer, " qn=%s\n", SNESNewtonTRQNTypes[tr->qn]));
920: if (tr->norm != NORM_2) PetscCall(PetscViewerASCIIPrintf(viewer, " norm=%s\n", NormTypes[tr->norm]));
921: }
922: PetscFunctionReturn(PETSC_SUCCESS);
923: }
925: /*@
926: SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
928: Logically Collective
930: Input Parameters:
931: + snes - the `SNES` context
932: - tol - tolerance
934: Level: deprecated
936: .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`, `SNESSetTolerances()`
937: @*/
938: PetscErrorCode SNESSetTrustRegionTolerance(SNES snes, PetscReal tol)
939: {
940: return SNESNewtonTRSetTolerances(snes, tol, PETSC_CURRENT, PETSC_CURRENT);
941: }
943: /*@
944: SNESNewtonTRSetTolerances - Sets the trust region parameter tolerances.
946: Logically Collective
948: Input Parameters:
949: + snes - the `SNES` context
950: . delta_min - minimum allowed trust region size
951: . delta_max - maximum allowed trust region size
952: - delta_0 - initial trust region size
954: Options Database Key:
955: + -snes_tr_deltamin <tol> - Set minimum size
956: . -snes_tr_deltamax <tol> - Set maximum size
957: - -snes_tr_delta0 <tol> - Set initial size
959: Note:
960: Use `PETSC_DETERMINE` to use the default value for the given `SNES`.
961: Use `PETSC_CURRENT` to retain a value.
963: Fortran Note:
964: Use `PETSC_DETERMINE_REAL`, `PETSC_CURRENT_REAL`
966: Level: intermediate
968: .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`, `SNESNewtonTRGetTolerances()`
969: @*/
970: PetscErrorCode SNESNewtonTRSetTolerances(SNES snes, PetscReal delta_min, PetscReal delta_max, PetscReal delta_0)
971: {
972: PetscFunctionBegin;
977: PetscTryMethod(snes, "SNESNewtonTRSetTolerances_C", (SNES, PetscReal, PetscReal, PetscReal), (snes, delta_min, delta_max, delta_0));
978: PetscFunctionReturn(PETSC_SUCCESS);
979: }
981: /*@
982: SNESNewtonTRGetTolerances - Gets the trust region parameter tolerances.
984: Not Collective
986: Input Parameter:
987: . snes - the `SNES` context
989: Output Parameters:
990: + delta_min - minimum allowed trust region size or `NULL`
991: . delta_max - maximum allowed trust region size or `NULL`
992: - delta_0 - initial trust region size or `NULL`
994: Level: intermediate
996: .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`, `SNESNewtonTRSetTolerances()`
997: @*/
998: PetscErrorCode SNESNewtonTRGetTolerances(SNES snes, PetscReal *delta_min, PetscReal *delta_max, PetscReal *delta_0)
999: {
1000: PetscFunctionBegin;
1002: if (delta_min) PetscAssertPointer(delta_min, 2);
1003: if (delta_max) PetscAssertPointer(delta_max, 3);
1004: if (delta_0) PetscAssertPointer(delta_0, 4);
1005: PetscUseMethod(snes, "SNESNewtonTRGetTolerances_C", (SNES, PetscReal *, PetscReal *, PetscReal *), (snes, delta_min, delta_max, delta_0));
1006: PetscFunctionReturn(PETSC_SUCCESS);
1007: }
1009: /*@
1010: SNESNewtonTRSetUpdateParameters - Sets the trust region update parameters.
1012: Logically Collective
1014: Input Parameters:
1015: + snes - the `SNES` context
1016: . eta1 - acceptance tolerance
1017: . eta2 - shrinking tolerance
1018: . eta3 - enlarging tolerance
1019: . t1 - shrink factor
1020: - t2 - enlarge factor
1022: Options Database Key:
1023: + -snes_tr_eta1 <tol> - Set `eta1`
1024: . -snes_tr_eta2 <tol> - Set `eta2`
1025: . -snes_tr_eta3 <tol> - Set `eta3`
1026: . -snes_tr_t1 <tol> - Set `t1`
1027: - -snes_tr_t2 <tol> - Set `t2`
1029: Notes:
1030: Given the ratio $\rho = \frac{f(x_k) - f(x_k+s_k)}{m(0) - m(s_k)}$, with $x_k$ the current iterate,
1031: $s_k$ the computed step, $f$ the objective function, and $m$ the quadratic model, the trust region
1032: radius is modified as follows
1034: $
1035: \delta =
1036: \begin{cases}
1037: \delta * t_1 ,& \rho < \eta_2 \\
1038: \delta * t_2 ,& \rho > \eta_3 \\
1039: \end{cases}
1040: $
1042: The step is accepted if $\rho > \eta_1$.
1043: Use `PETSC_DETERMINE` to use the default value for the given `SNES`.
1044: Use `PETSC_CURRENT` to retain a value.
1046: Fortran Note:
1047: Use `PETSC_DETERMINE_REAL`, `PETSC_CURRENT_REAL`
1049: Level: intermediate
1051: .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`, `SNESSetObjective()`, `SNESNewtonTRGetUpdateParameters()`
1052: @*/
1053: PetscErrorCode SNESNewtonTRSetUpdateParameters(SNES snes, PetscReal eta1, PetscReal eta2, PetscReal eta3, PetscReal t1, PetscReal t2)
1054: {
1055: PetscBool flg;
1057: PetscFunctionBegin;
1064: PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
1065: if (flg) {
1066: SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
1068: if (eta1 == PETSC_DETERMINE) eta1 = tr->default_eta1;
1069: if (eta2 == PETSC_DETERMINE) eta2 = tr->default_eta2;
1070: if (eta3 == PETSC_DETERMINE) eta3 = tr->default_eta3;
1071: if (t1 == PETSC_DETERMINE) t1 = tr->default_t1;
1072: if (t2 == PETSC_DETERMINE) t2 = tr->default_t2;
1073: if (eta1 != PETSC_CURRENT) tr->eta1 = eta1;
1074: if (eta2 != PETSC_CURRENT) tr->eta2 = eta2;
1075: if (eta3 != PETSC_CURRENT) tr->eta3 = eta3;
1076: if (t1 != PETSC_CURRENT) tr->t1 = t1;
1077: if (t2 != PETSC_CURRENT) tr->t2 = t2;
1078: }
1079: PetscFunctionReturn(PETSC_SUCCESS);
1080: }
1082: /*@
1083: SNESNewtonTRGetUpdateParameters - Gets the trust region update parameters.
1085: Not Collective
1087: Input Parameter:
1088: . snes - the `SNES` context
1090: Output Parameters:
1091: + eta1 - acceptance tolerance
1092: . eta2 - shrinking tolerance
1093: . eta3 - enlarging tolerance
1094: . t1 - shrink factor
1095: - t2 - enlarge factor
1097: Level: intermediate
1099: .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`, `SNESNewtonTRSetUpdateParameters()`
1100: @*/
1101: PetscErrorCode SNESNewtonTRGetUpdateParameters(SNES snes, PetscReal *eta1, PetscReal *eta2, PetscReal *eta3, PetscReal *t1, PetscReal *t2)
1102: {
1103: SNES_NEWTONTR *tr;
1104: PetscBool flg;
1106: PetscFunctionBegin;
1108: if (eta1) PetscAssertPointer(eta1, 2);
1109: if (eta2) PetscAssertPointer(eta2, 3);
1110: if (eta3) PetscAssertPointer(eta3, 4);
1111: if (t1) PetscAssertPointer(t1, 5);
1112: if (t2) PetscAssertPointer(t2, 6);
1113: PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
1114: PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name);
1115: tr = (SNES_NEWTONTR *)snes->data;
1116: if (eta1) *eta1 = tr->eta1;
1117: if (eta2) *eta2 = tr->eta2;
1118: if (eta3) *eta3 = tr->eta3;
1119: if (t1) *t1 = tr->t1;
1120: if (t2) *t2 = tr->t2;
1121: PetscFunctionReturn(PETSC_SUCCESS);
1122: }
1124: /*MC
1125: SNESNEWTONTR - Newton based nonlinear solver that uses a trust-region strategy
1127: Options Database Keys:
1128: + -snes_tr_deltamin <deltamin> - trust region parameter, minimum size of trust region
1129: . -snes_tr_deltamax <deltamax> - trust region parameter, max size of trust region (default: 1e10)
1130: . -snes_tr_delta0 <delta0> - trust region parameter, initial size of trust region (default: 0.2)
1131: . -snes_tr_eta1 <eta1> - trust region parameter $eta1 \le eta2$, $rho > eta1$ breaks out of the inner iteration (default: 0.001)
1132: . -snes_tr_eta2 <eta2> - trust region parameter, $rho \le eta2$ shrinks the trust region (default: 0.25)
1133: . -snes_tr_eta3 <eta3> - trust region parameter $eta3 > eta2$, $rho \ge eta3$ expands the trust region (default: 0.75)
1134: . -snes_tr_t1 <t1> - trust region parameter, shrinking factor of trust region (default: 0.25)
1135: . -snes_tr_t2 <t2> - trust region parameter, expanding factor of trust region (default: 2.0)
1136: . -snes_tr_norm_type <1,2,infinity> - Type of norm for trust region bounds (default: 2)
1137: - -snes_tr_fallback_type <newton,cauchy,dogleg> - Solution strategy to test reduction when step is outside of trust region. Can use scaled Newton direction, Cauchy point (Steepest Descent direction) or dogleg method.
1139: Level: beginner
1141: Notes:
1142: The code is largely based on the book {cite}`nocedal2006numerical` and supports minimizing objective functions using a quadratic model.
1143: Quasi-Newton models are also supported.
1145: Default step computation uses the Newton direction, but a dogleg type update is also supported.
1146: The 1- and infinity-norms are also supported when computing the trust region bounds.
1148: .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESSetObjective()`,
1149: `SNESNewtonTRSetTolerances()`, `SNESNewtonTRSetUpdateParameters()`
1150: `SNESNewtonTRSetNormType()`, `SNESNewtonTRSetFallbackType()`, `SNESNewtonTRSetQNType()`
1151: `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRSetPreCheck()`,
1152: M*/
1153: PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes)
1154: {
1155: SNES_NEWTONTR *neP;
1157: PetscFunctionBegin;
1158: snes->ops->setup = SNESSetUp_NEWTONTR;
1159: snes->ops->solve = SNESSolve_NEWTONTR;
1160: snes->ops->reset = SNESReset_NEWTONTR;
1161: snes->ops->destroy = SNESDestroy_NEWTONTR;
1162: snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR;
1163: snes->ops->view = SNESView_NEWTONTR;
1165: PetscCall(SNESParametersInitialize(snes));
1166: PetscObjectParameterSetDefault(snes, stol, 0.0);
1168: snes->usesksp = PETSC_TRUE;
1169: snes->npcside = PC_RIGHT;
1170: snes->usesnpc = PETSC_TRUE;
1172: snes->alwayscomputesfinalresidual = PETSC_TRUE;
1174: PetscCall(PetscNew(&neP));
1175: snes->data = (void *)neP;
1177: PetscObjectParameterSetDefault(neP, eta1, 0.001);
1178: PetscObjectParameterSetDefault(neP, eta2, 0.25);
1179: PetscObjectParameterSetDefault(neP, eta3, 0.75);
1180: PetscObjectParameterSetDefault(neP, t1, 0.25);
1181: PetscObjectParameterSetDefault(neP, t2, 2.0);
1182: PetscObjectParameterSetDefault(neP, deltam, PetscDefined(USE_REAL_SINGLE) ? 1.e-6 : 1.e-12);
1183: PetscObjectParameterSetDefault(neP, delta0, 0.2);
1184: PetscObjectParameterSetDefault(neP, deltaM, 1.e10);
1186: neP->norm = NORM_2;
1187: neP->fallback = SNES_TR_FALLBACK_NEWTON;
1188: neP->kmdc = 0.0; /* by default do not use sufficient decrease */
1190: PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonTRSetTolerances_C", SNESNewtonTRSetTolerances_TR));
1191: PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonTRGetTolerances_C", SNESNewtonTRGetTolerances_TR));
1192: PetscFunctionReturn(PETSC_SUCCESS);
1193: }