141ba4c6cSHeeho Park 241ba4c6cSHeeho Park #include <../src/snes/impls/ntrdc/ntrdcimpl.h> /*I "petscsnes.h" I*/ 341ba4c6cSHeeho Park 441ba4c6cSHeeho Park typedef struct { 541ba4c6cSHeeho Park SNES snes; 641ba4c6cSHeeho Park /* Information on the regular SNES convergence test; which may have been user provided 741ba4c6cSHeeho Park Copied from tr.c (maybe able to disposed, but this is a private function) - Heeho 841ba4c6cSHeeho Park Same with SNESTR_KSPConverged_Private, SNESTR_KSPConverged_Destroy, and SNESTR_Converged_Private 941ba4c6cSHeeho Park */ 1041ba4c6cSHeeho Park 1141ba4c6cSHeeho Park PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 1241ba4c6cSHeeho Park PetscErrorCode (*convdestroy)(void *); 1341ba4c6cSHeeho Park void *convctx; 1441ba4c6cSHeeho Park } SNES_TRDC_KSPConverged_Ctx; 1541ba4c6cSHeeho Park 169371c9d4SSatish Balay static PetscErrorCode SNESTRDC_KSPConverged_Private(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *cctx) { 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 } 3341ba4c6cSHeeho Park PetscFunctionReturn(0); 3441ba4c6cSHeeho Park } 3541ba4c6cSHeeho Park 369371c9d4SSatish Balay static PetscErrorCode SNESTRDC_KSPConverged_Destroy(void *cctx) { 3741ba4c6cSHeeho Park SNES_TRDC_KSPConverged_Ctx *ctx = (SNES_TRDC_KSPConverged_Ctx *)cctx; 3841ba4c6cSHeeho Park 3941ba4c6cSHeeho Park PetscFunctionBegin; 409566063dSJacob Faibussowitsch PetscCall((*ctx->convdestroy)(ctx->convctx)); 419566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 4241ba4c6cSHeeho Park 4341ba4c6cSHeeho Park PetscFunctionReturn(0); 4441ba4c6cSHeeho Park } 4541ba4c6cSHeeho Park 4641ba4c6cSHeeho Park /* 4741ba4c6cSHeeho Park SNESTRDC_Converged_Private -test convergence JUST for 4841ba4c6cSHeeho Park the trust region tolerance. 4941ba4c6cSHeeho Park 5041ba4c6cSHeeho Park */ 519371c9d4SSatish Balay static PetscErrorCode SNESTRDC_Converged_Private(SNES snes, PetscInt it, PetscReal xnorm, PetscReal pnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy) { 5241ba4c6cSHeeho Park SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC *)snes->data; 5341ba4c6cSHeeho Park 5441ba4c6cSHeeho Park PetscFunctionBegin; 5541ba4c6cSHeeho Park *reason = SNES_CONVERGED_ITERATING; 5641ba4c6cSHeeho Park if (neP->delta < xnorm * snes->deltatol) { 579566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Diverged due to too small a trust region %g<%g*%g\n", (double)neP->delta, (double)xnorm, (double)snes->deltatol)); 5841ba4c6cSHeeho Park *reason = SNES_DIVERGED_TR_DELTA; 5941ba4c6cSHeeho Park } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 609566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n", snes->max_funcs)); 6141ba4c6cSHeeho Park *reason = SNES_DIVERGED_FUNCTION_COUNT; 6241ba4c6cSHeeho Park } 6341ba4c6cSHeeho Park PetscFunctionReturn(0); 6441ba4c6cSHeeho Park } 6541ba4c6cSHeeho Park 6641ba4c6cSHeeho Park /*@ 67f6dfbefdSBarry Smith SNESNewtonTRDCGetRhoFlag - Get whether the current solution update is within the trust-region. 6841ba4c6cSHeeho Park 69f6dfbefdSBarry Smith Input Parameter: 7041ba4c6cSHeeho Park . snes - the nonlinear solver object 7141ba4c6cSHeeho Park 72f6dfbefdSBarry Smith Output Parameter: 73f6dfbefdSBarry Smith . rho_flag: `PETSC_TRUE` if the solution update is in the trust-region; otherwise, `PETSC_FALSE` 7441ba4c6cSHeeho Park 7541ba4c6cSHeeho Park Level: developer 7641ba4c6cSHeeho Park 77f6dfbefdSBarry Smith .seealso: `SNESNEWTONTRDC`, `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, , `SNESNewtonTRDCSetPreCheck()`, 78f6dfbefdSBarry Smith `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()` 7941ba4c6cSHeeho Park @*/ 809371c9d4SSatish Balay PetscErrorCode SNESNewtonTRDCGetRhoFlag(SNES snes, PetscBool *rho_flag) { 8141ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 8241ba4c6cSHeeho Park 8341ba4c6cSHeeho Park PetscFunctionBegin; 8441ba4c6cSHeeho Park PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 8541ba4c6cSHeeho Park PetscValidBoolPointer(rho_flag, 2); 8641ba4c6cSHeeho Park *rho_flag = tr->rho_satisfied; 8741ba4c6cSHeeho Park PetscFunctionReturn(0); 8841ba4c6cSHeeho Park } 8941ba4c6cSHeeho Park 9041ba4c6cSHeeho Park /*@C 9141ba4c6cSHeeho Park SNESNewtonTRDCSetPreCheck - Sets a user function that is called before the search step has been determined. 9241ba4c6cSHeeho Park Allows the user a chance to change or override the trust region decision. 9341ba4c6cSHeeho Park 9441ba4c6cSHeeho Park Logically Collective on snes 9541ba4c6cSHeeho Park 9641ba4c6cSHeeho Park Input Parameters: 9741ba4c6cSHeeho Park + snes - the nonlinear solver object 98f6dfbefdSBarry Smith . func - [optional] function evaluation routine, see `SNESNewtonTRDCPreCheck()` for the calling sequence 9941ba4c6cSHeeho Park - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 10041ba4c6cSHeeho Park 10141ba4c6cSHeeho Park Level: intermediate 10241ba4c6cSHeeho Park 103f6dfbefdSBarry Smith Note: 104f6dfbefdSBarry Smith This function is called BEFORE the function evaluation within the `SNESNEWTONTRDC` solver. 10541ba4c6cSHeeho Park 106f6dfbefdSBarry Smith .seealso: `SNESNEWTONTRDC`, `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`, 107f6dfbefdSBarry Smith `SNESNewtonTRDCGetRhoFlag()` 10841ba4c6cSHeeho Park @*/ 1099371c9d4SSatish Balay PetscErrorCode SNESNewtonTRDCSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, PetscBool *, void *), void *ctx) { 11041ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 11141ba4c6cSHeeho Park 11241ba4c6cSHeeho Park PetscFunctionBegin; 11341ba4c6cSHeeho Park PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 11441ba4c6cSHeeho Park if (func) tr->precheck = func; 11541ba4c6cSHeeho Park if (ctx) tr->precheckctx = ctx; 11641ba4c6cSHeeho Park PetscFunctionReturn(0); 11741ba4c6cSHeeho Park } 11841ba4c6cSHeeho Park 11941ba4c6cSHeeho Park /*@C 12041ba4c6cSHeeho Park SNESNewtonTRDCGetPreCheck - Gets the pre-check function 12141ba4c6cSHeeho Park 12241ba4c6cSHeeho Park Not collective 12341ba4c6cSHeeho Park 12441ba4c6cSHeeho Park Input Parameter: 12541ba4c6cSHeeho Park . snes - the nonlinear solver context 12641ba4c6cSHeeho Park 12741ba4c6cSHeeho Park Output Parameters: 128f6dfbefdSBarry Smith + func - [optional] function evaluation routine, see for the calling sequence `SNESNewtonTRDCPreCheck()` 12941ba4c6cSHeeho Park - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 13041ba4c6cSHeeho Park 13141ba4c6cSHeeho Park Level: intermediate 13241ba4c6cSHeeho Park 133f6dfbefdSBarry Smith .seealso: `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCPreCheck()` 13441ba4c6cSHeeho Park @*/ 1359371c9d4SSatish Balay PetscErrorCode SNESNewtonTRDCGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, PetscBool *, void *), void **ctx) { 13641ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 13741ba4c6cSHeeho Park 13841ba4c6cSHeeho Park PetscFunctionBegin; 13941ba4c6cSHeeho Park PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 14041ba4c6cSHeeho Park if (func) *func = tr->precheck; 14141ba4c6cSHeeho Park if (ctx) *ctx = tr->precheckctx; 14241ba4c6cSHeeho Park PetscFunctionReturn(0); 14341ba4c6cSHeeho Park } 14441ba4c6cSHeeho Park 14541ba4c6cSHeeho Park /*@C 14641ba4c6cSHeeho Park SNESNewtonTRDCSetPostCheck - Sets a user function that is called after the search step has been determined but before the next 14741ba4c6cSHeeho Park function evaluation. Allows the user a chance to change or override the decision of the line search routine 14841ba4c6cSHeeho Park 14941ba4c6cSHeeho Park Logically Collective on snes 15041ba4c6cSHeeho Park 15141ba4c6cSHeeho Park Input Parameters: 15241ba4c6cSHeeho Park + snes - the nonlinear solver object 153f6dfbefdSBarry Smith . func - [optional] function evaluation routine, see `SNESNewtonTRDCPostCheck()` for the calling sequence 15441ba4c6cSHeeho Park - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 15541ba4c6cSHeeho Park 15641ba4c6cSHeeho Park Level: intermediate 15741ba4c6cSHeeho Park 158f6dfbefdSBarry Smith Note: 159f6dfbefdSBarry Smith This function is called BEFORE the function evaluation within the `SNESNEWTONTRDC` solver while the function set in 160f6dfbefdSBarry Smith `SNESLineSearchSetPostCheck()` is called AFTER the function evaluation. 16141ba4c6cSHeeho Park 162f6dfbefdSBarry Smith .seealso: `SNESNEWTONTRDC`, `SNESNewtonTRDCPostCheck()`, `SNESNewtonTRDCGetPostCheck()`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCGetPreCheck()` 16341ba4c6cSHeeho Park @*/ 1649371c9d4SSatish Balay PetscErrorCode SNESNewtonTRDCSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx) { 16541ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 16641ba4c6cSHeeho Park 16741ba4c6cSHeeho Park PetscFunctionBegin; 16841ba4c6cSHeeho Park PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 16941ba4c6cSHeeho Park if (func) tr->postcheck = func; 17041ba4c6cSHeeho Park if (ctx) tr->postcheckctx = ctx; 17141ba4c6cSHeeho Park PetscFunctionReturn(0); 17241ba4c6cSHeeho Park } 17341ba4c6cSHeeho Park 17441ba4c6cSHeeho Park /*@C 17541ba4c6cSHeeho Park SNESNewtonTRDCGetPostCheck - Gets the post-check function 17641ba4c6cSHeeho Park 17741ba4c6cSHeeho Park Not collective 17841ba4c6cSHeeho Park 17941ba4c6cSHeeho Park Input Parameter: 18041ba4c6cSHeeho Park . snes - the nonlinear solver context 18141ba4c6cSHeeho Park 18241ba4c6cSHeeho Park Output Parameters: 18341ba4c6cSHeeho Park + func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRDCPostCheck() 18441ba4c6cSHeeho Park - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 18541ba4c6cSHeeho Park 18641ba4c6cSHeeho Park Level: intermediate 18741ba4c6cSHeeho Park 188f6dfbefdSBarry Smith .seealso: `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCPostCheck()` 18941ba4c6cSHeeho Park @*/ 1909371c9d4SSatish Balay PetscErrorCode SNESNewtonTRDCGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) { 19141ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 19241ba4c6cSHeeho Park 19341ba4c6cSHeeho Park PetscFunctionBegin; 19441ba4c6cSHeeho Park PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 19541ba4c6cSHeeho Park if (func) *func = tr->postcheck; 19641ba4c6cSHeeho Park if (ctx) *ctx = tr->postcheckctx; 19741ba4c6cSHeeho Park PetscFunctionReturn(0); 19841ba4c6cSHeeho Park } 19941ba4c6cSHeeho Park 20041ba4c6cSHeeho Park /*@C 201f6dfbefdSBarry Smith SNESNewtonTRDCPreCheck - Called before the step has been determined in `SNESNEWTONTRDC` 20241ba4c6cSHeeho Park 20341ba4c6cSHeeho Park Logically Collective on snes 20441ba4c6cSHeeho Park 20541ba4c6cSHeeho Park Input Parameters: 20641ba4c6cSHeeho Park + snes - the solver 20741ba4c6cSHeeho Park . X - The last solution 20841ba4c6cSHeeho Park - Y - The step direction 20941ba4c6cSHeeho Park 21041ba4c6cSHeeho Park Output Parameters: 21141ba4c6cSHeeho Park . changed_Y - Indicator that the step direction Y has been changed. 21241ba4c6cSHeeho Park 21341ba4c6cSHeeho Park Level: developer 21441ba4c6cSHeeho Park 215f6dfbefdSBarry Smith .seealso: `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCPostCheck()` 21641ba4c6cSHeeho Park @*/ 2179371c9d4SSatish Balay static PetscErrorCode SNESNewtonTRDCPreCheck(SNES snes, Vec X, Vec Y, PetscBool *changed_Y) { 21841ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 21941ba4c6cSHeeho Park 22041ba4c6cSHeeho Park PetscFunctionBegin; 22141ba4c6cSHeeho Park *changed_Y = PETSC_FALSE; 22241ba4c6cSHeeho Park if (tr->precheck) { 2239566063dSJacob Faibussowitsch PetscCall((*tr->precheck)(snes, X, Y, changed_Y, tr->precheckctx)); 22441ba4c6cSHeeho Park PetscValidLogicalCollectiveBool(snes, *changed_Y, 4); 22541ba4c6cSHeeho Park } 22641ba4c6cSHeeho Park PetscFunctionReturn(0); 22741ba4c6cSHeeho Park } 22841ba4c6cSHeeho Park 22941ba4c6cSHeeho Park /*@C 230f6dfbefdSBarry Smith SNESNewtonTRDCPostCheck - Called after the step has been determined in `SNESNEWTONTRDC` but before the function evaluation at that step 23141ba4c6cSHeeho Park 23241ba4c6cSHeeho Park Logically Collective on snes 23341ba4c6cSHeeho Park 23441ba4c6cSHeeho Park Input Parameters: 235f1a722f8SMatthew G. Knepley + snes - the solver 236f1a722f8SMatthew G. Knepley . X - The last solution 23741ba4c6cSHeeho Park . Y - The full step direction 23841ba4c6cSHeeho Park - W - The updated solution, W = X - Y 23941ba4c6cSHeeho Park 24041ba4c6cSHeeho Park Output Parameters: 24141ba4c6cSHeeho Park + changed_Y - indicator if step has been changed 24241ba4c6cSHeeho Park - changed_W - Indicator if the new candidate solution W has been changed. 24341ba4c6cSHeeho Park 244f6dfbefdSBarry Smith Note: 24541ba4c6cSHeeho Park If Y is changed then W is recomputed as X - Y 24641ba4c6cSHeeho Park 24741ba4c6cSHeeho Park Level: developer 24841ba4c6cSHeeho Park 249f6dfbefdSBarry Smith .seealso: `SNESNEWTONTRDC`, `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`, `SNESNewtonTRDCPreCheck() 25041ba4c6cSHeeho Park @*/ 2519371c9d4SSatish Balay static PetscErrorCode SNESNewtonTRDCPostCheck(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) { 25241ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 25341ba4c6cSHeeho Park 25441ba4c6cSHeeho Park PetscFunctionBegin; 25541ba4c6cSHeeho Park *changed_Y = PETSC_FALSE; 25641ba4c6cSHeeho Park *changed_W = PETSC_FALSE; 25741ba4c6cSHeeho Park if (tr->postcheck) { 2589566063dSJacob Faibussowitsch PetscCall((*tr->postcheck)(snes, X, Y, W, changed_Y, changed_W, tr->postcheckctx)); 25941ba4c6cSHeeho Park PetscValidLogicalCollectiveBool(snes, *changed_Y, 5); 26041ba4c6cSHeeho Park PetscValidLogicalCollectiveBool(snes, *changed_W, 6); 26141ba4c6cSHeeho Park } 26241ba4c6cSHeeho Park PetscFunctionReturn(0); 26341ba4c6cSHeeho Park } 26441ba4c6cSHeeho Park 26541ba4c6cSHeeho Park /* 26641ba4c6cSHeeho Park SNESSolve_NEWTONTRDC - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy 26741ba4c6cSHeeho Park (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of 26841ba4c6cSHeeho Park nonlinear equations 26941ba4c6cSHeeho Park 27041ba4c6cSHeeho Park */ 2719371c9d4SSatish Balay static PetscErrorCode SNESSolve_NEWTONTRDC(SNES snes) { 27241ba4c6cSHeeho Park SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC *)snes->data; 27341ba4c6cSHeeho Park Vec X, F, Y, G, W, GradF, YNtmp; 27441ba4c6cSHeeho Park Vec YCtmp; 27541ba4c6cSHeeho Park Mat jac; 27641ba4c6cSHeeho Park PetscInt maxits, i, j, lits, inner_count, bs; 27741ba4c6cSHeeho Park PetscReal rho, fnorm, gnorm, xnorm = 0, delta, ynorm, temp_xnorm, temp_ynorm; /* TRDC inner iteration */ 27841ba4c6cSHeeho Park PetscReal inorms[99]; /* need to make it dynamic eventually, fixed max block size of 99 for now */ 27941ba4c6cSHeeho Park PetscReal deltaM, ynnorm, f0, mp, gTy, g, yTHy; /* rho calculation */ 28041ba4c6cSHeeho Park PetscReal auk, gfnorm, ycnorm, c0, c1, c2, tau, tau_pos, tau_neg, gTBg; /* Cauchy Point */ 28141ba4c6cSHeeho Park KSP ksp; 28241ba4c6cSHeeho Park SNESConvergedReason reason = SNES_CONVERGED_ITERATING; 28341ba4c6cSHeeho Park PetscBool breakout = PETSC_FALSE; 28441ba4c6cSHeeho Park SNES_TRDC_KSPConverged_Ctx *ctx; 28541ba4c6cSHeeho Park PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), (*convdestroy)(void *); 28641ba4c6cSHeeho Park void *convctx; 28741ba4c6cSHeeho Park 28841ba4c6cSHeeho Park PetscFunctionBegin; 28941ba4c6cSHeeho Park maxits = snes->max_its; /* maximum number of iterations */ 29041ba4c6cSHeeho Park X = snes->vec_sol; /* solution vector */ 29141ba4c6cSHeeho Park F = snes->vec_func; /* residual vector */ 29241ba4c6cSHeeho Park Y = snes->work[0]; /* update vector */ 29341ba4c6cSHeeho Park G = snes->work[1]; /* updated residual */ 29441ba4c6cSHeeho Park W = snes->work[2]; /* temporary vector */ 29541ba4c6cSHeeho Park GradF = snes->work[3]; /* grad f = J^T F */ 29641ba4c6cSHeeho Park YNtmp = snes->work[4]; /* Newton solution */ 29741ba4c6cSHeeho Park YCtmp = snes->work[5]; /* Cauchy solution */ 29841ba4c6cSHeeho Park 2990b121fc5SBarry 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); 30041ba4c6cSHeeho Park 3019566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(YNtmp, &bs)); 30241ba4c6cSHeeho Park 3039566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 30441ba4c6cSHeeho Park snes->iter = 0; 3059566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 30641ba4c6cSHeeho Park 30741ba4c6cSHeeho Park /* Set the linear stopping criteria to use the More' trick. From tr.c */ 3089566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 3099566063dSJacob Faibussowitsch PetscCall(KSPGetConvergenceTest(ksp, &convtest, &convctx, &convdestroy)); 31041ba4c6cSHeeho Park if (convtest != SNESTRDC_KSPConverged_Private) { 3119566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 31241ba4c6cSHeeho Park ctx->snes = snes; 3139566063dSJacob Faibussowitsch PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 3149566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(ksp, SNESTRDC_KSPConverged_Private, ctx, SNESTRDC_KSPConverged_Destroy)); 3159566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTRDC_KSPConverged_Private\n")); 31641ba4c6cSHeeho Park } 31741ba4c6cSHeeho Park 31841ba4c6cSHeeho Park if (!snes->vec_func_init_set) { 3199566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */ 32041ba4c6cSHeeho Park } else snes->vec_func_init_set = PETSC_FALSE; 32141ba4c6cSHeeho Park 3229566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */ 32341ba4c6cSHeeho Park SNESCheckFunctionNorm(snes, fnorm); 3249566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */ 32541ba4c6cSHeeho Park 3269566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 32741ba4c6cSHeeho Park snes->norm = fnorm; 3289566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 32941ba4c6cSHeeho Park delta = xnorm ? neP->delta0 * xnorm : neP->delta0; /* initial trust region size scaled by xnorm */ 33041ba4c6cSHeeho Park deltaM = xnorm ? neP->deltaM * xnorm : neP->deltaM; /* maximum trust region size scaled by xnorm */ 33141ba4c6cSHeeho Park neP->delta = delta; 3329566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 3339566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, fnorm)); 33441ba4c6cSHeeho Park 33541ba4c6cSHeeho Park neP->rho_satisfied = PETSC_FALSE; 33641ba4c6cSHeeho Park 33741ba4c6cSHeeho Park /* test convergence */ 338dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 33941ba4c6cSHeeho Park if (snes->reason) PetscFunctionReturn(0); 34041ba4c6cSHeeho Park 34141ba4c6cSHeeho Park for (i = 0; i < maxits; i++) { 34241ba4c6cSHeeho Park PetscBool changed_y; 34341ba4c6cSHeeho Park PetscBool changed_w; 34441ba4c6cSHeeho Park 34541ba4c6cSHeeho Park /* dogleg method */ 3469566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 34741ba4c6cSHeeho Park SNESCheckJacobianDomainerror(snes); 3489566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian)); 3499566063dSJacob Faibussowitsch PetscCall(KSPSolve(snes->ksp, F, YNtmp)); /* Quasi Newton Solution */ 35041ba4c6cSHeeho Park SNESCheckKSPSolve(snes); /* this is necessary but old tr.c did not have it*/ 3519566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 3529566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(snes, &jac, NULL, NULL, NULL)); 35341ba4c6cSHeeho Park 35441ba4c6cSHeeho Park /* rescale Jacobian, Newton solution update, and re-calculate delta for multiphase (multivariable) 35541ba4c6cSHeeho Park for inner iteration and Cauchy direction calculation 35641ba4c6cSHeeho Park */ 35741ba4c6cSHeeho Park if (bs > 1 && neP->auto_scale_multiphase) { 3589566063dSJacob Faibussowitsch PetscCall(VecStrideNormAll(YNtmp, NORM_INFINITY, inorms)); 35941ba4c6cSHeeho Park for (j = 0; j < bs; j++) { 36041ba4c6cSHeeho Park if (neP->auto_scale_max > 1.0) { 361ad540459SPierre Jolivet if (inorms[j] < 1.0 / neP->auto_scale_max) inorms[j] = 1.0 / neP->auto_scale_max; 36241ba4c6cSHeeho Park } 3639566063dSJacob Faibussowitsch PetscCall(VecStrideSet(W, j, inorms[j])); 3649566063dSJacob Faibussowitsch PetscCall(VecStrideScale(YNtmp, j, 1.0 / inorms[j])); 3659566063dSJacob Faibussowitsch PetscCall(VecStrideScale(X, j, 1.0 / inorms[j])); 36641ba4c6cSHeeho Park } 3679566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, &xnorm)); 36841ba4c6cSHeeho Park if (i == 0) { 36941ba4c6cSHeeho Park delta = neP->delta0 * xnorm; 37041ba4c6cSHeeho Park } else { 37141ba4c6cSHeeho Park delta = neP->delta * xnorm; 37241ba4c6cSHeeho Park } 37341ba4c6cSHeeho Park deltaM = neP->deltaM * xnorm; 3749566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(jac, PETSC_NULL, W)); 37541ba4c6cSHeeho Park } 37641ba4c6cSHeeho Park 37741ba4c6cSHeeho Park /* calculating GradF of minimization function */ 3789566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(jac, F, GradF)); /* grad f = J^T F */ 3799566063dSJacob Faibussowitsch PetscCall(VecNorm(YNtmp, NORM_2, &ynnorm)); /* ynnorm <- || Y_newton || */ 38041ba4c6cSHeeho Park 38141ba4c6cSHeeho Park inner_count = 0; 38241ba4c6cSHeeho Park neP->rho_satisfied = PETSC_FALSE; 38341ba4c6cSHeeho Park while (1) { 38441ba4c6cSHeeho Park if (ynnorm <= delta) { /* see if the Newton solution is within the trust region */ 3859566063dSJacob Faibussowitsch PetscCall(VecCopy(YNtmp, Y)); 38641ba4c6cSHeeho Park } else if (neP->use_cauchy) { /* use Cauchy direction if enabled */ 3879566063dSJacob Faibussowitsch PetscCall(MatMult(jac, GradF, W)); 3889566063dSJacob Faibussowitsch PetscCall(VecDotRealPart(W, W, &gTBg)); /* completes GradF^T J^T J GradF */ 3899566063dSJacob Faibussowitsch PetscCall(VecNorm(GradF, NORM_2, &gfnorm)); /* grad f norm <- || grad f || */ 39041ba4c6cSHeeho Park if (gTBg <= 0.0) { 39141ba4c6cSHeeho Park auk = PETSC_MAX_REAL; 39241ba4c6cSHeeho Park } else { 39341ba4c6cSHeeho Park auk = PetscSqr(gfnorm) / gTBg; 39441ba4c6cSHeeho Park } 39541ba4c6cSHeeho Park auk = PetscMin(delta / gfnorm, auk); 3969566063dSJacob Faibussowitsch PetscCall(VecCopy(GradF, YCtmp)); /* this could be improved */ 3979566063dSJacob Faibussowitsch PetscCall(VecScale(YCtmp, auk)); /* YCtmp, Cauchy solution*/ 3989566063dSJacob Faibussowitsch PetscCall(VecNorm(YCtmp, NORM_2, &ycnorm)); /* ycnorm <- || Y_cauchy || */ 39941ba4c6cSHeeho Park if (ycnorm >= delta) { /* see if the Cauchy solution meets the criteria */ 4009566063dSJacob Faibussowitsch PetscCall(VecCopy(YCtmp, Y)); 4019566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "DL evaluated. delta: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n", (double)delta, (double)ynnorm, (double)ycnorm)); 40241ba4c6cSHeeho Park } else { /* take ratio, tau, of Cauchy and Newton direction and step */ 4039566063dSJacob Faibussowitsch PetscCall(VecAXPY(YNtmp, -1.0, YCtmp)); /* YCtmp = A, YNtmp = B */ 4049566063dSJacob Faibussowitsch PetscCall(VecNorm(YNtmp, NORM_2, &c0)); /* this could be improved */ 40541ba4c6cSHeeho Park c0 = PetscSqr(c0); 4069566063dSJacob Faibussowitsch PetscCall(VecDotRealPart(YCtmp, YNtmp, &c1)); 40741ba4c6cSHeeho Park c1 = 2.0 * c1; 4089566063dSJacob Faibussowitsch PetscCall(VecNorm(YCtmp, NORM_2, &c2)); /* this could be improved */ 40941ba4c6cSHeeho Park c2 = PetscSqr(c2) - PetscSqr(delta); 41041ba4c6cSHeeho Park tau_pos = (c1 + PetscSqrtReal(PetscSqr(c1) - 4. * c0 * c2)) / (2. * c0); /* quadratic formula */ 41141ba4c6cSHeeho Park tau_neg = (c1 - PetscSqrtReal(PetscSqr(c1) - 4. * c0 * c2)) / (2. * c0); 41241ba4c6cSHeeho Park tau = PetscMax(tau_pos, tau_neg); /* can tau_neg > tau_pos? I don't think so, but just in case. */ 4139566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "DL evaluated. tau: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n", (double)tau, (double)ynnorm, (double)ycnorm)); 4149566063dSJacob Faibussowitsch PetscCall(VecWAXPY(W, tau, YNtmp, YCtmp)); 4159566063dSJacob Faibussowitsch PetscCall(VecAXPY(W, -tau, YCtmp)); 4169566063dSJacob Faibussowitsch PetscCall(VecCopy(W, Y)); /* this could be improved */ 41741ba4c6cSHeeho Park } 41841ba4c6cSHeeho Park } else { 41941ba4c6cSHeeho Park /* if Cauchy is disabled, only use Newton direction */ 42041ba4c6cSHeeho Park auk = delta / ynnorm; 4219566063dSJacob Faibussowitsch PetscCall(VecScale(YNtmp, auk)); 4229566063dSJacob Faibussowitsch PetscCall(VecCopy(YNtmp, Y)); /* this could be improved (many VecCopy, VecNorm)*/ 42341ba4c6cSHeeho Park } 42441ba4c6cSHeeho Park 4259566063dSJacob Faibussowitsch PetscCall(VecNorm(Y, NORM_2, &ynorm)); /* compute the final ynorm */ 42641ba4c6cSHeeho Park f0 = 0.5 * PetscSqr(fnorm); /* minimizing function f(X) */ 4279566063dSJacob Faibussowitsch PetscCall(MatMult(jac, Y, W)); 4289566063dSJacob Faibussowitsch PetscCall(VecDotRealPart(W, W, &yTHy)); /* completes GradY^T J^T J GradY */ 4299566063dSJacob Faibussowitsch PetscCall(VecDotRealPart(GradF, Y, &gTy)); 43041ba4c6cSHeeho Park mp = f0 - gTy + 0.5 * yTHy; /* quadratic model to satisfy, -gTy because our update is X-Y*/ 43141ba4c6cSHeeho Park 43241ba4c6cSHeeho Park /* scale back solution update */ 43341ba4c6cSHeeho Park if (bs > 1 && neP->auto_scale_multiphase) { 43441ba4c6cSHeeho Park for (j = 0; j < bs; j++) { 4359566063dSJacob Faibussowitsch PetscCall(VecStrideScale(Y, j, inorms[j])); 43641ba4c6cSHeeho Park if (inner_count == 0) { 43741ba4c6cSHeeho Park /* TRDC inner algorithm does not need scaled X after calculating delta in the outer iteration */ 43841ba4c6cSHeeho Park /* need to scale back X to match Y and provide proper update to the external code */ 4399566063dSJacob Faibussowitsch PetscCall(VecStrideScale(X, j, inorms[j])); 44041ba4c6cSHeeho Park } 44141ba4c6cSHeeho Park } 4429566063dSJacob Faibussowitsch if (inner_count == 0) PetscCall(VecNorm(X, NORM_2, &temp_xnorm)); /* only in the first iteration */ 4439566063dSJacob Faibussowitsch PetscCall(VecNorm(Y, NORM_2, &temp_ynorm)); 44441ba4c6cSHeeho Park } else { 44541ba4c6cSHeeho Park temp_xnorm = xnorm; 44641ba4c6cSHeeho Park temp_ynorm = ynorm; 44741ba4c6cSHeeho Park } 44841ba4c6cSHeeho Park inner_count++; 44941ba4c6cSHeeho Park 45041ba4c6cSHeeho Park /* Evaluate the solution to meet the improvement ratio criteria */ 4519566063dSJacob Faibussowitsch PetscCall(SNESNewtonTRDCPreCheck(snes, X, Y, &changed_y)); 4529566063dSJacob Faibussowitsch PetscCall(VecWAXPY(W, -1.0, Y, X)); 4539566063dSJacob Faibussowitsch PetscCall(SNESNewtonTRDCPostCheck(snes, X, Y, W, &changed_y, &changed_w)); 4549566063dSJacob Faibussowitsch if (changed_y) PetscCall(VecWAXPY(W, -1.0, Y, X)); 4559566063dSJacob Faibussowitsch PetscCall(VecCopy(Y, snes->vec_sol_update)); 4569566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, W, G)); /* F(X-Y) = G */ 4579566063dSJacob Faibussowitsch PetscCall(VecNorm(G, NORM_2, &gnorm)); /* gnorm <- || g || */ 45841ba4c6cSHeeho Park SNESCheckFunctionNorm(snes, gnorm); 45941ba4c6cSHeeho Park g = 0.5 * PetscSqr(gnorm); /* minimizing function g(W) */ 46041ba4c6cSHeeho Park if (f0 == mp) rho = 0.0; 46141ba4c6cSHeeho Park else rho = (f0 - g) / (f0 - mp); /* actual improvement over predicted improvement */ 46241ba4c6cSHeeho Park 46341ba4c6cSHeeho Park if (rho < neP->eta2) { 46441ba4c6cSHeeho Park delta *= neP->t1; /* shrink the region */ 46541ba4c6cSHeeho Park } else if (rho > neP->eta3) { 46641ba4c6cSHeeho Park delta = PetscMin(neP->t2 * delta, deltaM); /* expand the region, but not greater than deltaM */ 46741ba4c6cSHeeho Park } 46841ba4c6cSHeeho Park 46941ba4c6cSHeeho Park neP->delta = delta; 47041ba4c6cSHeeho Park if (rho >= neP->eta1) { 47141ba4c6cSHeeho Park /* unscale delta and xnorm before going to the next outer iteration */ 47241ba4c6cSHeeho Park if (bs > 1 && neP->auto_scale_multiphase) { 47341ba4c6cSHeeho Park neP->delta = delta / xnorm; 47441ba4c6cSHeeho Park xnorm = temp_xnorm; 47541ba4c6cSHeeho Park ynorm = temp_ynorm; 47641ba4c6cSHeeho Park } 47741ba4c6cSHeeho Park neP->rho_satisfied = PETSC_TRUE; 47841ba4c6cSHeeho Park break; /* the improvement ratio is satisfactory */ 47941ba4c6cSHeeho Park } 4809566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Trying again in smaller region\n")); 48141ba4c6cSHeeho Park 48241ba4c6cSHeeho Park /* check to see if progress is hopeless */ 48341ba4c6cSHeeho Park neP->itflag = PETSC_FALSE; 48441ba4c6cSHeeho Park /* both delta, ynorm, and xnorm are either scaled or unscaled */ 4859566063dSJacob Faibussowitsch PetscCall(SNESTRDC_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP)); 48641ba4c6cSHeeho Park if (!reason) { 48741ba4c6cSHeeho Park /* temp_xnorm, temp_ynorm is always unscaled */ 48841ba4c6cSHeeho Park /* also the inner iteration already calculated the Jacobian and solved the matrix */ 48941ba4c6cSHeeho Park /* therefore, it should be passing iteration number of iter+1 instead of iter+0 in the first iteration and after */ 4909566063dSJacob Faibussowitsch PetscCall((*snes->ops->converged)(snes, snes->iter + 1, temp_xnorm, temp_ynorm, fnorm, &reason, snes->cnvP)); 49141ba4c6cSHeeho Park } 49241ba4c6cSHeeho Park /* if multiphase state changes, break out inner iteration */ 49341ba4c6cSHeeho Park if (reason == SNES_BREAKOUT_INNER_ITER) { 49441ba4c6cSHeeho Park if (bs > 1 && neP->auto_scale_multiphase) { 49541ba4c6cSHeeho Park /* unscale delta and xnorm before going to the next outer iteration */ 49641ba4c6cSHeeho Park neP->delta = delta / xnorm; 49741ba4c6cSHeeho Park xnorm = temp_xnorm; 49841ba4c6cSHeeho Park ynorm = temp_ynorm; 49941ba4c6cSHeeho Park } 50041ba4c6cSHeeho Park reason = SNES_CONVERGED_ITERATING; 50141ba4c6cSHeeho Park break; 50241ba4c6cSHeeho Park } 50341ba4c6cSHeeho Park if (reason == SNES_CONVERGED_SNORM_RELATIVE) reason = SNES_DIVERGED_INNER; 50441ba4c6cSHeeho Park if (reason) { 50541ba4c6cSHeeho Park if (reason < 0) { 50641ba4c6cSHeeho Park /* We're not progressing, so return with the current iterate */ 5079566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, i + 1, fnorm)); 50841ba4c6cSHeeho Park breakout = PETSC_TRUE; 50941ba4c6cSHeeho Park break; 51041ba4c6cSHeeho Park } else if (reason > 0) { 51141ba4c6cSHeeho Park /* We're converged, so return with the current iterate and update solution */ 5129566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, i + 1, fnorm)); 51341ba4c6cSHeeho Park breakout = PETSC_FALSE; 51441ba4c6cSHeeho Park break; 51541ba4c6cSHeeho Park } 51641ba4c6cSHeeho Park } 51741ba4c6cSHeeho Park snes->numFailures++; 51841ba4c6cSHeeho Park } 51941ba4c6cSHeeho Park if (!breakout) { 52041ba4c6cSHeeho Park /* Update function and solution vectors */ 52141ba4c6cSHeeho Park fnorm = gnorm; 5229566063dSJacob Faibussowitsch PetscCall(VecCopy(G, F)); 5239566063dSJacob Faibussowitsch PetscCall(VecCopy(W, X)); 52441ba4c6cSHeeho Park /* Monitor convergence */ 5259566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 52641ba4c6cSHeeho Park snes->iter = i + 1; 52741ba4c6cSHeeho Park snes->norm = fnorm; 52841ba4c6cSHeeho Park snes->xnorm = xnorm; 52941ba4c6cSHeeho Park snes->ynorm = ynorm; 5309566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 5319566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 5329566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 53341ba4c6cSHeeho Park /* Test for convergence, xnorm = || X || */ 53441ba4c6cSHeeho Park neP->itflag = PETSC_TRUE; 5359566063dSJacob Faibussowitsch if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm)); 536dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP); 53741ba4c6cSHeeho Park if (reason) break; 53841ba4c6cSHeeho Park } else break; 53941ba4c6cSHeeho Park } 54041ba4c6cSHeeho Park 5419566063dSJacob Faibussowitsch /* PetscCall(PetscFree(inorms)); */ 54241ba4c6cSHeeho Park if (i == maxits) { 5439566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits)); 54441ba4c6cSHeeho Park if (!reason) reason = SNES_DIVERGED_MAX_IT; 54541ba4c6cSHeeho Park } 5469566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 54741ba4c6cSHeeho Park snes->reason = reason; 5489566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 54941ba4c6cSHeeho Park if (convtest != SNESTRDC_KSPConverged_Private) { 5509566063dSJacob Faibussowitsch PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 5519566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 5529566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(ksp, convtest, convctx, convdestroy)); 55341ba4c6cSHeeho Park } 55441ba4c6cSHeeho Park PetscFunctionReturn(0); 55541ba4c6cSHeeho Park } 55641ba4c6cSHeeho Park 5579371c9d4SSatish Balay static PetscErrorCode SNESSetUp_NEWTONTRDC(SNES snes) { 55841ba4c6cSHeeho Park PetscFunctionBegin; 5599566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes, 6)); 5609566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 56141ba4c6cSHeeho Park PetscFunctionReturn(0); 56241ba4c6cSHeeho Park } 56341ba4c6cSHeeho Park 5649371c9d4SSatish Balay PetscErrorCode SNESReset_NEWTONTRDC(SNES snes) { 56541ba4c6cSHeeho Park PetscFunctionBegin; 56641ba4c6cSHeeho Park PetscFunctionReturn(0); 56741ba4c6cSHeeho Park } 56841ba4c6cSHeeho Park 5699371c9d4SSatish Balay static PetscErrorCode SNESDestroy_NEWTONTRDC(SNES snes) { 57041ba4c6cSHeeho Park PetscFunctionBegin; 5719566063dSJacob Faibussowitsch PetscCall(SNESReset_NEWTONTRDC(snes)); 5729566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 57341ba4c6cSHeeho Park PetscFunctionReturn(0); 57441ba4c6cSHeeho Park } 57541ba4c6cSHeeho Park 5769371c9d4SSatish Balay static PetscErrorCode SNESSetFromOptions_NEWTONTRDC(SNES snes, PetscOptionItems *PetscOptionsObject) { 57741ba4c6cSHeeho Park SNES_NEWTONTRDC *ctx = (SNES_NEWTONTRDC *)snes->data; 57841ba4c6cSHeeho Park 57941ba4c6cSHeeho Park PetscFunctionBegin; 580d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations"); 5819566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_tol", "Trust region tolerance", "SNESSetTrustRegionTolerance", snes->deltatol, &snes->deltatol, NULL)); 5829566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_eta1", "eta1", "None", ctx->eta1, &ctx->eta1, NULL)); 5839566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_eta2", "eta2", "None", ctx->eta2, &ctx->eta2, NULL)); 5849566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_eta3", "eta3", "None", ctx->eta3, &ctx->eta3, NULL)); 5859566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_t1", "t1", "None", ctx->t1, &ctx->t1, NULL)); 5869566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_t2", "t2", "None", ctx->t2, &ctx->t2, NULL)); 5879566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_deltaM", "deltaM", "None", ctx->deltaM, &ctx->deltaM, NULL)); 5889566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_delta0", "delta0", "None", ctx->delta0, &ctx->delta0, NULL)); 5899566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_trdc_auto_scale_max", "auto_scale_max", "None", ctx->auto_scale_max, &ctx->auto_scale_max, NULL)); 5909566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_trdc_use_cauchy", "use_cauchy", "use Cauchy step and direction", ctx->use_cauchy, &ctx->use_cauchy, NULL)); 5919566063dSJacob 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)); 592d0609cedSBarry Smith PetscOptionsHeadEnd(); 59341ba4c6cSHeeho Park PetscFunctionReturn(0); 59441ba4c6cSHeeho Park } 59541ba4c6cSHeeho Park 5969371c9d4SSatish Balay static PetscErrorCode SNESView_NEWTONTRDC(SNES snes, PetscViewer viewer) { 59741ba4c6cSHeeho Park SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 59841ba4c6cSHeeho Park PetscBool iascii; 59941ba4c6cSHeeho Park 60041ba4c6cSHeeho Park PetscFunctionBegin; 6019566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 60241ba4c6cSHeeho Park if (iascii) { 60363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Trust region tolerance %g (-snes_trtol)\n", (double)snes->deltatol)); 6049566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " eta1=%g, eta2=%g, eta3=%g\n", (double)tr->eta1, (double)tr->eta2, (double)tr->eta3)); 6059566063dSJacob 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)); 60641ba4c6cSHeeho Park } 60741ba4c6cSHeeho Park PetscFunctionReturn(0); 60841ba4c6cSHeeho Park } 609f6dfbefdSBarry Smith 61041ba4c6cSHeeho Park /*MC 61141ba4c6cSHeeho Park SNESNEWTONTRDC - Newton based nonlinear solver that uses trust-region dogleg method with Cauchy direction 61241ba4c6cSHeeho Park 613f6dfbefdSBarry Smith Options Database Keys: 61441ba4c6cSHeeho Park + -snes_trdc_tol <tol> - trust region tolerance 61541ba4c6cSHeeho Park . -snes_trdc_eta1 <eta1> - trust region parameter 0.0 <= eta1 <= eta2, rho >= eta1 breaks out of the inner iteration (default: eta1=0.001) 61641ba4c6cSHeeho Park . -snes_trdc_eta2 <eta2> - trust region parameter 0.0 <= eta1 <= eta2, rho <= eta2 shrinks the trust region (default: eta2=0.25) 61741ba4c6cSHeeho Park . -snes_trdc_eta3 <eta3> - trust region parameter eta3 > eta2, rho >= eta3 expands the trust region (default: eta3=0.75) 61841ba4c6cSHeeho Park . -snes_trdc_t1 <t1> - trust region parameter, shrinking factor of trust region (default: 0.25) 61941ba4c6cSHeeho Park . -snes_trdc_t2 <t2> - trust region parameter, expanding factor of trust region (default: 2.0) 62041ba4c6cSHeeho Park . -snes_trdc_deltaM <deltaM> - trust region parameter, max size of trust region, deltaM*norm2(x) (default: 0.5) 62141ba4c6cSHeeho Park . -snes_trdc_delta0 <delta0> - trust region parameter, initial size of trust region, delta0*norm2(x) (default: 0.1) 62241ba4c6cSHeeho Park . -snes_trdc_auto_scale_max <auto_scale_max> - used with auto_scale_multiphase, caps the maximum auto-scaling factor 62341ba4c6cSHeeho Park . -snes_trdc_use_cauchy <use_cauchy> - True uses dogleg Cauchy (Steepest Descent direction) step & direction in the trust region algorithm 62441ba4c6cSHeeho Park - -snes_trdc_auto_scale_multiphase <auto_scale_multiphase> - True turns on auto-scaling for multivariable block matrix for Cauchy and trust region 62541ba4c6cSHeeho Park 626f6dfbefdSBarry Smith Reference: 627f6dfbefdSBarry Smith . - * "Linear and Nonlinear Solvers for Simulating Multiphase Flow 62841ba4c6cSHeeho Park within Large-Scale Engineered Subsurface Systems" by Heeho D. Park, Glenn E. Hammond, 62941ba4c6cSHeeho Park Albert J. Valocchi, Tara LaForce. 63041ba4c6cSHeeho Park 63141ba4c6cSHeeho Park Level: intermediate 63241ba4c6cSHeeho Park 633f6dfbefdSBarry Smith .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESSetTrustRegionTolerance()`, 634f6dfbefdSBarry Smith `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`, 635f6dfbefdSBarry Smith `SNESNewtonTRDCGetRhoFlag()`, `SNESNewtonTRDCSetPreCheck()` 63641ba4c6cSHeeho Park M*/ 6379371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTRDC(SNES snes) { 63841ba4c6cSHeeho Park SNES_NEWTONTRDC *neP; 63941ba4c6cSHeeho Park 64041ba4c6cSHeeho Park PetscFunctionBegin; 64141ba4c6cSHeeho Park snes->ops->setup = SNESSetUp_NEWTONTRDC; 64241ba4c6cSHeeho Park snes->ops->solve = SNESSolve_NEWTONTRDC; 64341ba4c6cSHeeho Park snes->ops->destroy = SNESDestroy_NEWTONTRDC; 64441ba4c6cSHeeho Park snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTRDC; 64541ba4c6cSHeeho Park snes->ops->view = SNESView_NEWTONTRDC; 64641ba4c6cSHeeho Park snes->ops->reset = SNESReset_NEWTONTRDC; 64741ba4c6cSHeeho Park 64841ba4c6cSHeeho Park snes->usesksp = PETSC_TRUE; 64941ba4c6cSHeeho Park snes->usesnpc = PETSC_FALSE; 65041ba4c6cSHeeho Park 65141ba4c6cSHeeho Park snes->alwayscomputesfinalresidual = PETSC_TRUE; 65241ba4c6cSHeeho Park 653*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&neP)); 65441ba4c6cSHeeho Park snes->data = (void *)neP; 65541ba4c6cSHeeho Park neP->delta = 0.0; 65641ba4c6cSHeeho Park neP->delta0 = 0.1; 65741ba4c6cSHeeho Park neP->eta1 = 0.001; 65841ba4c6cSHeeho Park neP->eta2 = 0.25; 65941ba4c6cSHeeho Park neP->eta3 = 0.75; 66041ba4c6cSHeeho Park neP->t1 = 0.25; 66141ba4c6cSHeeho Park neP->t2 = 2.0; 66241ba4c6cSHeeho Park neP->deltaM = 0.5; 66341ba4c6cSHeeho Park neP->sigma = 0.0001; 66441ba4c6cSHeeho Park neP->itflag = PETSC_FALSE; 66541ba4c6cSHeeho Park neP->rnorm0 = 0.0; 66641ba4c6cSHeeho Park neP->ttol = 0.0; 66741ba4c6cSHeeho Park neP->use_cauchy = PETSC_TRUE; 66841ba4c6cSHeeho Park neP->auto_scale_multiphase = PETSC_FALSE; 66941ba4c6cSHeeho Park neP->auto_scale_max = -1.0; 67041ba4c6cSHeeho Park neP->rho_satisfied = PETSC_FALSE; 67141ba4c6cSHeeho Park snes->deltatol = 1.e-12; 67241ba4c6cSHeeho Park 67341ba4c6cSHeeho Park /* for multiphase (multivariable) scaling */ 67441ba4c6cSHeeho Park /* may be used for dynamic allocation of inorms, but it fails snes_tutorials-ex3_13 67541ba4c6cSHeeho Park on test forced DIVERGED_JACOBIAN_DOMAIN test. I will use static array for now. 6769566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(snes->work[0],&neP->bs)); 6779566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(neP->bs,&neP->inorms)); 67841ba4c6cSHeeho Park */ 67941ba4c6cSHeeho Park 68041ba4c6cSHeeho Park PetscFunctionReturn(0); 68141ba4c6cSHeeho Park } 682