141ba4c6cSHeeho Park #include <../src/snes/impls/ntrdc/ntrdcimpl.h> /*I "petscsnes.h" I*/ 241ba4c6cSHeeho Park 341ba4c6cSHeeho Park typedef struct { 441ba4c6cSHeeho Park SNES snes; 541ba4c6cSHeeho Park /* Information on the regular SNES convergence test; which may have been user provided 641ba4c6cSHeeho Park Copied from tr.c (maybe able to disposed, but this is a private function) - Heeho 741ba4c6cSHeeho Park Same with SNESTR_KSPConverged_Private, SNESTR_KSPConverged_Destroy, and SNESTR_Converged_Private 841ba4c6cSHeeho Park */ 941ba4c6cSHeeho Park 1041ba4c6cSHeeho Park PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 1141ba4c6cSHeeho Park PetscErrorCode (*convdestroy)(void *); 1241ba4c6cSHeeho Park void *convctx; 1341ba4c6cSHeeho Park } SNES_TRDC_KSPConverged_Ctx; 1441ba4c6cSHeeho Park 15*3201ab8dSStefano Zampini static PetscErrorCode SNESNewtonTRSetTolerances_TRDC(SNES snes, PetscReal delta_min, PetscReal delta_max, PetscReal delta_0) 16*3201ab8dSStefano Zampini { 17*3201ab8dSStefano Zampini SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 18*3201ab8dSStefano Zampini 19*3201ab8dSStefano Zampini PetscFunctionBegin; 20*3201ab8dSStefano Zampini if (delta_min == PETSC_DETERMINE) delta_min = 1.e-12; 21*3201ab8dSStefano Zampini if (delta_max == PETSC_DETERMINE) delta_max = 0.5; 22*3201ab8dSStefano Zampini if (delta_0 == PETSC_DETERMINE) delta_0 = 0.1; 23*3201ab8dSStefano Zampini if (delta_min != PETSC_CURRENT) tr->deltatol = delta_min; 24*3201ab8dSStefano Zampini if (delta_max != PETSC_CURRENT) tr->deltaM = delta_max; 25*3201ab8dSStefano Zampini if (delta_0 != PETSC_CURRENT) tr->delta0 = delta_0; 26*3201ab8dSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 27*3201ab8dSStefano Zampini } 28*3201ab8dSStefano Zampini 29d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTRDC_KSPConverged_Private(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *cctx) 30d71ae5a4SJacob Faibussowitsch { 3141ba4c6cSHeeho Park SNES_TRDC_KSPConverged_Ctx *ctx = (SNES_TRDC_KSPConverged_Ctx *)cctx; 3241ba4c6cSHeeho Park SNES snes = ctx->snes; 3341ba4c6cSHeeho Park SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC *)snes->data; 3441ba4c6cSHeeho Park Vec x; 3541ba4c6cSHeeho Park PetscReal nrm; 3641ba4c6cSHeeho Park 3741ba4c6cSHeeho Park PetscFunctionBegin; 389566063dSJacob Faibussowitsch PetscCall((*ctx->convtest)(ksp, n, rnorm, reason, ctx->convctx)); 3948a46eb9SPierre Jolivet if (*reason) PetscCall(PetscInfo(snes, "Default or user provided convergence test KSP iterations=%" PetscInt_FMT ", rnorm=%g\n", n, (double)rnorm)); 4041ba4c6cSHeeho Park /* Determine norm of solution */ 419566063dSJacob Faibussowitsch PetscCall(KSPBuildSolution(ksp, NULL, &x)); 429566063dSJacob Faibussowitsch PetscCall(VecNorm(x, NORM_2, &nrm)); 4341ba4c6cSHeeho Park if (nrm >= neP->delta) { 449566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Ending linear iteration early, delta=%g, length=%g\n", (double)neP->delta, (double)nrm)); 4541ba4c6cSHeeho Park *reason = KSP_CONVERGED_STEP_LENGTH; 4641ba4c6cSHeeho Park } 473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4841ba4c6cSHeeho Park } 4941ba4c6cSHeeho Park 50d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTRDC_KSPConverged_Destroy(void *cctx) 51d71ae5a4SJacob Faibussowitsch { 5241ba4c6cSHeeho Park SNES_TRDC_KSPConverged_Ctx *ctx = (SNES_TRDC_KSPConverged_Ctx *)cctx; 5341ba4c6cSHeeho Park 5441ba4c6cSHeeho Park PetscFunctionBegin; 559566063dSJacob Faibussowitsch PetscCall((*ctx->convdestroy)(ctx->convctx)); 569566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5841ba4c6cSHeeho Park } 5941ba4c6cSHeeho Park 6041ba4c6cSHeeho Park /* 6141ba4c6cSHeeho Park SNESTRDC_Converged_Private -test convergence JUST for 6241ba4c6cSHeeho Park the trust region tolerance. 6341ba4c6cSHeeho Park 6441ba4c6cSHeeho Park */ 65d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTRDC_Converged_Private(SNES snes, PetscInt it, PetscReal xnorm, PetscReal pnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy) 66d71ae5a4SJacob Faibussowitsch { 6741ba4c6cSHeeho Park SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC *)snes->data; 6841ba4c6cSHeeho Park 6941ba4c6cSHeeho Park PetscFunctionBegin; 7041ba4c6cSHeeho Park *reason = SNES_CONVERGED_ITERATING; 71*3201ab8dSStefano Zampini if (neP->delta < xnorm * neP->deltatol) { 72*3201ab8dSStefano Zampini PetscCall(PetscInfo(snes, "Diverged due to too small a trust region %g<%g*%g\n", (double)neP->delta, (double)xnorm, (double)neP->deltatol)); 7341ba4c6cSHeeho Park *reason = SNES_DIVERGED_TR_DELTA; 7441ba4c6cSHeeho Park } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 759566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n", snes->max_funcs)); 7641ba4c6cSHeeho Park *reason = SNES_DIVERGED_FUNCTION_COUNT; 7741ba4c6cSHeeho Park } 783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7941ba4c6cSHeeho Park } 8041ba4c6cSHeeho Park 8141ba4c6cSHeeho Park /*@ 82f6dfbefdSBarry Smith SNESNewtonTRDCGetRhoFlag - Get whether the current solution update is within the trust-region. 8341ba4c6cSHeeho Park 8420f4b53cSBarry Smith Logically Collective 8520f4b53cSBarry Smith 86f6dfbefdSBarry Smith Input Parameter: 8741ba4c6cSHeeho Park . snes - the nonlinear solver object 8841ba4c6cSHeeho Park 89f6dfbefdSBarry Smith Output Parameter: 90420bcc1bSBarry Smith . rho_flag - `PETSC_FALSE` or `PETSC_TRUE` 9141ba4c6cSHeeho Park 9241ba4c6cSHeeho Park Level: developer 9341ba4c6cSHeeho Park 94420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCSetPreCheck()`, 95f6dfbefdSBarry Smith `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()` 9641ba4c6cSHeeho Park @*/ 97d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRDCGetRhoFlag(SNES snes, PetscBool *rho_flag) 98d71ae5a4SJacob Faibussowitsch { 9941ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 10041ba4c6cSHeeho Park 10141ba4c6cSHeeho Park PetscFunctionBegin; 10241ba4c6cSHeeho Park PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1034f572ea9SToby Isaac PetscAssertPointer(rho_flag, 2); 10441ba4c6cSHeeho Park *rho_flag = tr->rho_satisfied; 1053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10641ba4c6cSHeeho Park } 10741ba4c6cSHeeho Park 10841ba4c6cSHeeho Park /*@C 10941ba4c6cSHeeho Park SNESNewtonTRDCSetPreCheck - Sets a user function that is called before the search step has been determined. 11041ba4c6cSHeeho Park Allows the user a chance to change or override the trust region decision. 11141ba4c6cSHeeho Park 112c3339decSBarry Smith Logically Collective 11341ba4c6cSHeeho Park 11441ba4c6cSHeeho Park Input Parameters: 11541ba4c6cSHeeho Park + snes - the nonlinear solver object 11620f4b53cSBarry Smith . func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRDCPreCheck()` 11720f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 11841ba4c6cSHeeho Park 11941ba4c6cSHeeho Park Level: intermediate 12041ba4c6cSHeeho Park 121f6dfbefdSBarry Smith Note: 122f6dfbefdSBarry Smith This function is called BEFORE the function evaluation within the `SNESNEWTONTRDC` solver. 12341ba4c6cSHeeho Park 124420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`, 125f6dfbefdSBarry Smith `SNESNewtonTRDCGetRhoFlag()` 12641ba4c6cSHeeho Park @*/ 127d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRDCSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, PetscBool *, void *), void *ctx) 128d71ae5a4SJacob Faibussowitsch { 12941ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 13041ba4c6cSHeeho Park 13141ba4c6cSHeeho Park PetscFunctionBegin; 13241ba4c6cSHeeho Park PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 13341ba4c6cSHeeho Park if (func) tr->precheck = func; 13441ba4c6cSHeeho Park if (ctx) tr->precheckctx = ctx; 1353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13641ba4c6cSHeeho Park } 13741ba4c6cSHeeho Park 13841ba4c6cSHeeho Park /*@C 139420bcc1bSBarry Smith SNESNewtonTRDCGetPreCheck - Gets the pre-check function optionally set with `SNESNewtonTRDCSetPreCheck()` 14041ba4c6cSHeeho Park 14120f4b53cSBarry Smith Not Collective 14241ba4c6cSHeeho Park 14341ba4c6cSHeeho Park Input Parameter: 14441ba4c6cSHeeho Park . snes - the nonlinear solver context 14541ba4c6cSHeeho Park 14641ba4c6cSHeeho Park Output Parameters: 14720f4b53cSBarry Smith + func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRDCPreCheck()` 14820f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 14941ba4c6cSHeeho Park 15041ba4c6cSHeeho Park Level: intermediate 15141ba4c6cSHeeho Park 152420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCPreCheck()` 15341ba4c6cSHeeho Park @*/ 154d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRDCGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, PetscBool *, void *), void **ctx) 155d71ae5a4SJacob Faibussowitsch { 15641ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 15741ba4c6cSHeeho Park 15841ba4c6cSHeeho Park PetscFunctionBegin; 15941ba4c6cSHeeho Park PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 16041ba4c6cSHeeho Park if (func) *func = tr->precheck; 16141ba4c6cSHeeho Park if (ctx) *ctx = tr->precheckctx; 1623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16341ba4c6cSHeeho Park } 16441ba4c6cSHeeho Park 16541ba4c6cSHeeho Park /*@C 16641ba4c6cSHeeho Park SNESNewtonTRDCSetPostCheck - Sets a user function that is called after the search step has been determined but before the next 16741ba4c6cSHeeho Park function evaluation. Allows the user a chance to change or override the decision of the line search routine 16841ba4c6cSHeeho Park 169c3339decSBarry Smith Logically Collective 17041ba4c6cSHeeho Park 17141ba4c6cSHeeho Park Input Parameters: 17241ba4c6cSHeeho Park + snes - the nonlinear solver object 17320f4b53cSBarry Smith . func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRDCPostCheck()` 17420f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 17541ba4c6cSHeeho Park 17641ba4c6cSHeeho Park Level: intermediate 17741ba4c6cSHeeho Park 178f6dfbefdSBarry Smith Note: 179f6dfbefdSBarry Smith This function is called BEFORE the function evaluation within the `SNESNEWTONTRDC` solver while the function set in 180f6dfbefdSBarry Smith `SNESLineSearchSetPostCheck()` is called AFTER the function evaluation. 18141ba4c6cSHeeho Park 182420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCPostCheck()`, `SNESNewtonTRDCGetPostCheck()`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCGetPreCheck()` 18341ba4c6cSHeeho Park @*/ 184d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRDCSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx) 185d71ae5a4SJacob Faibussowitsch { 18641ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 18741ba4c6cSHeeho Park 18841ba4c6cSHeeho Park PetscFunctionBegin; 18941ba4c6cSHeeho Park PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 19041ba4c6cSHeeho Park if (func) tr->postcheck = func; 19141ba4c6cSHeeho Park if (ctx) tr->postcheckctx = ctx; 1923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19341ba4c6cSHeeho Park } 19441ba4c6cSHeeho Park 19541ba4c6cSHeeho Park /*@C 196420bcc1bSBarry Smith SNESNewtonTRDCGetPostCheck - Gets the post-check function optionally set with `SNESNewtonTRDCSetPostCheck()` 19741ba4c6cSHeeho Park 19820f4b53cSBarry Smith Not Collective 19941ba4c6cSHeeho Park 20041ba4c6cSHeeho Park Input Parameter: 20141ba4c6cSHeeho Park . snes - the nonlinear solver context 20241ba4c6cSHeeho Park 20341ba4c6cSHeeho Park Output Parameters: 20420f4b53cSBarry Smith + func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRDCPostCheck()` 20520f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 20641ba4c6cSHeeho Park 20741ba4c6cSHeeho Park Level: intermediate 20841ba4c6cSHeeho Park 209420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCPostCheck()`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCGetPreCheck()` 21041ba4c6cSHeeho Park @*/ 211d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRDCGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) 212d71ae5a4SJacob Faibussowitsch { 21341ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 21441ba4c6cSHeeho Park 21541ba4c6cSHeeho Park PetscFunctionBegin; 21641ba4c6cSHeeho Park PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 21741ba4c6cSHeeho Park if (func) *func = tr->postcheck; 21841ba4c6cSHeeho Park if (ctx) *ctx = tr->postcheckctx; 2193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 22041ba4c6cSHeeho Park } 22141ba4c6cSHeeho Park 222ceaaa498SBarry Smith // PetscClangLinter pragma disable: -fdoc-internal-linkage 22341ba4c6cSHeeho Park /*@C 224f6dfbefdSBarry Smith SNESNewtonTRDCPreCheck - Called before the step has been determined in `SNESNEWTONTRDC` 22541ba4c6cSHeeho Park 226c3339decSBarry Smith Logically Collective 22741ba4c6cSHeeho Park 22841ba4c6cSHeeho Park Input Parameters: 22941ba4c6cSHeeho Park + snes - the solver 23041ba4c6cSHeeho Park . X - The last solution 23141ba4c6cSHeeho Park - Y - The step direction 23241ba4c6cSHeeho Park 2332fe279fdSBarry Smith Output Parameter: 23420f4b53cSBarry Smith . changed_Y - Indicator that the step direction `Y` has been changed. 23541ba4c6cSHeeho Park 23641ba4c6cSHeeho Park Level: developer 23741ba4c6cSHeeho Park 238420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCPostCheck()` 23941ba4c6cSHeeho Park @*/ 240d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNewtonTRDCPreCheck(SNES snes, Vec X, Vec Y, PetscBool *changed_Y) 241d71ae5a4SJacob Faibussowitsch { 24241ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 24341ba4c6cSHeeho Park 24441ba4c6cSHeeho Park PetscFunctionBegin; 24541ba4c6cSHeeho Park *changed_Y = PETSC_FALSE; 24641ba4c6cSHeeho Park if (tr->precheck) { 2479566063dSJacob Faibussowitsch PetscCall((*tr->precheck)(snes, X, Y, changed_Y, tr->precheckctx)); 24841ba4c6cSHeeho Park PetscValidLogicalCollectiveBool(snes, *changed_Y, 4); 24941ba4c6cSHeeho Park } 2503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 25141ba4c6cSHeeho Park } 25241ba4c6cSHeeho Park 253ceaaa498SBarry Smith // PetscClangLinter pragma disable: -fdoc-internal-linkage 25441ba4c6cSHeeho Park /*@C 255f6dfbefdSBarry Smith SNESNewtonTRDCPostCheck - Called after the step has been determined in `SNESNEWTONTRDC` but before the function evaluation at that step 25641ba4c6cSHeeho Park 257c3339decSBarry Smith Logically Collective 25841ba4c6cSHeeho Park 25941ba4c6cSHeeho Park Input Parameters: 260f1a722f8SMatthew G. Knepley + snes - the solver 261f1a722f8SMatthew G. Knepley . X - The last solution 26241ba4c6cSHeeho Park . Y - The full step direction 26341ba4c6cSHeeho Park - W - The updated solution, W = X - Y 26441ba4c6cSHeeho Park 26541ba4c6cSHeeho Park Output Parameters: 26641ba4c6cSHeeho Park + changed_Y - indicator if step has been changed 26720f4b53cSBarry Smith - changed_W - Indicator if the new candidate solution `W` has been changed. 26841ba4c6cSHeeho Park 2692fe279fdSBarry Smith Level: developer 2702fe279fdSBarry Smith 271f6dfbefdSBarry Smith Note: 27220f4b53cSBarry Smith If `Y` is changed then `W` is recomputed as `X` - `Y` 27341ba4c6cSHeeho Park 274420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`, `SNESNewtonTRDCPreCheck() 27541ba4c6cSHeeho Park @*/ 276d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNewtonTRDCPostCheck(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) 277d71ae5a4SJacob Faibussowitsch { 27841ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 27941ba4c6cSHeeho Park 28041ba4c6cSHeeho Park PetscFunctionBegin; 28141ba4c6cSHeeho Park *changed_Y = PETSC_FALSE; 28241ba4c6cSHeeho Park *changed_W = PETSC_FALSE; 28341ba4c6cSHeeho Park if (tr->postcheck) { 2849566063dSJacob Faibussowitsch PetscCall((*tr->postcheck)(snes, X, Y, W, changed_Y, changed_W, tr->postcheckctx)); 28541ba4c6cSHeeho Park PetscValidLogicalCollectiveBool(snes, *changed_Y, 5); 28641ba4c6cSHeeho Park PetscValidLogicalCollectiveBool(snes, *changed_W, 6); 28741ba4c6cSHeeho Park } 2883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 28941ba4c6cSHeeho Park } 29041ba4c6cSHeeho Park 29141ba4c6cSHeeho Park /* 29241ba4c6cSHeeho Park SNESSolve_NEWTONTRDC - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy 29341ba4c6cSHeeho Park (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of 29441ba4c6cSHeeho Park nonlinear equations 29541ba4c6cSHeeho Park 29641ba4c6cSHeeho Park */ 297d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NEWTONTRDC(SNES snes) 298d71ae5a4SJacob Faibussowitsch { 29941ba4c6cSHeeho Park SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC *)snes->data; 30041ba4c6cSHeeho Park Vec X, F, Y, G, W, GradF, YNtmp; 30141ba4c6cSHeeho Park Vec YCtmp; 30241ba4c6cSHeeho Park Mat jac; 30341ba4c6cSHeeho Park PetscInt maxits, i, j, lits, inner_count, bs; 30441ba4c6cSHeeho Park PetscReal rho, fnorm, gnorm, xnorm = 0, delta, ynorm, temp_xnorm, temp_ynorm; /* TRDC inner iteration */ 30541ba4c6cSHeeho Park PetscReal inorms[99]; /* need to make it dynamic eventually, fixed max block size of 99 for now */ 30641ba4c6cSHeeho Park PetscReal deltaM, ynnorm, f0, mp, gTy, g, yTHy; /* rho calculation */ 30741ba4c6cSHeeho Park PetscReal auk, gfnorm, ycnorm, c0, c1, c2, tau, tau_pos, tau_neg, gTBg; /* Cauchy Point */ 30841ba4c6cSHeeho Park KSP ksp; 30941ba4c6cSHeeho Park SNESConvergedReason reason = SNES_CONVERGED_ITERATING; 31041ba4c6cSHeeho Park PetscBool breakout = PETSC_FALSE; 31141ba4c6cSHeeho Park SNES_TRDC_KSPConverged_Ctx *ctx; 31241ba4c6cSHeeho Park PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), (*convdestroy)(void *); 31341ba4c6cSHeeho Park void *convctx; 31441ba4c6cSHeeho Park 31541ba4c6cSHeeho Park PetscFunctionBegin; 31641ba4c6cSHeeho Park maxits = snes->max_its; /* maximum number of iterations */ 31741ba4c6cSHeeho Park X = snes->vec_sol; /* solution vector */ 31841ba4c6cSHeeho Park F = snes->vec_func; /* residual vector */ 31941ba4c6cSHeeho Park Y = snes->work[0]; /* update vector */ 32041ba4c6cSHeeho Park G = snes->work[1]; /* updated residual */ 32141ba4c6cSHeeho Park W = snes->work[2]; /* temporary vector */ 32241ba4c6cSHeeho Park GradF = snes->work[3]; /* grad f = J^T F */ 32341ba4c6cSHeeho Park YNtmp = snes->work[4]; /* Newton solution */ 32441ba4c6cSHeeho Park YCtmp = snes->work[5]; /* Cauchy solution */ 32541ba4c6cSHeeho Park 3260b121fc5SBarry Smith 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); 32741ba4c6cSHeeho Park 3289566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(YNtmp, &bs)); 32941ba4c6cSHeeho Park 3309566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 33141ba4c6cSHeeho Park snes->iter = 0; 3329566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 33341ba4c6cSHeeho Park 33441ba4c6cSHeeho Park /* Set the linear stopping criteria to use the More' trick. From tr.c */ 3359566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 3369566063dSJacob Faibussowitsch PetscCall(KSPGetConvergenceTest(ksp, &convtest, &convctx, &convdestroy)); 33741ba4c6cSHeeho Park if (convtest != SNESTRDC_KSPConverged_Private) { 3389566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 33941ba4c6cSHeeho Park ctx->snes = snes; 3409566063dSJacob Faibussowitsch PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 3419566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(ksp, SNESTRDC_KSPConverged_Private, ctx, SNESTRDC_KSPConverged_Destroy)); 3429566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTRDC_KSPConverged_Private\n")); 34341ba4c6cSHeeho Park } 34441ba4c6cSHeeho Park 34541ba4c6cSHeeho Park if (!snes->vec_func_init_set) { 3469566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */ 34741ba4c6cSHeeho Park } else snes->vec_func_init_set = PETSC_FALSE; 34841ba4c6cSHeeho Park 3499566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */ 35041ba4c6cSHeeho Park SNESCheckFunctionNorm(snes, fnorm); 3519566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */ 35241ba4c6cSHeeho Park 3539566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 35441ba4c6cSHeeho Park snes->norm = fnorm; 3559566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 35641ba4c6cSHeeho Park delta = xnorm ? neP->delta0 * xnorm : neP->delta0; /* initial trust region size scaled by xnorm */ 35741ba4c6cSHeeho Park deltaM = xnorm ? neP->deltaM * xnorm : neP->deltaM; /* maximum trust region size scaled by xnorm */ 35841ba4c6cSHeeho Park neP->delta = delta; 3599566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 3609566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, fnorm)); 36141ba4c6cSHeeho Park 36241ba4c6cSHeeho Park neP->rho_satisfied = PETSC_FALSE; 36341ba4c6cSHeeho Park 36441ba4c6cSHeeho Park /* test convergence */ 365dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 3663ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 36741ba4c6cSHeeho Park 36841ba4c6cSHeeho Park for (i = 0; i < maxits; i++) { 36941ba4c6cSHeeho Park PetscBool changed_y; 37041ba4c6cSHeeho Park PetscBool changed_w; 37141ba4c6cSHeeho Park 37241ba4c6cSHeeho Park /* dogleg method */ 3739566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 37441ba4c6cSHeeho Park SNESCheckJacobianDomainerror(snes); 3759566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian)); 3769566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, F, YNtmp)); /* Quasi Newton Solution */ 37741ba4c6cSHeeho Park SNESCheckKSPSolve(snes); /* this is necessary but old tr.c did not have it*/ 3789566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 3799566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(snes, &jac, NULL, NULL, NULL)); 38041ba4c6cSHeeho Park 38141ba4c6cSHeeho Park /* rescale Jacobian, Newton solution update, and re-calculate delta for multiphase (multivariable) 38241ba4c6cSHeeho Park for inner iteration and Cauchy direction calculation 38341ba4c6cSHeeho Park */ 38441ba4c6cSHeeho Park if (bs > 1 && neP->auto_scale_multiphase) { 3859566063dSJacob Faibussowitsch PetscCall(VecStrideNormAll(YNtmp, NORM_INFINITY, inorms)); 38641ba4c6cSHeeho Park for (j = 0; j < bs; j++) { 38741ba4c6cSHeeho Park if (neP->auto_scale_max > 1.0) { 388ad540459SPierre Jolivet if (inorms[j] < 1.0 / neP->auto_scale_max) inorms[j] = 1.0 / neP->auto_scale_max; 38941ba4c6cSHeeho Park } 3909566063dSJacob Faibussowitsch PetscCall(VecStrideSet(W, j, inorms[j])); 3919566063dSJacob Faibussowitsch PetscCall(VecStrideScale(YNtmp, j, 1.0 / inorms[j])); 3929566063dSJacob Faibussowitsch PetscCall(VecStrideScale(X, j, 1.0 / inorms[j])); 39341ba4c6cSHeeho Park } 3949566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, &xnorm)); 39541ba4c6cSHeeho Park if (i == 0) { 39641ba4c6cSHeeho Park delta = neP->delta0 * xnorm; 39741ba4c6cSHeeho Park } else { 39841ba4c6cSHeeho Park delta = neP->delta * xnorm; 39941ba4c6cSHeeho Park } 40041ba4c6cSHeeho Park deltaM = neP->deltaM * xnorm; 401f3fa974cSJacob Faibussowitsch PetscCall(MatDiagonalScale(jac, NULL, W)); 40241ba4c6cSHeeho Park } 40341ba4c6cSHeeho Park 40441ba4c6cSHeeho Park /* calculating GradF of minimization function */ 4059566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(jac, F, GradF)); /* grad f = J^T F */ 4069566063dSJacob Faibussowitsch PetscCall(VecNorm(YNtmp, NORM_2, &ynnorm)); /* ynnorm <- || Y_newton || */ 40741ba4c6cSHeeho Park 40841ba4c6cSHeeho Park inner_count = 0; 40941ba4c6cSHeeho Park neP->rho_satisfied = PETSC_FALSE; 41041ba4c6cSHeeho Park while (1) { 41141ba4c6cSHeeho Park if (ynnorm <= delta) { /* see if the Newton solution is within the trust region */ 4129566063dSJacob Faibussowitsch PetscCall(VecCopy(YNtmp, Y)); 41341ba4c6cSHeeho Park } else if (neP->use_cauchy) { /* use Cauchy direction if enabled */ 4149566063dSJacob Faibussowitsch PetscCall(MatMult(jac, GradF, W)); 4159566063dSJacob Faibussowitsch PetscCall(VecDotRealPart(W, W, &gTBg)); /* completes GradF^T J^T J GradF */ 4169566063dSJacob Faibussowitsch PetscCall(VecNorm(GradF, NORM_2, &gfnorm)); /* grad f norm <- || grad f || */ 41741ba4c6cSHeeho Park if (gTBg <= 0.0) { 41841ba4c6cSHeeho Park auk = PETSC_MAX_REAL; 41941ba4c6cSHeeho Park } else { 42041ba4c6cSHeeho Park auk = PetscSqr(gfnorm) / gTBg; 42141ba4c6cSHeeho Park } 42241ba4c6cSHeeho Park auk = PetscMin(delta / gfnorm, auk); 4239566063dSJacob Faibussowitsch PetscCall(VecCopy(GradF, YCtmp)); /* this could be improved */ 4249566063dSJacob Faibussowitsch PetscCall(VecScale(YCtmp, auk)); /* YCtmp, Cauchy solution*/ 4259566063dSJacob Faibussowitsch PetscCall(VecNorm(YCtmp, NORM_2, &ycnorm)); /* ycnorm <- || Y_cauchy || */ 42641ba4c6cSHeeho Park if (ycnorm >= delta) { /* see if the Cauchy solution meets the criteria */ 4279566063dSJacob Faibussowitsch PetscCall(VecCopy(YCtmp, Y)); 4289566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "DL evaluated. delta: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n", (double)delta, (double)ynnorm, (double)ycnorm)); 42941ba4c6cSHeeho Park } else { /* take ratio, tau, of Cauchy and Newton direction and step */ 4309566063dSJacob Faibussowitsch PetscCall(VecAXPY(YNtmp, -1.0, YCtmp)); /* YCtmp = A, YNtmp = B */ 4319566063dSJacob Faibussowitsch PetscCall(VecNorm(YNtmp, NORM_2, &c0)); /* this could be improved */ 43241ba4c6cSHeeho Park c0 = PetscSqr(c0); 4339566063dSJacob Faibussowitsch PetscCall(VecDotRealPart(YCtmp, YNtmp, &c1)); 43441ba4c6cSHeeho Park c1 = 2.0 * c1; 4359566063dSJacob Faibussowitsch PetscCall(VecNorm(YCtmp, NORM_2, &c2)); /* this could be improved */ 43641ba4c6cSHeeho Park c2 = PetscSqr(c2) - PetscSqr(delta); 43741ba4c6cSHeeho Park tau_pos = (c1 + PetscSqrtReal(PetscSqr(c1) - 4. * c0 * c2)) / (2. * c0); /* quadratic formula */ 43841ba4c6cSHeeho Park tau_neg = (c1 - PetscSqrtReal(PetscSqr(c1) - 4. * c0 * c2)) / (2. * c0); 43941ba4c6cSHeeho Park tau = PetscMax(tau_pos, tau_neg); /* can tau_neg > tau_pos? I don't think so, but just in case. */ 4409566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "DL evaluated. tau: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n", (double)tau, (double)ynnorm, (double)ycnorm)); 4419566063dSJacob Faibussowitsch PetscCall(VecWAXPY(W, tau, YNtmp, YCtmp)); 4429566063dSJacob Faibussowitsch PetscCall(VecAXPY(W, -tau, YCtmp)); 4439566063dSJacob Faibussowitsch PetscCall(VecCopy(W, Y)); /* this could be improved */ 44441ba4c6cSHeeho Park } 44541ba4c6cSHeeho Park } else { 44641ba4c6cSHeeho Park /* if Cauchy is disabled, only use Newton direction */ 44741ba4c6cSHeeho Park auk = delta / ynnorm; 4489566063dSJacob Faibussowitsch PetscCall(VecScale(YNtmp, auk)); 4499566063dSJacob Faibussowitsch PetscCall(VecCopy(YNtmp, Y)); /* this could be improved (many VecCopy, VecNorm)*/ 45041ba4c6cSHeeho Park } 45141ba4c6cSHeeho Park 4529566063dSJacob Faibussowitsch PetscCall(VecNorm(Y, NORM_2, &ynorm)); /* compute the final ynorm */ 45341ba4c6cSHeeho Park f0 = 0.5 * PetscSqr(fnorm); /* minimizing function f(X) */ 4549566063dSJacob Faibussowitsch PetscCall(MatMult(jac, Y, W)); 4559566063dSJacob Faibussowitsch PetscCall(VecDotRealPart(W, W, &yTHy)); /* completes GradY^T J^T J GradY */ 4569566063dSJacob Faibussowitsch PetscCall(VecDotRealPart(GradF, Y, &gTy)); 45741ba4c6cSHeeho Park mp = f0 - gTy + 0.5 * yTHy; /* quadratic model to satisfy, -gTy because our update is X-Y*/ 45841ba4c6cSHeeho Park 45941ba4c6cSHeeho Park /* scale back solution update */ 46041ba4c6cSHeeho Park if (bs > 1 && neP->auto_scale_multiphase) { 46141ba4c6cSHeeho Park for (j = 0; j < bs; j++) { 4629566063dSJacob Faibussowitsch PetscCall(VecStrideScale(Y, j, inorms[j])); 46341ba4c6cSHeeho Park if (inner_count == 0) { 46441ba4c6cSHeeho Park /* TRDC inner algorithm does not need scaled X after calculating delta in the outer iteration */ 46541ba4c6cSHeeho Park /* need to scale back X to match Y and provide proper update to the external code */ 4669566063dSJacob Faibussowitsch PetscCall(VecStrideScale(X, j, inorms[j])); 46741ba4c6cSHeeho Park } 46841ba4c6cSHeeho Park } 4699566063dSJacob Faibussowitsch if (inner_count == 0) PetscCall(VecNorm(X, NORM_2, &temp_xnorm)); /* only in the first iteration */ 4709566063dSJacob Faibussowitsch PetscCall(VecNorm(Y, NORM_2, &temp_ynorm)); 47141ba4c6cSHeeho Park } else { 47241ba4c6cSHeeho Park temp_xnorm = xnorm; 47341ba4c6cSHeeho Park temp_ynorm = ynorm; 47441ba4c6cSHeeho Park } 47541ba4c6cSHeeho Park inner_count++; 47641ba4c6cSHeeho Park 47741ba4c6cSHeeho Park /* Evaluate the solution to meet the improvement ratio criteria */ 4789566063dSJacob Faibussowitsch PetscCall(SNESNewtonTRDCPreCheck(snes, X, Y, &changed_y)); 4799566063dSJacob Faibussowitsch PetscCall(VecWAXPY(W, -1.0, Y, X)); 4809566063dSJacob Faibussowitsch PetscCall(SNESNewtonTRDCPostCheck(snes, X, Y, W, &changed_y, &changed_w)); 4819566063dSJacob Faibussowitsch if (changed_y) PetscCall(VecWAXPY(W, -1.0, Y, X)); 4829566063dSJacob Faibussowitsch PetscCall(VecCopy(Y, snes->vec_sol_update)); 4839566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, W, G)); /* F(X-Y) = G */ 4849566063dSJacob Faibussowitsch PetscCall(VecNorm(G, NORM_2, &gnorm)); /* gnorm <- || g || */ 48541ba4c6cSHeeho Park SNESCheckFunctionNorm(snes, gnorm); 48641ba4c6cSHeeho Park g = 0.5 * PetscSqr(gnorm); /* minimizing function g(W) */ 48741ba4c6cSHeeho Park if (f0 == mp) rho = 0.0; 48841ba4c6cSHeeho Park else rho = (f0 - g) / (f0 - mp); /* actual improvement over predicted improvement */ 48941ba4c6cSHeeho Park 49041ba4c6cSHeeho Park if (rho < neP->eta2) { 49141ba4c6cSHeeho Park delta *= neP->t1; /* shrink the region */ 49241ba4c6cSHeeho Park } else if (rho > neP->eta3) { 49341ba4c6cSHeeho Park delta = PetscMin(neP->t2 * delta, deltaM); /* expand the region, but not greater than deltaM */ 49441ba4c6cSHeeho Park } 49541ba4c6cSHeeho Park 49641ba4c6cSHeeho Park neP->delta = delta; 49741ba4c6cSHeeho Park if (rho >= neP->eta1) { 49841ba4c6cSHeeho Park /* unscale delta and xnorm before going to the next outer iteration */ 49941ba4c6cSHeeho Park if (bs > 1 && neP->auto_scale_multiphase) { 50041ba4c6cSHeeho Park neP->delta = delta / xnorm; 50141ba4c6cSHeeho Park xnorm = temp_xnorm; 50241ba4c6cSHeeho Park ynorm = temp_ynorm; 50341ba4c6cSHeeho Park } 50441ba4c6cSHeeho Park neP->rho_satisfied = PETSC_TRUE; 50541ba4c6cSHeeho Park break; /* the improvement ratio is satisfactory */ 50641ba4c6cSHeeho Park } 5079566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Trying again in smaller region\n")); 50841ba4c6cSHeeho Park 50941ba4c6cSHeeho Park /* check to see if progress is hopeless */ 51041ba4c6cSHeeho Park neP->itflag = PETSC_FALSE; 51141ba4c6cSHeeho Park /* both delta, ynorm, and xnorm are either scaled or unscaled */ 5129566063dSJacob Faibussowitsch PetscCall(SNESTRDC_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP)); 51341ba4c6cSHeeho Park /* if multiphase state changes, break out inner iteration */ 51441ba4c6cSHeeho Park if (reason == SNES_BREAKOUT_INNER_ITER) { 51541ba4c6cSHeeho Park if (bs > 1 && neP->auto_scale_multiphase) { 51641ba4c6cSHeeho Park /* unscale delta and xnorm before going to the next outer iteration */ 51741ba4c6cSHeeho Park neP->delta = delta / xnorm; 51841ba4c6cSHeeho Park xnorm = temp_xnorm; 51941ba4c6cSHeeho Park ynorm = temp_ynorm; 52041ba4c6cSHeeho Park } 52141ba4c6cSHeeho Park reason = SNES_CONVERGED_ITERATING; 52241ba4c6cSHeeho Park break; 52341ba4c6cSHeeho Park } 52441ba4c6cSHeeho Park if (reason == SNES_CONVERGED_SNORM_RELATIVE) reason = SNES_DIVERGED_INNER; 52541ba4c6cSHeeho Park if (reason) { 52641ba4c6cSHeeho Park if (reason < 0) { 52741ba4c6cSHeeho Park /* We're not progressing, so return with the current iterate */ 5289566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, i + 1, fnorm)); 52941ba4c6cSHeeho Park breakout = PETSC_TRUE; 53041ba4c6cSHeeho Park break; 53141ba4c6cSHeeho Park } else if (reason > 0) { 53241ba4c6cSHeeho Park /* We're converged, so return with the current iterate and update solution */ 5339566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, i + 1, fnorm)); 53441ba4c6cSHeeho Park breakout = PETSC_FALSE; 53541ba4c6cSHeeho Park break; 53641ba4c6cSHeeho Park } 53741ba4c6cSHeeho Park } 53841ba4c6cSHeeho Park snes->numFailures++; 53941ba4c6cSHeeho Park } 54041ba4c6cSHeeho Park if (!breakout) { 54141ba4c6cSHeeho Park /* Update function and solution vectors */ 54241ba4c6cSHeeho Park fnorm = gnorm; 5439566063dSJacob Faibussowitsch PetscCall(VecCopy(G, F)); 5449566063dSJacob Faibussowitsch PetscCall(VecCopy(W, X)); 54541ba4c6cSHeeho Park /* Monitor convergence */ 5469566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 54741ba4c6cSHeeho Park snes->iter = i + 1; 54841ba4c6cSHeeho Park snes->norm = fnorm; 54941ba4c6cSHeeho Park snes->xnorm = xnorm; 55041ba4c6cSHeeho Park snes->ynorm = ynorm; 5519566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 5529566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 5539566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 55441ba4c6cSHeeho Park /* Test for convergence, xnorm = || X || */ 55541ba4c6cSHeeho Park neP->itflag = PETSC_TRUE; 5569566063dSJacob Faibussowitsch if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm)); 557dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP); 55841ba4c6cSHeeho Park if (reason) break; 55941ba4c6cSHeeho Park } else break; 56041ba4c6cSHeeho Park } 56141ba4c6cSHeeho Park 5629566063dSJacob Faibussowitsch /* PetscCall(PetscFree(inorms)); */ 56341ba4c6cSHeeho Park if (i == maxits) { 5649566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits)); 56541ba4c6cSHeeho Park if (!reason) reason = SNES_DIVERGED_MAX_IT; 56641ba4c6cSHeeho Park } 5679566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 56841ba4c6cSHeeho Park snes->reason = reason; 5699566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 57041ba4c6cSHeeho Park if (convtest != SNESTRDC_KSPConverged_Private) { 5719566063dSJacob Faibussowitsch PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 5729566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 5739566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(ksp, convtest, convctx, convdestroy)); 57441ba4c6cSHeeho Park } 5753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 57641ba4c6cSHeeho Park } 57741ba4c6cSHeeho Park 578d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_NEWTONTRDC(SNES snes) 579d71ae5a4SJacob Faibussowitsch { 58041ba4c6cSHeeho Park PetscFunctionBegin; 5819566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes, 6)); 5829566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 5833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 58441ba4c6cSHeeho Park } 58541ba4c6cSHeeho Park 58666976f2fSJacob Faibussowitsch static PetscErrorCode SNESReset_NEWTONTRDC(SNES snes) 587d71ae5a4SJacob Faibussowitsch { 58841ba4c6cSHeeho Park PetscFunctionBegin; 5893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 59041ba4c6cSHeeho Park } 59141ba4c6cSHeeho Park 592d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_NEWTONTRDC(SNES snes) 593d71ae5a4SJacob Faibussowitsch { 59441ba4c6cSHeeho Park PetscFunctionBegin; 5959566063dSJacob Faibussowitsch PetscCall(SNESReset_NEWTONTRDC(snes)); 596*3201ab8dSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonTRSetTolerances_C", NULL)); 5979566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 5983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 59941ba4c6cSHeeho Park } 60041ba4c6cSHeeho Park 601d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NEWTONTRDC(SNES snes, PetscOptionItems *PetscOptionsObject) 602d71ae5a4SJacob Faibussowitsch { 60341ba4c6cSHeeho Park SNES_NEWTONTRDC *ctx = (SNES_NEWTONTRDC *)snes->data; 60441ba4c6cSHeeho Park 60541ba4c6cSHeeho Park PetscFunctionBegin; 606d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations"); 607*3201ab8dSStefano Zampini PetscCall(PetscOptionsReal("-snes_trdc_tol", "Trust region tolerance", "SNESNewtonTRSetTolerances", ctx->deltatol, &ctx->deltatol, NULL)); 6089566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_eta1", "eta1", "None", ctx->eta1, &ctx->eta1, NULL)); 6099566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_eta2", "eta2", "None", ctx->eta2, &ctx->eta2, NULL)); 6109566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_eta3", "eta3", "None", ctx->eta3, &ctx->eta3, NULL)); 6119566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_t1", "t1", "None", ctx->t1, &ctx->t1, NULL)); 6129566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_t2", "t2", "None", ctx->t2, &ctx->t2, NULL)); 6139566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_deltaM", "deltaM", "None", ctx->deltaM, &ctx->deltaM, NULL)); 6149566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_delta0", "delta0", "None", ctx->delta0, &ctx->delta0, NULL)); 6159566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_auto_scale_max", "auto_scale_max", "None", ctx->auto_scale_max, &ctx->auto_scale_max, NULL)); 6169566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_trdc_use_cauchy", "use_cauchy", "use Cauchy step and direction", ctx->use_cauchy, &ctx->use_cauchy, NULL)); 6179566063dSJacob Faibussowitsch 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)); 618d0609cedSBarry Smith PetscOptionsHeadEnd(); 6193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 62041ba4c6cSHeeho Park } 62141ba4c6cSHeeho Park 622d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NEWTONTRDC(SNES snes, PetscViewer viewer) 623d71ae5a4SJacob Faibussowitsch { 62441ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 62541ba4c6cSHeeho Park PetscBool iascii; 62641ba4c6cSHeeho Park 62741ba4c6cSHeeho Park PetscFunctionBegin; 6289566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 62941ba4c6cSHeeho Park if (iascii) { 630*3201ab8dSStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " Trust region tolerance %g\n", (double)tr->deltatol)); 6319566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " eta1=%g, eta2=%g, eta3=%g\n", (double)tr->eta1, (double)tr->eta2, (double)tr->eta3)); 6329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " delta0=%g, t1=%g, t2=%g, deltaM=%g\n", (double)tr->delta0, (double)tr->t1, (double)tr->t2, (double)tr->deltaM)); 63341ba4c6cSHeeho Park } 6343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 63541ba4c6cSHeeho Park } 636f6dfbefdSBarry Smith 63741ba4c6cSHeeho Park /*MC 63841ba4c6cSHeeho Park SNESNEWTONTRDC - Newton based nonlinear solver that uses trust-region dogleg method with Cauchy direction 63941ba4c6cSHeeho Park 640f6dfbefdSBarry Smith Options Database Keys: 64141ba4c6cSHeeho Park + -snes_trdc_tol <tol> - trust region tolerance 64241ba4c6cSHeeho Park . -snes_trdc_eta1 <eta1> - trust region parameter 0.0 <= eta1 <= eta2, rho >= eta1 breaks out of the inner iteration (default: eta1=0.001) 64341ba4c6cSHeeho Park . -snes_trdc_eta2 <eta2> - trust region parameter 0.0 <= eta1 <= eta2, rho <= eta2 shrinks the trust region (default: eta2=0.25) 64441ba4c6cSHeeho Park . -snes_trdc_eta3 <eta3> - trust region parameter eta3 > eta2, rho >= eta3 expands the trust region (default: eta3=0.75) 64541ba4c6cSHeeho Park . -snes_trdc_t1 <t1> - trust region parameter, shrinking factor of trust region (default: 0.25) 64641ba4c6cSHeeho Park . -snes_trdc_t2 <t2> - trust region parameter, expanding factor of trust region (default: 2.0) 6471d27aa22SBarry Smith . -snes_trdc_deltaM <deltaM> - trust region parameter, max size of trust region, $deltaM*norm2(x)$ (default: 0.5) 6481d27aa22SBarry Smith . -snes_trdc_delta0 <delta0> - trust region parameter, initial size of trust region, $delta0*norm2(x)$ (default: 0.1) 64941ba4c6cSHeeho Park . -snes_trdc_auto_scale_max <auto_scale_max> - used with auto_scale_multiphase, caps the maximum auto-scaling factor 65041ba4c6cSHeeho Park . -snes_trdc_use_cauchy <use_cauchy> - True uses dogleg Cauchy (Steepest Descent direction) step & direction in the trust region algorithm 65141ba4c6cSHeeho Park - -snes_trdc_auto_scale_multiphase <auto_scale_multiphase> - True turns on auto-scaling for multivariable block matrix for Cauchy and trust region 65241ba4c6cSHeeho Park 653*3201ab8dSStefano Zampini Level: advanced 65420f4b53cSBarry Smith 655*3201ab8dSStefano Zampini Notes: 656*3201ab8dSStefano Zampini `SNESNEWTONTRDC` only works for root-finding problems and does not support objective functions. 657*3201ab8dSStefano Zampini The main difference with respect to `SNESNEWTONTR` is that `SNESNEWTONTRDC` scales the trust region by the norm of the current linearization point. 658*3201ab8dSStefano Zampini Future version may extend the `SNESNEWTONTR` code and deprecate `SNESNEWTONTRDC`. 65941ba4c6cSHeeho Park 660*3201ab8dSStefano Zampini For details, see {cite}`park2021linear` 661*3201ab8dSStefano Zampini 662*3201ab8dSStefano Zampini .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNewtonTRSetTolerances()`, 663f6dfbefdSBarry Smith `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`, 664f6dfbefdSBarry Smith `SNESNewtonTRDCGetRhoFlag()`, `SNESNewtonTRDCSetPreCheck()` 66541ba4c6cSHeeho Park M*/ 666d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTRDC(SNES snes) 667d71ae5a4SJacob Faibussowitsch { 66841ba4c6cSHeeho Park SNES_NEWTONTRDC *neP; 66941ba4c6cSHeeho Park 67041ba4c6cSHeeho Park PetscFunctionBegin; 67141ba4c6cSHeeho Park snes->ops->setup = SNESSetUp_NEWTONTRDC; 67241ba4c6cSHeeho Park snes->ops->solve = SNESSolve_NEWTONTRDC; 67341ba4c6cSHeeho Park snes->ops->destroy = SNESDestroy_NEWTONTRDC; 67441ba4c6cSHeeho Park snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTRDC; 67541ba4c6cSHeeho Park snes->ops->view = SNESView_NEWTONTRDC; 67641ba4c6cSHeeho Park snes->ops->reset = SNESReset_NEWTONTRDC; 67741ba4c6cSHeeho Park 67841ba4c6cSHeeho Park snes->usesksp = PETSC_TRUE; 67941ba4c6cSHeeho Park snes->usesnpc = PETSC_FALSE; 68041ba4c6cSHeeho Park 68141ba4c6cSHeeho Park snes->alwayscomputesfinalresidual = PETSC_TRUE; 68241ba4c6cSHeeho Park 68377e5a1f9SBarry Smith PetscCall(SNESParametersInitialize(snes)); 68477e5a1f9SBarry Smith 6854dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&neP)); 68641ba4c6cSHeeho Park snes->data = (void *)neP; 68741ba4c6cSHeeho Park neP->eta1 = 0.001; 68841ba4c6cSHeeho Park neP->eta2 = 0.25; 68941ba4c6cSHeeho Park neP->eta3 = 0.75; 69041ba4c6cSHeeho Park neP->t1 = 0.25; 69141ba4c6cSHeeho Park neP->t2 = 2.0; 69241ba4c6cSHeeho Park neP->sigma = 0.0001; 69341ba4c6cSHeeho Park neP->itflag = PETSC_FALSE; 69441ba4c6cSHeeho Park neP->rnorm0 = 0.0; 69541ba4c6cSHeeho Park neP->ttol = 0.0; 69641ba4c6cSHeeho Park neP->use_cauchy = PETSC_TRUE; 69741ba4c6cSHeeho Park neP->auto_scale_multiphase = PETSC_FALSE; 69841ba4c6cSHeeho Park neP->auto_scale_max = -1.0; 69941ba4c6cSHeeho Park neP->rho_satisfied = PETSC_FALSE; 700*3201ab8dSStefano Zampini neP->delta = 0.0; 701*3201ab8dSStefano Zampini neP->deltaM = 0.5; 702*3201ab8dSStefano Zampini neP->delta0 = 0.1; 703*3201ab8dSStefano Zampini neP->deltatol = 1.e-12; 70441ba4c6cSHeeho Park 70541ba4c6cSHeeho Park /* for multiphase (multivariable) scaling */ 70641ba4c6cSHeeho Park /* may be used for dynamic allocation of inorms, but it fails snes_tutorials-ex3_13 70741ba4c6cSHeeho Park on test forced DIVERGED_JACOBIAN_DOMAIN test. I will use static array for now. 7089566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(snes->work[0],&neP->bs)); 7099566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(neP->bs,&neP->inorms)); 71041ba4c6cSHeeho Park */ 711*3201ab8dSStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonTRSetTolerances_C", SNESNewtonTRSetTolerances_TRDC)); 7123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 71341ba4c6cSHeeho Park } 714