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 15d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTRDC_KSPConverged_Private(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *cctx) 16d71ae5a4SJacob Faibussowitsch { 1741ba4c6cSHeeho Park SNES_TRDC_KSPConverged_Ctx *ctx = (SNES_TRDC_KSPConverged_Ctx *)cctx; 1841ba4c6cSHeeho Park SNES snes = ctx->snes; 1941ba4c6cSHeeho Park SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC *)snes->data; 2041ba4c6cSHeeho Park Vec x; 2141ba4c6cSHeeho Park PetscReal nrm; 2241ba4c6cSHeeho Park 2341ba4c6cSHeeho Park PetscFunctionBegin; 249566063dSJacob Faibussowitsch PetscCall((*ctx->convtest)(ksp, n, rnorm, reason, ctx->convctx)); 2548a46eb9SPierre Jolivet if (*reason) PetscCall(PetscInfo(snes, "Default or user provided convergence test KSP iterations=%" PetscInt_FMT ", rnorm=%g\n", n, (double)rnorm)); 2641ba4c6cSHeeho Park /* Determine norm of solution */ 279566063dSJacob Faibussowitsch PetscCall(KSPBuildSolution(ksp, NULL, &x)); 289566063dSJacob Faibussowitsch PetscCall(VecNorm(x, NORM_2, &nrm)); 2941ba4c6cSHeeho Park if (nrm >= neP->delta) { 309566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Ending linear iteration early, delta=%g, length=%g\n", (double)neP->delta, (double)nrm)); 3141ba4c6cSHeeho Park *reason = KSP_CONVERGED_STEP_LENGTH; 3241ba4c6cSHeeho Park } 333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3441ba4c6cSHeeho Park } 3541ba4c6cSHeeho Park 36d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTRDC_KSPConverged_Destroy(void *cctx) 37d71ae5a4SJacob Faibussowitsch { 3841ba4c6cSHeeho Park SNES_TRDC_KSPConverged_Ctx *ctx = (SNES_TRDC_KSPConverged_Ctx *)cctx; 3941ba4c6cSHeeho Park 4041ba4c6cSHeeho Park PetscFunctionBegin; 419566063dSJacob Faibussowitsch PetscCall((*ctx->convdestroy)(ctx->convctx)); 429566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 4341ba4c6cSHeeho Park 443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4541ba4c6cSHeeho Park } 4641ba4c6cSHeeho Park 4741ba4c6cSHeeho Park /* 4841ba4c6cSHeeho Park SNESTRDC_Converged_Private -test convergence JUST for 4941ba4c6cSHeeho Park the trust region tolerance. 5041ba4c6cSHeeho Park 5141ba4c6cSHeeho Park */ 52d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTRDC_Converged_Private(SNES snes, PetscInt it, PetscReal xnorm, PetscReal pnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy) 53d71ae5a4SJacob Faibussowitsch { 5441ba4c6cSHeeho Park SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC *)snes->data; 5541ba4c6cSHeeho Park 5641ba4c6cSHeeho Park PetscFunctionBegin; 5741ba4c6cSHeeho Park *reason = SNES_CONVERGED_ITERATING; 5841ba4c6cSHeeho Park if (neP->delta < xnorm * snes->deltatol) { 599566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Diverged due to too small a trust region %g<%g*%g\n", (double)neP->delta, (double)xnorm, (double)snes->deltatol)); 6041ba4c6cSHeeho Park *reason = SNES_DIVERGED_TR_DELTA; 6141ba4c6cSHeeho Park } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 629566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n", snes->max_funcs)); 6341ba4c6cSHeeho Park *reason = SNES_DIVERGED_FUNCTION_COUNT; 6441ba4c6cSHeeho Park } 653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6641ba4c6cSHeeho Park } 6741ba4c6cSHeeho Park 6841ba4c6cSHeeho Park /*@ 69f6dfbefdSBarry Smith SNESNewtonTRDCGetRhoFlag - Get whether the current solution update is within the trust-region. 7041ba4c6cSHeeho Park 7120f4b53cSBarry Smith Logically Collective 7220f4b53cSBarry Smith 73f6dfbefdSBarry Smith Input Parameter: 7441ba4c6cSHeeho Park . snes - the nonlinear solver object 7541ba4c6cSHeeho Park 76f6dfbefdSBarry Smith Output Parameter: 77420bcc1bSBarry Smith . rho_flag - `PETSC_FALSE` or `PETSC_TRUE` 7841ba4c6cSHeeho Park 7941ba4c6cSHeeho Park Level: developer 8041ba4c6cSHeeho Park 81420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCSetPreCheck()`, 82f6dfbefdSBarry Smith `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()` 8341ba4c6cSHeeho Park @*/ 84d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRDCGetRhoFlag(SNES snes, PetscBool *rho_flag) 85d71ae5a4SJacob Faibussowitsch { 8641ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 8741ba4c6cSHeeho Park 8841ba4c6cSHeeho Park PetscFunctionBegin; 8941ba4c6cSHeeho Park PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 904f572ea9SToby Isaac PetscAssertPointer(rho_flag, 2); 9141ba4c6cSHeeho Park *rho_flag = tr->rho_satisfied; 923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9341ba4c6cSHeeho Park } 9441ba4c6cSHeeho Park 9541ba4c6cSHeeho Park /*@C 9641ba4c6cSHeeho Park SNESNewtonTRDCSetPreCheck - Sets a user function that is called before the search step has been determined. 9741ba4c6cSHeeho Park Allows the user a chance to change or override the trust region decision. 9841ba4c6cSHeeho Park 99c3339decSBarry Smith Logically Collective 10041ba4c6cSHeeho Park 10141ba4c6cSHeeho Park Input Parameters: 10241ba4c6cSHeeho Park + snes - the nonlinear solver object 10320f4b53cSBarry Smith . func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRDCPreCheck()` 10420f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 10541ba4c6cSHeeho Park 10641ba4c6cSHeeho Park Level: intermediate 10741ba4c6cSHeeho Park 108f6dfbefdSBarry Smith Note: 109f6dfbefdSBarry Smith This function is called BEFORE the function evaluation within the `SNESNEWTONTRDC` solver. 11041ba4c6cSHeeho Park 111420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`, 112f6dfbefdSBarry Smith `SNESNewtonTRDCGetRhoFlag()` 11341ba4c6cSHeeho Park @*/ 114d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRDCSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, PetscBool *, void *), void *ctx) 115d71ae5a4SJacob Faibussowitsch { 11641ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 11741ba4c6cSHeeho Park 11841ba4c6cSHeeho Park PetscFunctionBegin; 11941ba4c6cSHeeho Park PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 12041ba4c6cSHeeho Park if (func) tr->precheck = func; 12141ba4c6cSHeeho Park if (ctx) tr->precheckctx = ctx; 1223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12341ba4c6cSHeeho Park } 12441ba4c6cSHeeho Park 12541ba4c6cSHeeho Park /*@C 126420bcc1bSBarry Smith SNESNewtonTRDCGetPreCheck - Gets the pre-check function optionally set with `SNESNewtonTRDCSetPreCheck()` 12741ba4c6cSHeeho Park 12820f4b53cSBarry Smith Not Collective 12941ba4c6cSHeeho Park 13041ba4c6cSHeeho Park Input Parameter: 13141ba4c6cSHeeho Park . snes - the nonlinear solver context 13241ba4c6cSHeeho Park 13341ba4c6cSHeeho Park Output Parameters: 13420f4b53cSBarry Smith + func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRDCPreCheck()` 13520f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 13641ba4c6cSHeeho Park 13741ba4c6cSHeeho Park Level: intermediate 13841ba4c6cSHeeho Park 139420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCPreCheck()` 14041ba4c6cSHeeho Park @*/ 141d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRDCGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, PetscBool *, void *), void **ctx) 142d71ae5a4SJacob Faibussowitsch { 14341ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 14441ba4c6cSHeeho Park 14541ba4c6cSHeeho Park PetscFunctionBegin; 14641ba4c6cSHeeho Park PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 14741ba4c6cSHeeho Park if (func) *func = tr->precheck; 14841ba4c6cSHeeho Park if (ctx) *ctx = tr->precheckctx; 1493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15041ba4c6cSHeeho Park } 15141ba4c6cSHeeho Park 15241ba4c6cSHeeho Park /*@C 15341ba4c6cSHeeho Park SNESNewtonTRDCSetPostCheck - Sets a user function that is called after the search step has been determined but before the next 15441ba4c6cSHeeho Park function evaluation. Allows the user a chance to change or override the decision of the line search routine 15541ba4c6cSHeeho Park 156c3339decSBarry Smith Logically Collective 15741ba4c6cSHeeho Park 15841ba4c6cSHeeho Park Input Parameters: 15941ba4c6cSHeeho Park + snes - the nonlinear solver object 16020f4b53cSBarry Smith . func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRDCPostCheck()` 16120f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 16241ba4c6cSHeeho Park 16341ba4c6cSHeeho Park Level: intermediate 16441ba4c6cSHeeho Park 165f6dfbefdSBarry Smith Note: 166f6dfbefdSBarry Smith This function is called BEFORE the function evaluation within the `SNESNEWTONTRDC` solver while the function set in 167f6dfbefdSBarry Smith `SNESLineSearchSetPostCheck()` is called AFTER the function evaluation. 16841ba4c6cSHeeho Park 169420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCPostCheck()`, `SNESNewtonTRDCGetPostCheck()`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCGetPreCheck()` 17041ba4c6cSHeeho Park @*/ 171d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRDCSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx) 172d71ae5a4SJacob Faibussowitsch { 17341ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 17441ba4c6cSHeeho Park 17541ba4c6cSHeeho Park PetscFunctionBegin; 17641ba4c6cSHeeho Park PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 17741ba4c6cSHeeho Park if (func) tr->postcheck = func; 17841ba4c6cSHeeho Park if (ctx) tr->postcheckctx = ctx; 1793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 18041ba4c6cSHeeho Park } 18141ba4c6cSHeeho Park 18241ba4c6cSHeeho Park /*@C 183420bcc1bSBarry Smith SNESNewtonTRDCGetPostCheck - Gets the post-check function optionally set with `SNESNewtonTRDCSetPostCheck()` 18441ba4c6cSHeeho Park 18520f4b53cSBarry Smith Not Collective 18641ba4c6cSHeeho Park 18741ba4c6cSHeeho Park Input Parameter: 18841ba4c6cSHeeho Park . snes - the nonlinear solver context 18941ba4c6cSHeeho Park 19041ba4c6cSHeeho Park Output Parameters: 19120f4b53cSBarry Smith + func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRDCPostCheck()` 19220f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 19341ba4c6cSHeeho Park 19441ba4c6cSHeeho Park Level: intermediate 19541ba4c6cSHeeho Park 196420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCPostCheck()`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCGetPreCheck()` 19741ba4c6cSHeeho Park @*/ 198d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRDCGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) 199d71ae5a4SJacob Faibussowitsch { 20041ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 20141ba4c6cSHeeho Park 20241ba4c6cSHeeho Park PetscFunctionBegin; 20341ba4c6cSHeeho Park PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 20441ba4c6cSHeeho Park if (func) *func = tr->postcheck; 20541ba4c6cSHeeho Park if (ctx) *ctx = tr->postcheckctx; 2063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 20741ba4c6cSHeeho Park } 20841ba4c6cSHeeho Park 209ceaaa498SBarry Smith // PetscClangLinter pragma disable: -fdoc-internal-linkage 21041ba4c6cSHeeho Park /*@C 211f6dfbefdSBarry Smith SNESNewtonTRDCPreCheck - Called before the step has been determined in `SNESNEWTONTRDC` 21241ba4c6cSHeeho Park 213c3339decSBarry Smith Logically Collective 21441ba4c6cSHeeho Park 21541ba4c6cSHeeho Park Input Parameters: 21641ba4c6cSHeeho Park + snes - the solver 21741ba4c6cSHeeho Park . X - The last solution 21841ba4c6cSHeeho Park - Y - The step direction 21941ba4c6cSHeeho Park 2202fe279fdSBarry Smith Output Parameter: 22120f4b53cSBarry Smith . changed_Y - Indicator that the step direction `Y` has been changed. 22241ba4c6cSHeeho Park 22341ba4c6cSHeeho Park Level: developer 22441ba4c6cSHeeho Park 225420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCPostCheck()` 22641ba4c6cSHeeho Park @*/ 227d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNewtonTRDCPreCheck(SNES snes, Vec X, Vec Y, PetscBool *changed_Y) 228d71ae5a4SJacob Faibussowitsch { 22941ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 23041ba4c6cSHeeho Park 23141ba4c6cSHeeho Park PetscFunctionBegin; 23241ba4c6cSHeeho Park *changed_Y = PETSC_FALSE; 23341ba4c6cSHeeho Park if (tr->precheck) { 2349566063dSJacob Faibussowitsch PetscCall((*tr->precheck)(snes, X, Y, changed_Y, tr->precheckctx)); 23541ba4c6cSHeeho Park PetscValidLogicalCollectiveBool(snes, *changed_Y, 4); 23641ba4c6cSHeeho Park } 2373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 23841ba4c6cSHeeho Park } 23941ba4c6cSHeeho Park 240ceaaa498SBarry Smith // PetscClangLinter pragma disable: -fdoc-internal-linkage 24141ba4c6cSHeeho Park /*@C 242f6dfbefdSBarry Smith SNESNewtonTRDCPostCheck - Called after the step has been determined in `SNESNEWTONTRDC` but before the function evaluation at that step 24341ba4c6cSHeeho Park 244c3339decSBarry Smith Logically Collective 24541ba4c6cSHeeho Park 24641ba4c6cSHeeho Park Input Parameters: 247f1a722f8SMatthew G. Knepley + snes - the solver 248f1a722f8SMatthew G. Knepley . X - The last solution 24941ba4c6cSHeeho Park . Y - The full step direction 25041ba4c6cSHeeho Park - W - The updated solution, W = X - Y 25141ba4c6cSHeeho Park 25241ba4c6cSHeeho Park Output Parameters: 25341ba4c6cSHeeho Park + changed_Y - indicator if step has been changed 25420f4b53cSBarry Smith - changed_W - Indicator if the new candidate solution `W` has been changed. 25541ba4c6cSHeeho Park 2562fe279fdSBarry Smith Level: developer 2572fe279fdSBarry Smith 258f6dfbefdSBarry Smith Note: 25920f4b53cSBarry Smith If `Y` is changed then `W` is recomputed as `X` - `Y` 26041ba4c6cSHeeho Park 261420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`, `SNESNewtonTRDCPreCheck() 26241ba4c6cSHeeho Park @*/ 263d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNewtonTRDCPostCheck(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) 264d71ae5a4SJacob Faibussowitsch { 26541ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 26641ba4c6cSHeeho Park 26741ba4c6cSHeeho Park PetscFunctionBegin; 26841ba4c6cSHeeho Park *changed_Y = PETSC_FALSE; 26941ba4c6cSHeeho Park *changed_W = PETSC_FALSE; 27041ba4c6cSHeeho Park if (tr->postcheck) { 2719566063dSJacob Faibussowitsch PetscCall((*tr->postcheck)(snes, X, Y, W, changed_Y, changed_W, tr->postcheckctx)); 27241ba4c6cSHeeho Park PetscValidLogicalCollectiveBool(snes, *changed_Y, 5); 27341ba4c6cSHeeho Park PetscValidLogicalCollectiveBool(snes, *changed_W, 6); 27441ba4c6cSHeeho Park } 2753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 27641ba4c6cSHeeho Park } 27741ba4c6cSHeeho Park 27841ba4c6cSHeeho Park /* 27941ba4c6cSHeeho Park SNESSolve_NEWTONTRDC - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy 28041ba4c6cSHeeho Park (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of 28141ba4c6cSHeeho Park nonlinear equations 28241ba4c6cSHeeho Park 28341ba4c6cSHeeho Park */ 284d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NEWTONTRDC(SNES snes) 285d71ae5a4SJacob Faibussowitsch { 28641ba4c6cSHeeho Park SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC *)snes->data; 28741ba4c6cSHeeho Park Vec X, F, Y, G, W, GradF, YNtmp; 28841ba4c6cSHeeho Park Vec YCtmp; 28941ba4c6cSHeeho Park Mat jac; 29041ba4c6cSHeeho Park PetscInt maxits, i, j, lits, inner_count, bs; 29141ba4c6cSHeeho Park PetscReal rho, fnorm, gnorm, xnorm = 0, delta, ynorm, temp_xnorm, temp_ynorm; /* TRDC inner iteration */ 29241ba4c6cSHeeho Park PetscReal inorms[99]; /* need to make it dynamic eventually, fixed max block size of 99 for now */ 29341ba4c6cSHeeho Park PetscReal deltaM, ynnorm, f0, mp, gTy, g, yTHy; /* rho calculation */ 29441ba4c6cSHeeho Park PetscReal auk, gfnorm, ycnorm, c0, c1, c2, tau, tau_pos, tau_neg, gTBg; /* Cauchy Point */ 29541ba4c6cSHeeho Park KSP ksp; 29641ba4c6cSHeeho Park SNESConvergedReason reason = SNES_CONVERGED_ITERATING; 29741ba4c6cSHeeho Park PetscBool breakout = PETSC_FALSE; 29841ba4c6cSHeeho Park SNES_TRDC_KSPConverged_Ctx *ctx; 29941ba4c6cSHeeho Park PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), (*convdestroy)(void *); 30041ba4c6cSHeeho Park void *convctx; 30141ba4c6cSHeeho Park 30241ba4c6cSHeeho Park PetscFunctionBegin; 30341ba4c6cSHeeho Park maxits = snes->max_its; /* maximum number of iterations */ 30441ba4c6cSHeeho Park X = snes->vec_sol; /* solution vector */ 30541ba4c6cSHeeho Park F = snes->vec_func; /* residual vector */ 30641ba4c6cSHeeho Park Y = snes->work[0]; /* update vector */ 30741ba4c6cSHeeho Park G = snes->work[1]; /* updated residual */ 30841ba4c6cSHeeho Park W = snes->work[2]; /* temporary vector */ 30941ba4c6cSHeeho Park GradF = snes->work[3]; /* grad f = J^T F */ 31041ba4c6cSHeeho Park YNtmp = snes->work[4]; /* Newton solution */ 31141ba4c6cSHeeho Park YCtmp = snes->work[5]; /* Cauchy solution */ 31241ba4c6cSHeeho Park 3130b121fc5SBarry 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); 31441ba4c6cSHeeho Park 3159566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(YNtmp, &bs)); 31641ba4c6cSHeeho Park 3179566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 31841ba4c6cSHeeho Park snes->iter = 0; 3199566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 32041ba4c6cSHeeho Park 32141ba4c6cSHeeho Park /* Set the linear stopping criteria to use the More' trick. From tr.c */ 3229566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 3239566063dSJacob Faibussowitsch PetscCall(KSPGetConvergenceTest(ksp, &convtest, &convctx, &convdestroy)); 32441ba4c6cSHeeho Park if (convtest != SNESTRDC_KSPConverged_Private) { 3259566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 32641ba4c6cSHeeho Park ctx->snes = snes; 3279566063dSJacob Faibussowitsch PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 3289566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(ksp, SNESTRDC_KSPConverged_Private, ctx, SNESTRDC_KSPConverged_Destroy)); 3299566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTRDC_KSPConverged_Private\n")); 33041ba4c6cSHeeho Park } 33141ba4c6cSHeeho Park 33241ba4c6cSHeeho Park if (!snes->vec_func_init_set) { 3339566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */ 33441ba4c6cSHeeho Park } else snes->vec_func_init_set = PETSC_FALSE; 33541ba4c6cSHeeho Park 3369566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */ 33741ba4c6cSHeeho Park SNESCheckFunctionNorm(snes, fnorm); 3389566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */ 33941ba4c6cSHeeho Park 3409566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 34141ba4c6cSHeeho Park snes->norm = fnorm; 3429566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 34341ba4c6cSHeeho Park delta = xnorm ? neP->delta0 * xnorm : neP->delta0; /* initial trust region size scaled by xnorm */ 34441ba4c6cSHeeho Park deltaM = xnorm ? neP->deltaM * xnorm : neP->deltaM; /* maximum trust region size scaled by xnorm */ 34541ba4c6cSHeeho Park neP->delta = delta; 3469566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 3479566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, fnorm)); 34841ba4c6cSHeeho Park 34941ba4c6cSHeeho Park neP->rho_satisfied = PETSC_FALSE; 35041ba4c6cSHeeho Park 35141ba4c6cSHeeho Park /* test convergence */ 352dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 3533ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 35441ba4c6cSHeeho Park 35541ba4c6cSHeeho Park for (i = 0; i < maxits; i++) { 35641ba4c6cSHeeho Park PetscBool changed_y; 35741ba4c6cSHeeho Park PetscBool changed_w; 35841ba4c6cSHeeho Park 35941ba4c6cSHeeho Park /* dogleg method */ 3609566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 36141ba4c6cSHeeho Park SNESCheckJacobianDomainerror(snes); 3629566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian)); 3639566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, F, YNtmp)); /* Quasi Newton Solution */ 36441ba4c6cSHeeho Park SNESCheckKSPSolve(snes); /* this is necessary but old tr.c did not have it*/ 3659566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 3669566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(snes, &jac, NULL, NULL, NULL)); 36741ba4c6cSHeeho Park 36841ba4c6cSHeeho Park /* rescale Jacobian, Newton solution update, and re-calculate delta for multiphase (multivariable) 36941ba4c6cSHeeho Park for inner iteration and Cauchy direction calculation 37041ba4c6cSHeeho Park */ 37141ba4c6cSHeeho Park if (bs > 1 && neP->auto_scale_multiphase) { 3729566063dSJacob Faibussowitsch PetscCall(VecStrideNormAll(YNtmp, NORM_INFINITY, inorms)); 37341ba4c6cSHeeho Park for (j = 0; j < bs; j++) { 37441ba4c6cSHeeho Park if (neP->auto_scale_max > 1.0) { 375ad540459SPierre Jolivet if (inorms[j] < 1.0 / neP->auto_scale_max) inorms[j] = 1.0 / neP->auto_scale_max; 37641ba4c6cSHeeho Park } 3779566063dSJacob Faibussowitsch PetscCall(VecStrideSet(W, j, inorms[j])); 3789566063dSJacob Faibussowitsch PetscCall(VecStrideScale(YNtmp, j, 1.0 / inorms[j])); 3799566063dSJacob Faibussowitsch PetscCall(VecStrideScale(X, j, 1.0 / inorms[j])); 38041ba4c6cSHeeho Park } 3819566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, &xnorm)); 38241ba4c6cSHeeho Park if (i == 0) { 38341ba4c6cSHeeho Park delta = neP->delta0 * xnorm; 38441ba4c6cSHeeho Park } else { 38541ba4c6cSHeeho Park delta = neP->delta * xnorm; 38641ba4c6cSHeeho Park } 38741ba4c6cSHeeho Park deltaM = neP->deltaM * xnorm; 388f3fa974cSJacob Faibussowitsch PetscCall(MatDiagonalScale(jac, NULL, W)); 38941ba4c6cSHeeho Park } 39041ba4c6cSHeeho Park 39141ba4c6cSHeeho Park /* calculating GradF of minimization function */ 3929566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(jac, F, GradF)); /* grad f = J^T F */ 3939566063dSJacob Faibussowitsch PetscCall(VecNorm(YNtmp, NORM_2, &ynnorm)); /* ynnorm <- || Y_newton || */ 39441ba4c6cSHeeho Park 39541ba4c6cSHeeho Park inner_count = 0; 39641ba4c6cSHeeho Park neP->rho_satisfied = PETSC_FALSE; 39741ba4c6cSHeeho Park while (1) { 39841ba4c6cSHeeho Park if (ynnorm <= delta) { /* see if the Newton solution is within the trust region */ 3999566063dSJacob Faibussowitsch PetscCall(VecCopy(YNtmp, Y)); 40041ba4c6cSHeeho Park } else if (neP->use_cauchy) { /* use Cauchy direction if enabled */ 4019566063dSJacob Faibussowitsch PetscCall(MatMult(jac, GradF, W)); 4029566063dSJacob Faibussowitsch PetscCall(VecDotRealPart(W, W, &gTBg)); /* completes GradF^T J^T J GradF */ 4039566063dSJacob Faibussowitsch PetscCall(VecNorm(GradF, NORM_2, &gfnorm)); /* grad f norm <- || grad f || */ 40441ba4c6cSHeeho Park if (gTBg <= 0.0) { 40541ba4c6cSHeeho Park auk = PETSC_MAX_REAL; 40641ba4c6cSHeeho Park } else { 40741ba4c6cSHeeho Park auk = PetscSqr(gfnorm) / gTBg; 40841ba4c6cSHeeho Park } 40941ba4c6cSHeeho Park auk = PetscMin(delta / gfnorm, auk); 4109566063dSJacob Faibussowitsch PetscCall(VecCopy(GradF, YCtmp)); /* this could be improved */ 4119566063dSJacob Faibussowitsch PetscCall(VecScale(YCtmp, auk)); /* YCtmp, Cauchy solution*/ 4129566063dSJacob Faibussowitsch PetscCall(VecNorm(YCtmp, NORM_2, &ycnorm)); /* ycnorm <- || Y_cauchy || */ 41341ba4c6cSHeeho Park if (ycnorm >= delta) { /* see if the Cauchy solution meets the criteria */ 4149566063dSJacob Faibussowitsch PetscCall(VecCopy(YCtmp, Y)); 4159566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "DL evaluated. delta: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n", (double)delta, (double)ynnorm, (double)ycnorm)); 41641ba4c6cSHeeho Park } else { /* take ratio, tau, of Cauchy and Newton direction and step */ 4179566063dSJacob Faibussowitsch PetscCall(VecAXPY(YNtmp, -1.0, YCtmp)); /* YCtmp = A, YNtmp = B */ 4189566063dSJacob Faibussowitsch PetscCall(VecNorm(YNtmp, NORM_2, &c0)); /* this could be improved */ 41941ba4c6cSHeeho Park c0 = PetscSqr(c0); 4209566063dSJacob Faibussowitsch PetscCall(VecDotRealPart(YCtmp, YNtmp, &c1)); 42141ba4c6cSHeeho Park c1 = 2.0 * c1; 4229566063dSJacob Faibussowitsch PetscCall(VecNorm(YCtmp, NORM_2, &c2)); /* this could be improved */ 42341ba4c6cSHeeho Park c2 = PetscSqr(c2) - PetscSqr(delta); 42441ba4c6cSHeeho Park tau_pos = (c1 + PetscSqrtReal(PetscSqr(c1) - 4. * c0 * c2)) / (2. * c0); /* quadratic formula */ 42541ba4c6cSHeeho Park tau_neg = (c1 - PetscSqrtReal(PetscSqr(c1) - 4. * c0 * c2)) / (2. * c0); 42641ba4c6cSHeeho Park tau = PetscMax(tau_pos, tau_neg); /* can tau_neg > tau_pos? I don't think so, but just in case. */ 4279566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "DL evaluated. tau: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n", (double)tau, (double)ynnorm, (double)ycnorm)); 4289566063dSJacob Faibussowitsch PetscCall(VecWAXPY(W, tau, YNtmp, YCtmp)); 4299566063dSJacob Faibussowitsch PetscCall(VecAXPY(W, -tau, YCtmp)); 4309566063dSJacob Faibussowitsch PetscCall(VecCopy(W, Y)); /* this could be improved */ 43141ba4c6cSHeeho Park } 43241ba4c6cSHeeho Park } else { 43341ba4c6cSHeeho Park /* if Cauchy is disabled, only use Newton direction */ 43441ba4c6cSHeeho Park auk = delta / ynnorm; 4359566063dSJacob Faibussowitsch PetscCall(VecScale(YNtmp, auk)); 4369566063dSJacob Faibussowitsch PetscCall(VecCopy(YNtmp, Y)); /* this could be improved (many VecCopy, VecNorm)*/ 43741ba4c6cSHeeho Park } 43841ba4c6cSHeeho Park 4399566063dSJacob Faibussowitsch PetscCall(VecNorm(Y, NORM_2, &ynorm)); /* compute the final ynorm */ 44041ba4c6cSHeeho Park f0 = 0.5 * PetscSqr(fnorm); /* minimizing function f(X) */ 4419566063dSJacob Faibussowitsch PetscCall(MatMult(jac, Y, W)); 4429566063dSJacob Faibussowitsch PetscCall(VecDotRealPart(W, W, &yTHy)); /* completes GradY^T J^T J GradY */ 4439566063dSJacob Faibussowitsch PetscCall(VecDotRealPart(GradF, Y, &gTy)); 44441ba4c6cSHeeho Park mp = f0 - gTy + 0.5 * yTHy; /* quadratic model to satisfy, -gTy because our update is X-Y*/ 44541ba4c6cSHeeho Park 44641ba4c6cSHeeho Park /* scale back solution update */ 44741ba4c6cSHeeho Park if (bs > 1 && neP->auto_scale_multiphase) { 44841ba4c6cSHeeho Park for (j = 0; j < bs; j++) { 4499566063dSJacob Faibussowitsch PetscCall(VecStrideScale(Y, j, inorms[j])); 45041ba4c6cSHeeho Park if (inner_count == 0) { 45141ba4c6cSHeeho Park /* TRDC inner algorithm does not need scaled X after calculating delta in the outer iteration */ 45241ba4c6cSHeeho Park /* need to scale back X to match Y and provide proper update to the external code */ 4539566063dSJacob Faibussowitsch PetscCall(VecStrideScale(X, j, inorms[j])); 45441ba4c6cSHeeho Park } 45541ba4c6cSHeeho Park } 4569566063dSJacob Faibussowitsch if (inner_count == 0) PetscCall(VecNorm(X, NORM_2, &temp_xnorm)); /* only in the first iteration */ 4579566063dSJacob Faibussowitsch PetscCall(VecNorm(Y, NORM_2, &temp_ynorm)); 45841ba4c6cSHeeho Park } else { 45941ba4c6cSHeeho Park temp_xnorm = xnorm; 46041ba4c6cSHeeho Park temp_ynorm = ynorm; 46141ba4c6cSHeeho Park } 46241ba4c6cSHeeho Park inner_count++; 46341ba4c6cSHeeho Park 46441ba4c6cSHeeho Park /* Evaluate the solution to meet the improvement ratio criteria */ 4659566063dSJacob Faibussowitsch PetscCall(SNESNewtonTRDCPreCheck(snes, X, Y, &changed_y)); 4669566063dSJacob Faibussowitsch PetscCall(VecWAXPY(W, -1.0, Y, X)); 4679566063dSJacob Faibussowitsch PetscCall(SNESNewtonTRDCPostCheck(snes, X, Y, W, &changed_y, &changed_w)); 4689566063dSJacob Faibussowitsch if (changed_y) PetscCall(VecWAXPY(W, -1.0, Y, X)); 4699566063dSJacob Faibussowitsch PetscCall(VecCopy(Y, snes->vec_sol_update)); 4709566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, W, G)); /* F(X-Y) = G */ 4719566063dSJacob Faibussowitsch PetscCall(VecNorm(G, NORM_2, &gnorm)); /* gnorm <- || g || */ 47241ba4c6cSHeeho Park SNESCheckFunctionNorm(snes, gnorm); 47341ba4c6cSHeeho Park g = 0.5 * PetscSqr(gnorm); /* minimizing function g(W) */ 47441ba4c6cSHeeho Park if (f0 == mp) rho = 0.0; 47541ba4c6cSHeeho Park else rho = (f0 - g) / (f0 - mp); /* actual improvement over predicted improvement */ 47641ba4c6cSHeeho Park 47741ba4c6cSHeeho Park if (rho < neP->eta2) { 47841ba4c6cSHeeho Park delta *= neP->t1; /* shrink the region */ 47941ba4c6cSHeeho Park } else if (rho > neP->eta3) { 48041ba4c6cSHeeho Park delta = PetscMin(neP->t2 * delta, deltaM); /* expand the region, but not greater than deltaM */ 48141ba4c6cSHeeho Park } 48241ba4c6cSHeeho Park 48341ba4c6cSHeeho Park neP->delta = delta; 48441ba4c6cSHeeho Park if (rho >= neP->eta1) { 48541ba4c6cSHeeho Park /* unscale delta and xnorm before going to the next outer iteration */ 48641ba4c6cSHeeho Park if (bs > 1 && neP->auto_scale_multiphase) { 48741ba4c6cSHeeho Park neP->delta = delta / xnorm; 48841ba4c6cSHeeho Park xnorm = temp_xnorm; 48941ba4c6cSHeeho Park ynorm = temp_ynorm; 49041ba4c6cSHeeho Park } 49141ba4c6cSHeeho Park neP->rho_satisfied = PETSC_TRUE; 49241ba4c6cSHeeho Park break; /* the improvement ratio is satisfactory */ 49341ba4c6cSHeeho Park } 4949566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Trying again in smaller region\n")); 49541ba4c6cSHeeho Park 49641ba4c6cSHeeho Park /* check to see if progress is hopeless */ 49741ba4c6cSHeeho Park neP->itflag = PETSC_FALSE; 49841ba4c6cSHeeho Park /* both delta, ynorm, and xnorm are either scaled or unscaled */ 4999566063dSJacob Faibussowitsch PetscCall(SNESTRDC_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP)); 50041ba4c6cSHeeho Park /* if multiphase state changes, break out inner iteration */ 50141ba4c6cSHeeho Park if (reason == SNES_BREAKOUT_INNER_ITER) { 50241ba4c6cSHeeho Park if (bs > 1 && neP->auto_scale_multiphase) { 50341ba4c6cSHeeho Park /* unscale delta and xnorm before going to the next outer iteration */ 50441ba4c6cSHeeho Park neP->delta = delta / xnorm; 50541ba4c6cSHeeho Park xnorm = temp_xnorm; 50641ba4c6cSHeeho Park ynorm = temp_ynorm; 50741ba4c6cSHeeho Park } 50841ba4c6cSHeeho Park reason = SNES_CONVERGED_ITERATING; 50941ba4c6cSHeeho Park break; 51041ba4c6cSHeeho Park } 51141ba4c6cSHeeho Park if (reason == SNES_CONVERGED_SNORM_RELATIVE) reason = SNES_DIVERGED_INNER; 51241ba4c6cSHeeho Park if (reason) { 51341ba4c6cSHeeho Park if (reason < 0) { 51441ba4c6cSHeeho Park /* We're not progressing, so return with the current iterate */ 5159566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, i + 1, fnorm)); 51641ba4c6cSHeeho Park breakout = PETSC_TRUE; 51741ba4c6cSHeeho Park break; 51841ba4c6cSHeeho Park } else if (reason > 0) { 51941ba4c6cSHeeho Park /* We're converged, so return with the current iterate and update solution */ 5209566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, i + 1, fnorm)); 52141ba4c6cSHeeho Park breakout = PETSC_FALSE; 52241ba4c6cSHeeho Park break; 52341ba4c6cSHeeho Park } 52441ba4c6cSHeeho Park } 52541ba4c6cSHeeho Park snes->numFailures++; 52641ba4c6cSHeeho Park } 52741ba4c6cSHeeho Park if (!breakout) { 52841ba4c6cSHeeho Park /* Update function and solution vectors */ 52941ba4c6cSHeeho Park fnorm = gnorm; 5309566063dSJacob Faibussowitsch PetscCall(VecCopy(G, F)); 5319566063dSJacob Faibussowitsch PetscCall(VecCopy(W, X)); 53241ba4c6cSHeeho Park /* Monitor convergence */ 5339566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 53441ba4c6cSHeeho Park snes->iter = i + 1; 53541ba4c6cSHeeho Park snes->norm = fnorm; 53641ba4c6cSHeeho Park snes->xnorm = xnorm; 53741ba4c6cSHeeho Park snes->ynorm = ynorm; 5389566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 5399566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 5409566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 54141ba4c6cSHeeho Park /* Test for convergence, xnorm = || X || */ 54241ba4c6cSHeeho Park neP->itflag = PETSC_TRUE; 5439566063dSJacob Faibussowitsch if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm)); 544dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP); 54541ba4c6cSHeeho Park if (reason) break; 54641ba4c6cSHeeho Park } else break; 54741ba4c6cSHeeho Park } 54841ba4c6cSHeeho Park 5499566063dSJacob Faibussowitsch /* PetscCall(PetscFree(inorms)); */ 55041ba4c6cSHeeho Park if (i == maxits) { 5519566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits)); 55241ba4c6cSHeeho Park if (!reason) reason = SNES_DIVERGED_MAX_IT; 55341ba4c6cSHeeho Park } 5549566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 55541ba4c6cSHeeho Park snes->reason = reason; 5569566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 55741ba4c6cSHeeho Park if (convtest != SNESTRDC_KSPConverged_Private) { 5589566063dSJacob Faibussowitsch PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 5599566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 5609566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(ksp, convtest, convctx, convdestroy)); 56141ba4c6cSHeeho Park } 5623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 56341ba4c6cSHeeho Park } 56441ba4c6cSHeeho Park 565d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_NEWTONTRDC(SNES snes) 566d71ae5a4SJacob Faibussowitsch { 56741ba4c6cSHeeho Park PetscFunctionBegin; 5689566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes, 6)); 5699566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 5703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 57141ba4c6cSHeeho Park } 57241ba4c6cSHeeho Park 57366976f2fSJacob Faibussowitsch static PetscErrorCode SNESReset_NEWTONTRDC(SNES snes) 574d71ae5a4SJacob Faibussowitsch { 57541ba4c6cSHeeho Park PetscFunctionBegin; 5763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 57741ba4c6cSHeeho Park } 57841ba4c6cSHeeho Park 579d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_NEWTONTRDC(SNES snes) 580d71ae5a4SJacob Faibussowitsch { 58141ba4c6cSHeeho Park PetscFunctionBegin; 5829566063dSJacob Faibussowitsch PetscCall(SNESReset_NEWTONTRDC(snes)); 5839566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 5843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 58541ba4c6cSHeeho Park } 58641ba4c6cSHeeho Park 587d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NEWTONTRDC(SNES snes, PetscOptionItems *PetscOptionsObject) 588d71ae5a4SJacob Faibussowitsch { 58941ba4c6cSHeeho Park SNES_NEWTONTRDC *ctx = (SNES_NEWTONTRDC *)snes->data; 59041ba4c6cSHeeho Park 59141ba4c6cSHeeho Park PetscFunctionBegin; 592d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations"); 5939566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_tol", "Trust region tolerance", "SNESSetTrustRegionTolerance", snes->deltatol, &snes->deltatol, NULL)); 5949566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_eta1", "eta1", "None", ctx->eta1, &ctx->eta1, NULL)); 5959566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_eta2", "eta2", "None", ctx->eta2, &ctx->eta2, NULL)); 5969566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_eta3", "eta3", "None", ctx->eta3, &ctx->eta3, NULL)); 5979566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_t1", "t1", "None", ctx->t1, &ctx->t1, NULL)); 5989566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_t2", "t2", "None", ctx->t2, &ctx->t2, NULL)); 5999566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_deltaM", "deltaM", "None", ctx->deltaM, &ctx->deltaM, NULL)); 6009566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_delta0", "delta0", "None", ctx->delta0, &ctx->delta0, NULL)); 6019566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_auto_scale_max", "auto_scale_max", "None", ctx->auto_scale_max, &ctx->auto_scale_max, NULL)); 6029566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_trdc_use_cauchy", "use_cauchy", "use Cauchy step and direction", ctx->use_cauchy, &ctx->use_cauchy, NULL)); 6039566063dSJacob 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)); 604d0609cedSBarry Smith PetscOptionsHeadEnd(); 6053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 60641ba4c6cSHeeho Park } 60741ba4c6cSHeeho Park 608d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NEWTONTRDC(SNES snes, PetscViewer viewer) 609d71ae5a4SJacob Faibussowitsch { 61041ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 61141ba4c6cSHeeho Park PetscBool iascii; 61241ba4c6cSHeeho Park 61341ba4c6cSHeeho Park PetscFunctionBegin; 6149566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 61541ba4c6cSHeeho Park if (iascii) { 61663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Trust region tolerance %g (-snes_trtol)\n", (double)snes->deltatol)); 6179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " eta1=%g, eta2=%g, eta3=%g\n", (double)tr->eta1, (double)tr->eta2, (double)tr->eta3)); 6189566063dSJacob 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)); 61941ba4c6cSHeeho Park } 6203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 62141ba4c6cSHeeho Park } 622f6dfbefdSBarry Smith 62341ba4c6cSHeeho Park /*MC 62441ba4c6cSHeeho Park SNESNEWTONTRDC - Newton based nonlinear solver that uses trust-region dogleg method with Cauchy direction 62541ba4c6cSHeeho Park 626f6dfbefdSBarry Smith Options Database Keys: 62741ba4c6cSHeeho Park + -snes_trdc_tol <tol> - trust region tolerance 62841ba4c6cSHeeho Park . -snes_trdc_eta1 <eta1> - trust region parameter 0.0 <= eta1 <= eta2, rho >= eta1 breaks out of the inner iteration (default: eta1=0.001) 62941ba4c6cSHeeho Park . -snes_trdc_eta2 <eta2> - trust region parameter 0.0 <= eta1 <= eta2, rho <= eta2 shrinks the trust region (default: eta2=0.25) 63041ba4c6cSHeeho Park . -snes_trdc_eta3 <eta3> - trust region parameter eta3 > eta2, rho >= eta3 expands the trust region (default: eta3=0.75) 63141ba4c6cSHeeho Park . -snes_trdc_t1 <t1> - trust region parameter, shrinking factor of trust region (default: 0.25) 63241ba4c6cSHeeho Park . -snes_trdc_t2 <t2> - trust region parameter, expanding factor of trust region (default: 2.0) 633*1d27aa22SBarry Smith . -snes_trdc_deltaM <deltaM> - trust region parameter, max size of trust region, $deltaM*norm2(x)$ (default: 0.5) 634*1d27aa22SBarry Smith . -snes_trdc_delta0 <delta0> - trust region parameter, initial size of trust region, $delta0*norm2(x)$ (default: 0.1) 63541ba4c6cSHeeho Park . -snes_trdc_auto_scale_max <auto_scale_max> - used with auto_scale_multiphase, caps the maximum auto-scaling factor 63641ba4c6cSHeeho Park . -snes_trdc_use_cauchy <use_cauchy> - True uses dogleg Cauchy (Steepest Descent direction) step & direction in the trust region algorithm 63741ba4c6cSHeeho Park - -snes_trdc_auto_scale_multiphase <auto_scale_multiphase> - True turns on auto-scaling for multivariable block matrix for Cauchy and trust region 63841ba4c6cSHeeho Park 63920f4b53cSBarry Smith Level: intermediate 64020f4b53cSBarry Smith 641*1d27aa22SBarry Smith Note: 642*1d27aa22SBarry Smith See {cite}`park2021linear` 64341ba4c6cSHeeho Park 644420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESSetTrustRegionTolerance()`, 645f6dfbefdSBarry Smith `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`, 646f6dfbefdSBarry Smith `SNESNewtonTRDCGetRhoFlag()`, `SNESNewtonTRDCSetPreCheck()` 64741ba4c6cSHeeho Park M*/ 648d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTRDC(SNES snes) 649d71ae5a4SJacob Faibussowitsch { 65041ba4c6cSHeeho Park SNES_NEWTONTRDC *neP; 65141ba4c6cSHeeho Park 65241ba4c6cSHeeho Park PetscFunctionBegin; 65341ba4c6cSHeeho Park snes->ops->setup = SNESSetUp_NEWTONTRDC; 65441ba4c6cSHeeho Park snes->ops->solve = SNESSolve_NEWTONTRDC; 65541ba4c6cSHeeho Park snes->ops->destroy = SNESDestroy_NEWTONTRDC; 65641ba4c6cSHeeho Park snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTRDC; 65741ba4c6cSHeeho Park snes->ops->view = SNESView_NEWTONTRDC; 65841ba4c6cSHeeho Park snes->ops->reset = SNESReset_NEWTONTRDC; 65941ba4c6cSHeeho Park 66041ba4c6cSHeeho Park snes->usesksp = PETSC_TRUE; 66141ba4c6cSHeeho Park snes->usesnpc = PETSC_FALSE; 66241ba4c6cSHeeho Park 66341ba4c6cSHeeho Park snes->alwayscomputesfinalresidual = PETSC_TRUE; 66441ba4c6cSHeeho Park 6654dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&neP)); 66641ba4c6cSHeeho Park snes->data = (void *)neP; 66741ba4c6cSHeeho Park neP->delta = 0.0; 66841ba4c6cSHeeho Park neP->delta0 = 0.1; 66941ba4c6cSHeeho Park neP->eta1 = 0.001; 67041ba4c6cSHeeho Park neP->eta2 = 0.25; 67141ba4c6cSHeeho Park neP->eta3 = 0.75; 67241ba4c6cSHeeho Park neP->t1 = 0.25; 67341ba4c6cSHeeho Park neP->t2 = 2.0; 67441ba4c6cSHeeho Park neP->deltaM = 0.5; 67541ba4c6cSHeeho Park neP->sigma = 0.0001; 67641ba4c6cSHeeho Park neP->itflag = PETSC_FALSE; 67741ba4c6cSHeeho Park neP->rnorm0 = 0.0; 67841ba4c6cSHeeho Park neP->ttol = 0.0; 67941ba4c6cSHeeho Park neP->use_cauchy = PETSC_TRUE; 68041ba4c6cSHeeho Park neP->auto_scale_multiphase = PETSC_FALSE; 68141ba4c6cSHeeho Park neP->auto_scale_max = -1.0; 68241ba4c6cSHeeho Park neP->rho_satisfied = PETSC_FALSE; 68341ba4c6cSHeeho Park snes->deltatol = 1.e-12; 68441ba4c6cSHeeho Park 68541ba4c6cSHeeho Park /* for multiphase (multivariable) scaling */ 68641ba4c6cSHeeho Park /* may be used for dynamic allocation of inorms, but it fails snes_tutorials-ex3_13 68741ba4c6cSHeeho Park on test forced DIVERGED_JACOBIAN_DOMAIN test. I will use static array for now. 6889566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(snes->work[0],&neP->bs)); 6899566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(neP->bs,&neP->inorms)); 69041ba4c6cSHeeho Park */ 69141ba4c6cSHeeho Park 6923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 69341ba4c6cSHeeho Park } 694