1c6db04a5SJed Brown #include <../src/snes/impls/tr/trimpl.h> /*I "petscsnes.h" I*/ 24800dd8cSBarry Smith 3971273eeSBarry Smith typedef struct { 4971273eeSBarry Smith SNES snes; 5df8705c3SBarry Smith PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 6df8705c3SBarry Smith PetscErrorCode (*convdestroy)(void *); 7df8705c3SBarry Smith void *convctx; 8971273eeSBarry Smith } SNES_TR_KSPConverged_Ctx; 9971273eeSBarry Smith 104a221d59SStefano Zampini const char *const SNESNewtonTRFallbackTypes[] = {"NEWTON", "CAUCHY", "DOGLEG", "SNESNewtonTRFallbackType", "SNES_TR_FALLBACK_", NULL}; 11*24fb275aSStefano Zampini const char *const SNESNewtonTRQNTypes[] = {"NONE", "SAME", "DIFFERENT", "SNESNewtonTRQNType", "SNES_TR_QN_", NULL}; 12*24fb275aSStefano Zampini 13*24fb275aSStefano Zampini static PetscErrorCode SNESComputeJacobian_MATLMVM(SNES snes, Vec X, Mat J, Mat B, void *dummy) 14*24fb275aSStefano Zampini { 15*24fb275aSStefano Zampini PetscFunctionBegin; 16*24fb275aSStefano Zampini // PetscCall(MatLMVMSymBroydenSetDelta(B, _some_delta)); 17*24fb275aSStefano Zampini PetscCall(MatLMVMUpdate(B, X, snes->vec_func)); 18*24fb275aSStefano Zampini PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 19*24fb275aSStefano Zampini PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 20*24fb275aSStefano Zampini if (J != B) { 21*24fb275aSStefano Zampini // PetscCall(MatLMVMSymBroydenSetDelta(J, _some_delta)); 22*24fb275aSStefano Zampini PetscCall(MatLMVMUpdate(J, X, snes->vec_func)); 23*24fb275aSStefano Zampini PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 24*24fb275aSStefano Zampini PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 25*24fb275aSStefano Zampini } 26*24fb275aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 27*24fb275aSStefano Zampini } 284a221d59SStefano Zampini 29d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTR_KSPConverged_Private(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *cctx) 30d71ae5a4SJacob Faibussowitsch { 31971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx *)cctx; 32971273eeSBarry Smith SNES snes = ctx->snes; 3304d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 34df60cc22SBarry Smith Vec x; 35064f8208SBarry Smith PetscReal nrm; 36df60cc22SBarry Smith 373a40ed3dSBarry Smith 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)); 40a935fc98SLois Curfman McInnes /* Determine norm of solution */ 419566063dSJacob Faibussowitsch PetscCall(KSPBuildSolution(ksp, NULL, &x)); 42*24fb275aSStefano Zampini PetscCall(VecNorm(x, neP->norm, &nrm)); 43064f8208SBarry Smith if (nrm >= neP->delta) { 449566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Ending linear iteration early, delta=%g, length=%g\n", (double)neP->delta, (double)nrm)); 45329f5518SBarry Smith *reason = KSP_CONVERGED_STEP_LENGTH; 46df60cc22SBarry Smith } 473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 48df60cc22SBarry Smith } 4982bf6240SBarry Smith 50d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTR_KSPConverged_Destroy(void *cctx) 51d71ae5a4SJacob Faibussowitsch { 52971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx *)cctx; 53971273eeSBarry Smith 54971273eeSBarry Smith PetscFunctionBegin; 559566063dSJacob Faibussowitsch PetscCall((*ctx->convdestroy)(ctx->convctx)); 569566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 58971273eeSBarry Smith } 59971273eeSBarry Smith 60d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTR_Converged_Private(SNES snes, PetscInt it, PetscReal xnorm, PetscReal pnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy) 61d71ae5a4SJacob Faibussowitsch { 6204d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 6385385478SLisandro Dalcin 6485385478SLisandro Dalcin PetscFunctionBegin; 6585385478SLisandro Dalcin *reason = SNES_CONVERGED_ITERATING; 664a221d59SStefano Zampini if (neP->delta < snes->deltatol) { 674a221d59SStefano Zampini PetscCall(PetscInfo(snes, "Diverged due to too small a trust region %g<%g\n", (double)neP->delta, (double)snes->deltatol)); 681c6b2ff8SBarry Smith *reason = SNES_DIVERGED_TR_DELTA; 69e71169deSBarry Smith } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 7063a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n", snes->max_funcs)); 7185385478SLisandro Dalcin *reason = SNES_DIVERGED_FUNCTION_COUNT; 7285385478SLisandro Dalcin } 733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7485385478SLisandro Dalcin } 7585385478SLisandro Dalcin 764a221d59SStefano Zampini /*@ 77*24fb275aSStefano Zampini SNESNewtonTRSetNormType - Specify the type of norm to use for the computation of the trust region. 78*24fb275aSStefano Zampini 79*24fb275aSStefano Zampini Input Parameters: 80*24fb275aSStefano Zampini + snes - the nonlinear solver object 81*24fb275aSStefano Zampini - norm - the norm type 82*24fb275aSStefano Zampini 83*24fb275aSStefano Zampini Level: intermediate 84*24fb275aSStefano Zampini 85*24fb275aSStefano Zampini .seealso: `SNESNEWTONTR`, `NormType` 86*24fb275aSStefano Zampini @*/ 87*24fb275aSStefano Zampini PetscErrorCode SNESNewtonTRSetNormType(SNES snes, NormType norm) 88*24fb275aSStefano Zampini { 89*24fb275aSStefano Zampini PetscBool flg; 90*24fb275aSStefano Zampini 91*24fb275aSStefano Zampini PetscFunctionBegin; 92*24fb275aSStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 93*24fb275aSStefano Zampini PetscValidLogicalCollectiveEnum(snes, norm, 2); 94*24fb275aSStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 95*24fb275aSStefano Zampini if (flg) { 96*24fb275aSStefano Zampini SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 97*24fb275aSStefano Zampini 98*24fb275aSStefano Zampini tr->norm = norm; 99*24fb275aSStefano Zampini } 100*24fb275aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 101*24fb275aSStefano Zampini } 102*24fb275aSStefano Zampini 103*24fb275aSStefano Zampini /*@ 104*24fb275aSStefano Zampini SNESNewtonTRSetQNType - Specify to use a quasi-Newton model. 105*24fb275aSStefano Zampini 106*24fb275aSStefano Zampini Input Parameters: 107*24fb275aSStefano Zampini + snes - the nonlinear solver object 108*24fb275aSStefano Zampini - use - the type of approximations to be used 109*24fb275aSStefano Zampini 110*24fb275aSStefano Zampini Level: intermediate 111*24fb275aSStefano Zampini 112*24fb275aSStefano Zampini Notes: 113*24fb275aSStefano Zampini Options for the approximations can be set with the snes_tr_qn_ and snes_tr_qn_pre_ prefixes. 114*24fb275aSStefano Zampini 115*24fb275aSStefano Zampini .seealso: `SNESNEWTONTR`, `SNESNewtonTRQNType`, `MATLMVM` 116*24fb275aSStefano Zampini @*/ 117*24fb275aSStefano Zampini PetscErrorCode SNESNewtonTRSetQNType(SNES snes, SNESNewtonTRQNType use) 118*24fb275aSStefano Zampini { 119*24fb275aSStefano Zampini PetscBool flg; 120*24fb275aSStefano Zampini 121*24fb275aSStefano Zampini PetscFunctionBegin; 122*24fb275aSStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 123*24fb275aSStefano Zampini PetscValidLogicalCollectiveEnum(snes, use, 2); 124*24fb275aSStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 125*24fb275aSStefano Zampini if (flg) { 126*24fb275aSStefano Zampini SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 127*24fb275aSStefano Zampini 128*24fb275aSStefano Zampini tr->qn = use; 129*24fb275aSStefano Zampini } 130*24fb275aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 131*24fb275aSStefano Zampini } 132*24fb275aSStefano Zampini 133*24fb275aSStefano Zampini /*@ 134420bcc1bSBarry Smith SNESNewtonTRSetFallbackType - Set the type of fallback to use if the solution of the trust region subproblem is outside the radius 1354a221d59SStefano Zampini 1364a221d59SStefano Zampini Input Parameters: 1374a221d59SStefano Zampini + snes - the nonlinear solver object 1384a221d59SStefano Zampini - ftype - the fallback type, see `SNESNewtonTRFallbackType` 1394a221d59SStefano Zampini 1404a221d59SStefano Zampini Level: intermediate 1414a221d59SStefano Zampini 142420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPreCheck()`, 1434a221d59SStefano Zampini `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()` 1444a221d59SStefano Zampini @*/ 1454a221d59SStefano Zampini PetscErrorCode SNESNewtonTRSetFallbackType(SNES snes, SNESNewtonTRFallbackType ftype) 1464a221d59SStefano Zampini { 1474a221d59SStefano Zampini SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 1484a221d59SStefano Zampini PetscBool flg; 1494a221d59SStefano Zampini 1504a221d59SStefano Zampini PetscFunctionBegin; 1514a221d59SStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1524a221d59SStefano Zampini PetscValidLogicalCollectiveEnum(snes, ftype, 2); 1534a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 1544a221d59SStefano Zampini if (flg) tr->fallback = ftype; 1554a221d59SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 1564a221d59SStefano Zampini } 1574a221d59SStefano Zampini 1587cb011f5SBarry Smith /*@C 159c9368356SGlenn Hammond SNESNewtonTRSetPreCheck - Sets a user function that is called before the search step has been determined. 1604a221d59SStefano Zampini Allows the user a chance to change or override the trust region decision. 161f6dfbefdSBarry Smith 162c3339decSBarry Smith Logically Collective 163c9368356SGlenn Hammond 164c9368356SGlenn Hammond Input Parameters: 165c9368356SGlenn Hammond + snes - the nonlinear solver object 16620f4b53cSBarry Smith . func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPreCheck()` 16720f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 168c9368356SGlenn Hammond 16920f4b53cSBarry Smith Level: deprecated (since 3.19) 170c9368356SGlenn Hammond 171f6dfbefdSBarry Smith Note: 1724a221d59SStefano Zampini This function is called BEFORE the function evaluation within the solver. 173c9368356SGlenn Hammond 174420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`, 175c9368356SGlenn Hammond @*/ 176d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, PetscBool *, void *), void *ctx) 177d71ae5a4SJacob Faibussowitsch { 178c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 1794a221d59SStefano Zampini PetscBool flg; 180c9368356SGlenn Hammond 181c9368356SGlenn Hammond PetscFunctionBegin; 182c9368356SGlenn Hammond PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1834a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 1844a221d59SStefano Zampini if (flg) { 185c9368356SGlenn Hammond if (func) tr->precheck = func; 186c9368356SGlenn Hammond if (ctx) tr->precheckctx = ctx; 1874a221d59SStefano Zampini } 1883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 189c9368356SGlenn Hammond } 190c9368356SGlenn Hammond 191c9368356SGlenn Hammond /*@C 192c9368356SGlenn Hammond SNESNewtonTRGetPreCheck - Gets the pre-check function 193c9368356SGlenn Hammond 19420f4b53cSBarry Smith Deprecated use `SNESNEWTONDCTRDC` 19520f4b53cSBarry Smith 19620f4b53cSBarry Smith Not Collective 197c9368356SGlenn Hammond 198c9368356SGlenn Hammond Input Parameter: 199c9368356SGlenn Hammond . snes - the nonlinear solver context 200c9368356SGlenn Hammond 201c9368356SGlenn Hammond Output Parameters: 20220f4b53cSBarry Smith + func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPreCheck()` 20320f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 204c9368356SGlenn Hammond 20520f4b53cSBarry Smith Level: deprecated (since 3.19) 206c9368356SGlenn Hammond 207420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRPreCheck()` 208c9368356SGlenn Hammond @*/ 209d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, PetscBool *, void *), void **ctx) 210d71ae5a4SJacob Faibussowitsch { 211c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 2124a221d59SStefano Zampini PetscBool flg; 213c9368356SGlenn Hammond 214c9368356SGlenn Hammond PetscFunctionBegin; 215c9368356SGlenn Hammond PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 2164a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 2174a221d59SStefano Zampini PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 218c9368356SGlenn Hammond if (func) *func = tr->precheck; 219c9368356SGlenn Hammond if (ctx) *ctx = tr->precheckctx; 2203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 221c9368356SGlenn Hammond } 222c9368356SGlenn Hammond 223c9368356SGlenn Hammond /*@C 2247cb011f5SBarry Smith SNESNewtonTRSetPostCheck - Sets a user function that is called after the search step has been determined but before the next 2254a221d59SStefano Zampini function evaluation. Allows the user a chance to change or override the internal decision of the solver 226f6dfbefdSBarry Smith 227c3339decSBarry Smith Logically Collective 2287cb011f5SBarry Smith 2297cb011f5SBarry Smith Input Parameters: 2307cb011f5SBarry Smith + snes - the nonlinear solver object 23120f4b53cSBarry Smith . func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPostCheck()` 23220f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 2337cb011f5SBarry Smith 23420f4b53cSBarry Smith Level: deprecated (since 3.19) 2357cb011f5SBarry Smith 236f6dfbefdSBarry Smith Note: 2374a221d59SStefano Zampini This function is called BEFORE the function evaluation within the solver while the function set in 238f6dfbefdSBarry Smith `SNESLineSearchSetPostCheck()` is called AFTER the function evaluation. 2397cb011f5SBarry Smith 240420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRPostCheck()`, `SNESNewtonTRGetPostCheck()`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()` 2417cb011f5SBarry Smith @*/ 242d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx) 243d71ae5a4SJacob Faibussowitsch { 2447cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 2454a221d59SStefano Zampini PetscBool flg; 2467cb011f5SBarry Smith 2477cb011f5SBarry Smith PetscFunctionBegin; 2487cb011f5SBarry Smith PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 2494a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 2504a221d59SStefano Zampini if (flg) { 2517cb011f5SBarry Smith if (func) tr->postcheck = func; 2527cb011f5SBarry Smith if (ctx) tr->postcheckctx = ctx; 2534a221d59SStefano Zampini } 2543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2557cb011f5SBarry Smith } 2567cb011f5SBarry Smith 2577cb011f5SBarry Smith /*@C 2587cb011f5SBarry Smith SNESNewtonTRGetPostCheck - Gets the post-check function 2597cb011f5SBarry Smith 26020f4b53cSBarry Smith Not Collective 2617cb011f5SBarry Smith 2627cb011f5SBarry Smith Input Parameter: 2637cb011f5SBarry Smith . snes - the nonlinear solver context 2647cb011f5SBarry Smith 2657cb011f5SBarry Smith Output Parameters: 26620f4b53cSBarry Smith + func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPostCheck()` 26720f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 2687cb011f5SBarry Smith 2697cb011f5SBarry Smith Level: intermediate 2707cb011f5SBarry Smith 271420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRPostCheck()` 2727cb011f5SBarry Smith @*/ 273d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) 274d71ae5a4SJacob Faibussowitsch { 2757cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 2764a221d59SStefano Zampini PetscBool flg; 2777cb011f5SBarry Smith 2787cb011f5SBarry Smith PetscFunctionBegin; 2797cb011f5SBarry Smith PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 2804a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 2814a221d59SStefano Zampini PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 2827cb011f5SBarry Smith if (func) *func = tr->postcheck; 2837cb011f5SBarry Smith if (ctx) *ctx = tr->postcheckctx; 2843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2857cb011f5SBarry Smith } 2867cb011f5SBarry Smith 2877cb011f5SBarry Smith /*@C 2884a221d59SStefano Zampini SNESNewtonTRPreCheck - Runs the precheck routine 289c9368356SGlenn Hammond 290c3339decSBarry Smith Logically Collective 291c9368356SGlenn Hammond 292c9368356SGlenn Hammond Input Parameters: 293c9368356SGlenn Hammond + snes - the solver 294c9368356SGlenn Hammond . X - The last solution 295c9368356SGlenn Hammond - Y - The step direction 296c9368356SGlenn Hammond 2972fe279fdSBarry Smith Output Parameter: 2982fe279fdSBarry Smith . changed_Y - Indicator that the step direction `Y` has been changed. 299c9368356SGlenn Hammond 3004a221d59SStefano Zampini Level: intermediate 301c9368356SGlenn Hammond 302420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRPostCheck()` 303c9368356SGlenn Hammond @*/ 3044a221d59SStefano Zampini PetscErrorCode SNESNewtonTRPreCheck(SNES snes, Vec X, Vec Y, PetscBool *changed_Y) 305d71ae5a4SJacob Faibussowitsch { 306c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 3074a221d59SStefano Zampini PetscBool flg; 308c9368356SGlenn Hammond 309c9368356SGlenn Hammond PetscFunctionBegin; 3104a221d59SStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 3114a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 3124a221d59SStefano Zampini PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 313c9368356SGlenn Hammond *changed_Y = PETSC_FALSE; 314c9368356SGlenn Hammond if (tr->precheck) { 3159566063dSJacob Faibussowitsch PetscCall((*tr->precheck)(snes, X, Y, changed_Y, tr->precheckctx)); 316c9368356SGlenn Hammond PetscValidLogicalCollectiveBool(snes, *changed_Y, 4); 317c9368356SGlenn Hammond } 3183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 319c9368356SGlenn Hammond } 320c9368356SGlenn Hammond 321c9368356SGlenn Hammond /*@C 3224a221d59SStefano Zampini SNESNewtonTRPostCheck - Runs the postcheck routine 3237cb011f5SBarry Smith 324c3339decSBarry Smith Logically Collective 3257cb011f5SBarry Smith 3267cb011f5SBarry Smith Input Parameters: 3276b867d5aSJose E. Roman + snes - the solver 3286b867d5aSJose E. Roman . X - The last solution 3297cb011f5SBarry Smith . Y - The full step direction 3303312a946SBarry Smith - W - The updated solution, W = X - Y 3317cb011f5SBarry Smith 3327cb011f5SBarry Smith Output Parameters: 3333312a946SBarry Smith + changed_Y - indicator if step has been changed 3343312a946SBarry Smith - changed_W - Indicator if the new candidate solution W has been changed. 3357cb011f5SBarry Smith 336f6dfbefdSBarry Smith Note: 3373312a946SBarry Smith If Y is changed then W is recomputed as X - Y 3387cb011f5SBarry Smith 3394a221d59SStefano Zampini Level: intermediate 3407cb011f5SBarry Smith 341420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`, `SNESNewtonTRPreCheck()` 3427cb011f5SBarry Smith @*/ 3434a221d59SStefano Zampini PetscErrorCode SNESNewtonTRPostCheck(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) 344d71ae5a4SJacob Faibussowitsch { 3457cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 3464a221d59SStefano Zampini PetscBool flg; 3477cb011f5SBarry Smith 3487cb011f5SBarry Smith PetscFunctionBegin; 3494a221d59SStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 3504a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 3514a221d59SStefano Zampini PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 352c9368356SGlenn Hammond *changed_Y = PETSC_FALSE; 3537cb011f5SBarry Smith *changed_W = PETSC_FALSE; 3547cb011f5SBarry Smith if (tr->postcheck) { 3559566063dSJacob Faibussowitsch PetscCall((*tr->postcheck)(snes, X, Y, W, changed_Y, changed_W, tr->postcheckctx)); 356c9368356SGlenn Hammond PetscValidLogicalCollectiveBool(snes, *changed_Y, 5); 3577cb011f5SBarry Smith PetscValidLogicalCollectiveBool(snes, *changed_W, 6); 3587cb011f5SBarry Smith } 3593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3607cb011f5SBarry Smith } 36185385478SLisandro Dalcin 362*24fb275aSStefano Zampini /* stable implementation of roots of a*x^2 + b*x + c = 0 */ 3634a221d59SStefano Zampini static inline void PetscQuadraticRoots(PetscReal a, PetscReal b, PetscReal c, PetscReal *xm, PetscReal *xp) 3644a221d59SStefano Zampini { 3654a221d59SStefano Zampini PetscReal temp = -0.5 * (b + PetscCopysignReal(1.0, b) * PetscSqrtReal(b * b - 4 * a * c)); 3664a221d59SStefano Zampini PetscReal x1 = temp / a; 3674a221d59SStefano Zampini PetscReal x2 = c / temp; 3684a221d59SStefano Zampini *xm = PetscMin(x1, x2); 3694a221d59SStefano Zampini *xp = PetscMax(x1, x2); 3704a221d59SStefano Zampini } 3714a221d59SStefano Zampini 3727aa289d8SStefano Zampini /* Computes the quadratic model difference */ 373*24fb275aSStefano Zampini static PetscErrorCode SNESNewtonTRQuadraticDelta(SNES snes, Mat J, PetscBool has_objective, Vec Y, Vec GradF, Vec W, PetscReal *yTHy_, PetscReal *gTy_, PetscReal *deltaqm) 3747aa289d8SStefano Zampini { 375*24fb275aSStefano Zampini PetscReal yTHy, gTy; 376*24fb275aSStefano Zampini 3777aa289d8SStefano Zampini PetscFunctionBegin; 378*24fb275aSStefano Zampini PetscCall(MatMult(J, Y, W)); 379*24fb275aSStefano Zampini if (has_objective) PetscCall(VecDotRealPart(Y, W, &yTHy)); 380*24fb275aSStefano Zampini else PetscCall(VecDotRealPart(W, W, &yTHy)); /* Gauss-Newton approximation J^t * J */ 381*24fb275aSStefano Zampini PetscCall(VecDotRealPart(GradF, Y, &gTy)); 382*24fb275aSStefano Zampini *deltaqm = -(-(gTy) + 0.5 * (yTHy)); /* difference in quadratic model, -gTy because SNES solves it this way */ 383*24fb275aSStefano Zampini if (yTHy_) *yTHy_ = yTHy; 384*24fb275aSStefano Zampini if (gTy_) *gTy_ = gTy; 385*24fb275aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 386*24fb275aSStefano Zampini } 387*24fb275aSStefano Zampini 388*24fb275aSStefano Zampini /* Computes the new objective given X = Xk, Y = direction 389*24fb275aSStefano Zampini W work vector, on output W = X - Y 390*24fb275aSStefano Zampini G work vector, on output G = SNESFunction(W) */ 391*24fb275aSStefano Zampini static PetscErrorCode SNESNewtonTRObjective(SNES snes, PetscBool has_objective, Vec X, Vec Y, Vec W, Vec G, PetscReal *gnorm, PetscReal *fkp1) 392*24fb275aSStefano Zampini { 393*24fb275aSStefano Zampini PetscBool changed_y, changed_w; 394*24fb275aSStefano Zampini 395*24fb275aSStefano Zampini PetscFunctionBegin; 396*24fb275aSStefano Zampini /* TODO: we can add a linesearch here */ 397*24fb275aSStefano Zampini PetscCall(SNESNewtonTRPreCheck(snes, X, Y, &changed_y)); 398*24fb275aSStefano Zampini PetscCall(VecWAXPY(W, -1.0, Y, X)); /* Xkp1 */ 399*24fb275aSStefano Zampini PetscCall(SNESNewtonTRPostCheck(snes, X, Y, W, &changed_y, &changed_w)); 400*24fb275aSStefano Zampini if (changed_y && !changed_w) PetscCall(VecWAXPY(W, -1.0, Y, X)); 401*24fb275aSStefano Zampini 402*24fb275aSStefano Zampini PetscCall(SNESComputeFunction(snes, W, G)); /* F(Xkp1) = G */ 403*24fb275aSStefano Zampini PetscCall(VecNorm(G, NORM_2, gnorm)); 404*24fb275aSStefano Zampini if (has_objective) PetscCall(SNESComputeObjective(snes, W, fkp1)); 405*24fb275aSStefano Zampini else *fkp1 = 0.5 * PetscSqr(*gnorm); 406*24fb275aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 407*24fb275aSStefano Zampini } 408*24fb275aSStefano Zampini 409*24fb275aSStefano Zampini static PetscErrorCode SNESSetUpQN_NEWTONTR(SNES snes) 410*24fb275aSStefano Zampini { 411*24fb275aSStefano Zampini SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 412*24fb275aSStefano Zampini 413*24fb275aSStefano Zampini PetscFunctionBegin; 414*24fb275aSStefano Zampini PetscCall(MatDestroy(&tr->qnB)); 415*24fb275aSStefano Zampini PetscCall(MatDestroy(&tr->qnB_pre)); 416*24fb275aSStefano Zampini if (tr->qn) { 417*24fb275aSStefano Zampini PetscInt n, N; 418*24fb275aSStefano Zampini const char *optionsprefix; 419*24fb275aSStefano Zampini Mat B; 420*24fb275aSStefano Zampini 421*24fb275aSStefano Zampini PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &B)); 422*24fb275aSStefano Zampini PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix)); 423*24fb275aSStefano Zampini PetscCall(MatSetOptionsPrefix(B, "snes_tr_qn_")); 424*24fb275aSStefano Zampini PetscCall(MatAppendOptionsPrefix(B, optionsprefix)); 425*24fb275aSStefano Zampini PetscCall(MatSetType(B, MATLMVMBFGS)); 426*24fb275aSStefano Zampini PetscCall(VecGetLocalSize(snes->vec_sol, &n)); 427*24fb275aSStefano Zampini PetscCall(VecGetSize(snes->vec_sol, &N)); 428*24fb275aSStefano Zampini PetscCall(MatSetSizes(B, n, n, N, N)); 429*24fb275aSStefano Zampini PetscCall(MatSetUp(B)); 430*24fb275aSStefano Zampini PetscCall(MatSetFromOptions(B)); 431*24fb275aSStefano Zampini PetscCall(MatLMVMAllocate(B, snes->vec_sol, snes->vec_func)); 432*24fb275aSStefano Zampini tr->qnB = B; 433*24fb275aSStefano Zampini if (tr->qn == SNES_TR_QN_DIFFERENT) { 434*24fb275aSStefano Zampini PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &B)); 435*24fb275aSStefano Zampini PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix)); 436*24fb275aSStefano Zampini PetscCall(MatSetOptionsPrefix(B, "snes_tr_qn_pre_")); 437*24fb275aSStefano Zampini PetscCall(MatAppendOptionsPrefix(B, optionsprefix)); 438*24fb275aSStefano Zampini PetscCall(MatSetType(B, MATLMVMBFGS)); 439*24fb275aSStefano Zampini PetscCall(MatSetSizes(B, n, n, N, N)); 440*24fb275aSStefano Zampini PetscCall(MatSetUp(B)); 441*24fb275aSStefano Zampini PetscCall(MatSetFromOptions(B)); 442*24fb275aSStefano Zampini PetscCall(MatLMVMAllocate(B, snes->vec_sol, snes->vec_func)); 443*24fb275aSStefano Zampini tr->qnB_pre = B; 444*24fb275aSStefano Zampini } else { 445*24fb275aSStefano Zampini PetscCall(PetscObjectReference((PetscObject)tr->qnB)); 446*24fb275aSStefano Zampini tr->qnB_pre = tr->qnB; 447*24fb275aSStefano Zampini } 448*24fb275aSStefano Zampini } 4497aa289d8SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 4507aa289d8SStefano Zampini } 4517aa289d8SStefano Zampini 452df60cc22SBarry Smith /* 4534a221d59SStefano Zampini SNESSolve_NEWTONTR - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy 4544a221d59SStefano Zampini (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of 4554a221d59SStefano Zampini nonlinear equations 4564800dd8cSBarry Smith 4574800dd8cSBarry Smith */ 458d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NEWTONTR(SNES snes) 459d71ae5a4SJacob Faibussowitsch { 46004d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 461*24fb275aSStefano Zampini Vec X, F, Y, G, W, GradF, YU, Yc; 4624a221d59SStefano Zampini PetscInt maxits, lits; 4634b0a5c37SStefano Zampini PetscReal rho, fnorm, gnorm = 0.0, xnorm = 0.0, delta, ynorm; 464*24fb275aSStefano Zampini PetscReal deltaM, fk, fkp1, deltaqm = 0.0, gTy = 0.0, yTHy = 0.0; 465*24fb275aSStefano Zampini PetscReal auk, tauk, gfnorm, gfnorm_k, ycnorm, gTBg, objmin = 0.0, beta_k = 1.0; 466a0254a93SStefano Zampini PC pc; 467*24fb275aSStefano Zampini Mat J, Jp; 468*24fb275aSStefano Zampini PetscBool already_done = PETSC_FALSE, on_boundary; 4697aa289d8SStefano Zampini PetscBool clear_converged_test, rho_satisfied, has_objective; 470df8705c3SBarry Smith SNES_TR_KSPConverged_Ctx *ctx; 4715e28bcb6Sprj- void *convctx; 472*24fb275aSStefano Zampini 4734a221d59SStefano Zampini PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), (*convdestroy)(void *); 4744a221d59SStefano Zampini PetscErrorCode (*objective)(SNES, Vec, PetscReal *, void *); 4754800dd8cSBarry Smith 4763a40ed3dSBarry Smith PetscFunctionBegin; 4774a221d59SStefano Zampini PetscCall(SNESGetObjective(snes, &objective, NULL)); 4787aa289d8SStefano Zampini has_objective = objective ? PETSC_TRUE : PETSC_FALSE; 479c579b300SPatrick Farrell 480fbe28522SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 481fbe28522SBarry Smith X = snes->vec_sol; /* solution vector */ 48239e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 4834a221d59SStefano Zampini Y = snes->vec_sol_update; /* update vector */ 4844a221d59SStefano Zampini G = snes->work[0]; /* updated residual */ 4854a221d59SStefano Zampini W = snes->work[1]; /* temporary vector */ 4867aa289d8SStefano Zampini GradF = !has_objective ? snes->work[2] : snes->vec_func; /* grad f = J^T F */ 4874a221d59SStefano Zampini YU = snes->work[3]; /* work vector for dogleg method */ 488*24fb275aSStefano Zampini Yc = snes->work[4]; /* Cauchy point */ 4894a221d59SStefano Zampini 4904a221d59SStefano Zampini 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); 4914800dd8cSBarry Smith 4929566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 4934c49b128SBarry Smith snes->iter = 0; 4949566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 4954800dd8cSBarry Smith 496*24fb275aSStefano Zampini /* setup QN matrices if needed */ 497*24fb275aSStefano Zampini PetscCall(SNESSetUpQN_NEWTONTR(snes)); 498*24fb275aSStefano Zampini 4994a221d59SStefano Zampini /* Set the linear stopping criteria to use the More' trick if needed */ 5004a221d59SStefano Zampini clear_converged_test = PETSC_FALSE; 501a0254a93SStefano Zampini PetscCall(SNESGetKSP(snes, &snes->ksp)); 502a0254a93SStefano Zampini PetscCall(KSPGetConvergenceTest(snes->ksp, &convtest, &convctx, &convdestroy)); 503fcc61681SStefano Zampini if (convtest != SNESTR_KSPConverged_Private) { 5044a221d59SStefano Zampini clear_converged_test = PETSC_TRUE; 5059566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 506df8705c3SBarry Smith ctx->snes = snes; 507a0254a93SStefano Zampini PetscCall(KSPGetAndClearConvergenceTest(snes->ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 508a0254a93SStefano Zampini PetscCall(KSPSetConvergenceTest(snes->ksp, SNESTR_KSPConverged_Private, ctx, SNESTR_KSPConverged_Destroy)); 5099566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTR_KSPConverged_Private\n")); 510df8705c3SBarry Smith } 511df8705c3SBarry Smith 512e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 5139566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */ 5141aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 515e4ed7901SPeter Brune 5169566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */ 517422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 5189566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */ 5194a221d59SStefano Zampini 5209566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 521fbe28522SBarry Smith snes->norm = fnorm; 5229566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 5234a221d59SStefano Zampini delta = neP->delta0; 5244a221d59SStefano Zampini deltaM = neP->deltaM; 5254800dd8cSBarry Smith neP->delta = delta; 5269566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 527b37302c6SLois Curfman McInnes 52885385478SLisandro Dalcin /* test convergence */ 5294a221d59SStefano Zampini rho_satisfied = PETSC_FALSE; 5302d157150SStefano Zampini PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm)); 5312d157150SStefano Zampini PetscCall(SNESMonitor(snes, 0, fnorm)); 5323ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 5333f149594SLisandro Dalcin 5347aa289d8SStefano Zampini if (has_objective) PetscCall(SNESComputeObjective(snes, X, &fk)); 5354a221d59SStefano Zampini else fk = 0.5 * PetscSqr(fnorm); /* obj(x) = 0.5 * ||F(x)||^2 */ 53699a96b7cSMatthew Knepley 537a0254a93SStefano Zampini /* hook state vector to BFGS preconditioner */ 538a0254a93SStefano Zampini PetscCall(KSPGetPC(snes->ksp, &pc)); 539a0254a93SStefano Zampini PetscCall(PCLMVMSetUpdateVec(pc, X)); 540a0254a93SStefano Zampini 541*24fb275aSStefano Zampini if (neP->kmdc) PetscCall(KSPSetComputeEigenvalues(snes->ksp, PETSC_TRUE)); 5426b5873e3SBarry Smith 543*24fb275aSStefano Zampini while (snes->iter < maxits) { 54412d0050eSStefano Zampini /* calculating Jacobian and GradF of minimization function only once */ 5454a221d59SStefano Zampini if (!already_done) { 54612d0050eSStefano Zampini /* Call general purpose update function */ 54712d0050eSStefano Zampini PetscTryTypeMethod(snes, update, snes->iter); 54812d0050eSStefano Zampini 5494b0a5c37SStefano Zampini /* apply the nonlinear preconditioner */ 5504b0a5c37SStefano Zampini if (snes->npc && snes->npcside == PC_RIGHT) { 5514b0a5c37SStefano Zampini SNESConvergedReason reason; 5524b0a5c37SStefano Zampini 5534b0a5c37SStefano Zampini PetscCall(SNESSetInitialFunction(snes->npc, F)); 5544b0a5c37SStefano Zampini PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0)); 5554b0a5c37SStefano Zampini PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X)); 5564b0a5c37SStefano Zampini PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0)); 5574b0a5c37SStefano Zampini PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 5584b0a5c37SStefano Zampini if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) { 5594b0a5c37SStefano Zampini snes->reason = SNES_DIVERGED_INNER; 5604b0a5c37SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 5614b0a5c37SStefano Zampini } 5624b0a5c37SStefano Zampini // XXX 5634b0a5c37SStefano Zampini PetscCall(SNESGetNPCFunction(snes, F, &fnorm)); 5644b0a5c37SStefano Zampini } else if (snes->ops->update) { /* if update is present, recompute objective function and function norm */ 56512d0050eSStefano Zampini PetscCall(SNESComputeFunction(snes, X, F)); 56612d0050eSStefano Zampini } 56712d0050eSStefano Zampini 56812d0050eSStefano Zampini /* Jacobian */ 569*24fb275aSStefano Zampini J = NULL; 570*24fb275aSStefano Zampini Jp = NULL; 571*24fb275aSStefano Zampini if (!neP->qnB) { 5724a221d59SStefano Zampini PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 573*24fb275aSStefano Zampini J = snes->jacobian; 574*24fb275aSStefano Zampini Jp = snes->jacobian_pre; 575*24fb275aSStefano Zampini } else { /* QN model */ 576*24fb275aSStefano Zampini PetscCall(SNESComputeJacobian_MATLMVM(snes, X, neP->qnB, neP->qnB_pre, NULL)); 577*24fb275aSStefano Zampini J = neP->qnB; 578*24fb275aSStefano Zampini Jp = neP->qnB_pre; 579*24fb275aSStefano Zampini } 5804a221d59SStefano Zampini SNESCheckJacobianDomainerror(snes); 58112d0050eSStefano Zampini 582*24fb275aSStefano Zampini /* objective function */ 583*24fb275aSStefano Zampini PetscCall(VecNorm(F, NORM_2, &fnorm)); 584*24fb275aSStefano Zampini if (has_objective) PetscCall(SNESComputeObjective(snes, X, &fk)); 585*24fb275aSStefano Zampini else fk = 0.5 * PetscSqr(fnorm); /* obj(x) = 0.5 * ||F(x)||^2 */ 586*24fb275aSStefano Zampini 58712d0050eSStefano Zampini /* GradF */ 5887aa289d8SStefano Zampini if (has_objective) gfnorm = fnorm; 5897aa289d8SStefano Zampini else { 590*24fb275aSStefano Zampini PetscCall(MatMultTranspose(J, F, GradF)); /* grad f = J^T F */ 5917aa289d8SStefano Zampini PetscCall(VecNorm(GradF, NORM_2, &gfnorm)); 592fbe28522SBarry Smith } 593*24fb275aSStefano Zampini PetscCall(VecNorm(GradF, neP->norm, &gfnorm_k)); 5947aa289d8SStefano Zampini } 5957aa289d8SStefano Zampini already_done = PETSC_TRUE; 5967aa289d8SStefano Zampini 5974b0a5c37SStefano Zampini /* solve trust-region subproblem */ 5984b0a5c37SStefano Zampini 599*24fb275aSStefano Zampini /* first compute Cauchy Point */ 600*24fb275aSStefano Zampini PetscCall(MatMult(J, GradF, W)); 601*24fb275aSStefano Zampini if (has_objective) PetscCall(VecDotRealPart(GradF, W, &gTBg)); 602*24fb275aSStefano Zampini else PetscCall(VecDotRealPart(W, W, &gTBg)); /* B = J^t * J */ 603*24fb275aSStefano Zampini /* Eqs 4.11 and 4.12 in Nocedal and Wright 2nd Edition (4.7 and 4.8 in 1st Edition) */ 604*24fb275aSStefano Zampini auk = delta / gfnorm_k; 605*24fb275aSStefano Zampini if (gTBg < 0.0) tauk = 1.0; 606*24fb275aSStefano Zampini else tauk = PetscMin(gfnorm * gfnorm * gfnorm_k / (delta * gTBg), 1); 607*24fb275aSStefano Zampini auk *= tauk; 608*24fb275aSStefano Zampini ycnorm = auk * gfnorm; 609*24fb275aSStefano Zampini PetscCall(VecAXPBY(Yc, auk, 0.0, GradF)); 610*24fb275aSStefano Zampini 611*24fb275aSStefano Zampini on_boundary = PETSC_FALSE; 612*24fb275aSStefano Zampini if (tauk != 1.0) { 613*24fb275aSStefano Zampini KSPConvergedReason reason; 614*24fb275aSStefano Zampini 6154b0a5c37SStefano Zampini /* sufficient decrease (see 6.3.27 in Conn, Gould, Toint "Trust Region Methods") 616*24fb275aSStefano Zampini beta_k the largest eigenvalue of the Hessian. Here we use the previous estimated value */ 617*24fb275aSStefano Zampini objmin = -neP->kmdc * gnorm * PetscMin(gnorm / beta_k, delta); 618fb01098fSStefano Zampini PetscCall(KSPCGSetObjectiveTarget(snes->ksp, objmin)); 6194b0a5c37SStefano Zampini 620*24fb275aSStefano Zampini /* specify radius if looking for Newton step and trust region norm is the l2 norm */ 621*24fb275aSStefano Zampini PetscCall(KSPCGSetRadius(snes->ksp, neP->fallback == SNES_TR_FALLBACK_NEWTON && neP->norm == NORM_2 ? delta : 0.0)); 622*24fb275aSStefano Zampini PetscCall(KSPSetOperators(snes->ksp, J, Jp)); 6234a221d59SStefano Zampini PetscCall(KSPSolve(snes->ksp, F, Y)); 6244a221d59SStefano Zampini SNESCheckKSPSolve(snes); 6254a221d59SStefano Zampini PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 626*24fb275aSStefano Zampini PetscCall(KSPGetConvergedReason(snes->ksp, &reason)); 627*24fb275aSStefano Zampini on_boundary = (PetscBool)(reason == KSP_CONVERGED_STEP_LENGTH); 6284b0a5c37SStefano Zampini PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits)); 629*24fb275aSStefano Zampini if (neP->kmdc) { /* update estimated Hessian largest eigenvalue */ 630*24fb275aSStefano Zampini PetscReal emax, emin; 631*24fb275aSStefano Zampini PetscCall(KSPComputeExtremeSingularValues(snes->ksp, &emax, &emin)); 632*24fb275aSStefano Zampini if (emax > 0.0) beta_k = emax + 1; 633*24fb275aSStefano Zampini } 634*24fb275aSStefano Zampini } else { /* Cauchy point is on the boundary, accept it */ 635*24fb275aSStefano Zampini on_boundary = PETSC_TRUE; 636*24fb275aSStefano Zampini PetscCall(VecCopy(Yc, Y)); 637*24fb275aSStefano Zampini PetscCall(PetscInfo(snes, "CP evaluated on boundary. delta: %g, ycnorm: %g, gTBg: %g\n", (double)delta, (double)ycnorm, (double)gTBg)); 638*24fb275aSStefano Zampini } 639*24fb275aSStefano Zampini PetscCall(VecNorm(Y, neP->norm, &ynorm)); 6404800dd8cSBarry Smith 6414a221d59SStefano Zampini /* decide what to do when the update is outside of trust region */ 6425ec2728bSStefano Zampini if (ynorm > delta || ynorm == 0.0) { 6435ec2728bSStefano Zampini SNESNewtonTRFallbackType fallback = ynorm > 0.0 ? neP->fallback : SNES_TR_FALLBACK_CAUCHY; 6445ec2728bSStefano Zampini 645*24fb275aSStefano Zampini PetscCheck(neP->norm == NORM_2 || fallback != SNES_TR_FALLBACK_DOGLEG, PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "DOGLEG without l2 norm not implemented"); 6465ec2728bSStefano Zampini switch (fallback) { 6474a221d59SStefano Zampini case SNES_TR_FALLBACK_NEWTON: 6484a221d59SStefano Zampini auk = delta / ynorm; 6494a221d59SStefano Zampini PetscCall(VecScale(Y, auk)); 6504b0a5c37SStefano Zampini PetscCall(PetscInfo(snes, "SN evaluated. delta: %g, ynorm: %g\n", (double)delta, (double)ynorm)); 6514a221d59SStefano Zampini break; 6524a221d59SStefano Zampini case SNES_TR_FALLBACK_CAUCHY: 6534a221d59SStefano Zampini case SNES_TR_FALLBACK_DOGLEG: 6545ec2728bSStefano Zampini if (fallback == SNES_TR_FALLBACK_CAUCHY || gTBg <= 0.0) { 655*24fb275aSStefano Zampini PetscCall(VecCopy(Yc, Y)); 6564a221d59SStefano Zampini PetscCall(PetscInfo(snes, "CP evaluated. delta: %g, ynorm: %g, ycnorm: %g, gTBg: %g\n", (double)delta, (double)ynorm, (double)ycnorm, (double)gTBg)); 6574a221d59SStefano Zampini } else { /* take linear combination of Cauchy and Newton direction and step */ 658*24fb275aSStefano Zampini auk = gfnorm * gfnorm / gTBg; 659*24fb275aSStefano Zampini if (gfnorm_k * auk >= delta) { /* first leg: Cauchy point outside of trust region */ 660*24fb275aSStefano Zampini PetscCall(VecAXPBY(Y, delta / gfnorm_k, 0.0, GradF)); 661*24fb275aSStefano Zampini PetscCall(PetscInfo(snes, "CP evaluated (outside region). delta: %g, ynorm: %g, ycnorm: %g\n", (double)delta, (double)ynorm, (double)ycnorm)); 662*24fb275aSStefano Zampini } else { /* second leg */ 6634a221d59SStefano Zampini PetscReal c0, c1, c2, tau = 0.0, tpos, tneg; 6644a221d59SStefano Zampini PetscBool noroots; 665284fb49fSHeeho Park 666*24fb275aSStefano Zampini /* Find solutions of (Eq. 4.16 in Nocedal and Wright) 667*24fb275aSStefano Zampini ||p_U + lambda * (p_B - p_U)||^2 - delta^2 = 0, 668*24fb275aSStefano Zampini where p_U the Cauchy direction, p_B the Newton direction */ 6694a221d59SStefano Zampini PetscCall(VecAXPBY(YU, auk, 0.0, GradF)); 6704a221d59SStefano Zampini PetscCall(VecAXPY(Y, -1.0, YU)); 6714a221d59SStefano Zampini PetscCall(VecNorm(Y, NORM_2, &c0)); 6724a221d59SStefano Zampini PetscCall(VecDotRealPart(YU, Y, &c1)); 6734a221d59SStefano Zampini c0 = PetscSqr(c0); 6744a221d59SStefano Zampini c2 = PetscSqr(ycnorm) - PetscSqr(delta); 675*24fb275aSStefano Zampini PetscQuadraticRoots(c0, 2 * c1, c2, &tneg, &tpos); 6764a221d59SStefano Zampini 677*24fb275aSStefano Zampini /* In principle the DL strategy as a unique solution in [0,1] 678*24fb275aSStefano Zampini here we check that for some reason we numerically failed 679*24fb275aSStefano Zampini to compute it. In that case, we use the Cauchy point */ 6804a221d59SStefano Zampini noroots = PetscIsInfOrNanReal(tneg); 681*24fb275aSStefano Zampini if (!noroots) { 682*24fb275aSStefano Zampini if (tpos > 1) { 683*24fb275aSStefano Zampini if (tneg >= 0 && tneg <= 1) { 684*24fb275aSStefano Zampini tau = tneg; 685*24fb275aSStefano Zampini } else noroots = PETSC_TRUE; 686*24fb275aSStefano Zampini } else if (tpos >= 0) { 687*24fb275aSStefano Zampini tau = tpos; 688*24fb275aSStefano Zampini } else noroots = PETSC_TRUE; 689*24fb275aSStefano Zampini } 6904a221d59SStefano Zampini if (noroots) { /* No roots, select Cauchy point */ 691*24fb275aSStefano Zampini PetscCall(VecCopy(Yc, Y)); 6924a221d59SStefano Zampini } else { 693*24fb275aSStefano Zampini PetscCall(VecAXPBY(Y, 1.0, tau, YU)); 6944a221d59SStefano Zampini } 695*24fb275aSStefano Zampini PetscCall(PetscInfo(snes, "%s evaluated. roots: (%g, %g), tau %g, ynorm: %g, ycnorm: %g, gTBg: %g\n", noroots ? "CP" : "DL", (double)tneg, (double)tpos, (double)tau, (double)ynorm, (double)ycnorm, (double)gTBg)); 6964a221d59SStefano Zampini } 6974a221d59SStefano Zampini } 6984a221d59SStefano Zampini break; 6994a221d59SStefano Zampini default: 7004a221d59SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "Unknown fallback mode"); 701454a90a3SBarry Smith break; 70252392280SLois Curfman McInnes } 7034800dd8cSBarry Smith } 7044a221d59SStefano Zampini 7057aa289d8SStefano Zampini /* compute the quadratic model difference */ 706*24fb275aSStefano Zampini PetscCall(SNESNewtonTRQuadraticDelta(snes, J, has_objective, Y, GradF, W, &yTHy, &gTy, &deltaqm)); 7074a221d59SStefano Zampini 7084a221d59SStefano Zampini /* Compute new objective function */ 709*24fb275aSStefano Zampini PetscCall(SNESNewtonTRObjective(snes, has_objective, X, Y, W, G, &gnorm, &fkp1)); 710*24fb275aSStefano Zampini if (PetscIsInfOrNanReal(fkp1)) rho = neP->eta1; 711*24fb275aSStefano Zampini else { 7124a221d59SStefano Zampini if (deltaqm > 0.0) rho = (fk - fkp1) / deltaqm; /* actual improvement over predicted improvement */ 713*24fb275aSStefano Zampini else rho = neP->eta1; /* no reduction in quadratic model, step must be rejected */ 714*24fb275aSStefano Zampini } 7154a221d59SStefano Zampini 716*24fb275aSStefano Zampini PetscCall(VecNorm(Y, neP->norm, &ynorm)); 717*24fb275aSStefano Zampini PetscCall(PetscInfo(snes, "rho=%g, delta=%g, fk=%g, fkp1=%g, deltaqm=%g, gTy=%g, yTHy=%g, ynormk=%g\n", (double)rho, (double)delta, (double)fk, (double)fkp1, (double)deltaqm, (double)gTy, (double)yTHy, (double)ynorm)); 718*24fb275aSStefano Zampini 719*24fb275aSStefano Zampini /* update the size of the trust region */ 7204a221d59SStefano Zampini if (rho < neP->eta2) delta *= neP->t1; /* shrink the region */ 721*24fb275aSStefano Zampini else if (rho > neP->eta3 && on_boundary) delta *= neP->t2; /* expand the region */ 7224a221d59SStefano Zampini delta = PetscMin(delta, deltaM); /* but not greater than deltaM */ 7234a221d59SStefano Zampini 724*24fb275aSStefano Zampini /* log 2-norm of update for moniroting routines */ 725*24fb275aSStefano Zampini PetscCall(VecNorm(Y, NORM_2, &ynorm)); 726*24fb275aSStefano Zampini 727*24fb275aSStefano Zampini /* decide on new step */ 7284a221d59SStefano Zampini neP->delta = delta; 729*24fb275aSStefano Zampini if (rho > neP->eta1) { 7304a221d59SStefano Zampini rho_satisfied = PETSC_TRUE; 7314a221d59SStefano Zampini } else { 7324a221d59SStefano Zampini rho_satisfied = PETSC_FALSE; 7334a221d59SStefano Zampini PetscCall(PetscInfo(snes, "Trying again in smaller region\n")); 7344a221d59SStefano Zampini /* check to see if progress is hopeless */ 7354a221d59SStefano Zampini PetscCall(SNESTR_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP)); 7362d157150SStefano Zampini if (!snes->reason) PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm)); 7374b0a5c37SStefano Zampini if (snes->reason == SNES_CONVERGED_SNORM_RELATIVE) snes->reason = SNES_DIVERGED_TR_DELTA; 7384a221d59SStefano Zampini snes->numFailures++; 7394a221d59SStefano Zampini /* We're not progressing, so return with the current iterate */ 7404a221d59SStefano Zampini if (snes->reason) break; 7414a221d59SStefano Zampini } 7424a221d59SStefano Zampini if (rho_satisfied) { 7434a221d59SStefano Zampini /* Update function values */ 7444a221d59SStefano Zampini already_done = PETSC_FALSE; 7454800dd8cSBarry Smith fnorm = gnorm; 7464a221d59SStefano Zampini fk = fkp1; 7474a221d59SStefano Zampini 7484a221d59SStefano Zampini /* New residual and linearization point */ 7499566063dSJacob Faibussowitsch PetscCall(VecCopy(G, F)); 7509566063dSJacob Faibussowitsch PetscCall(VecCopy(W, X)); 7514a221d59SStefano Zampini 75285385478SLisandro Dalcin /* Monitor convergence */ 7539566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 7544a221d59SStefano Zampini snes->iter++; 755fbe28522SBarry Smith snes->norm = fnorm; 756c1e67a49SFande Kong snes->xnorm = xnorm; 757c1e67a49SFande Kong snes->ynorm = ynorm; 7589566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 7599566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 7604a221d59SStefano Zampini 76185385478SLisandro Dalcin /* Test for convergence, xnorm = || X || */ 7624a221d59SStefano Zampini PetscCall(VecNorm(X, NORM_2, &xnorm)); 7632d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm)); 7642d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 7654a221d59SStefano Zampini if (snes->reason) break; 7664a221d59SStefano Zampini } 76738442cffSBarry Smith } 768284fb49fSHeeho Park 7694a221d59SStefano Zampini if (clear_converged_test) { 770a0254a93SStefano Zampini PetscCall(KSPGetAndClearConvergenceTest(snes->ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 7719566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 772a0254a93SStefano Zampini PetscCall(KSPSetConvergenceTest(snes->ksp, convtest, convctx, convdestroy)); 7735e28bcb6Sprj- } 7743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7754800dd8cSBarry Smith } 776284fb49fSHeeho Park 777d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes) 778d71ae5a4SJacob Faibussowitsch { 7793a40ed3dSBarry Smith PetscFunctionBegin; 780*24fb275aSStefano Zampini PetscCall(SNESSetWorkVecs(snes, 5)); 7819566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 7823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7834800dd8cSBarry Smith } 7846b8b9a38SLisandro Dalcin 785*24fb275aSStefano Zampini static PetscErrorCode SNESReset_NEWTONTR(SNES snes) 786*24fb275aSStefano Zampini { 787*24fb275aSStefano Zampini SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 788*24fb275aSStefano Zampini 789*24fb275aSStefano Zampini PetscFunctionBegin; 790*24fb275aSStefano Zampini PetscCall(MatDestroy(&tr->qnB)); 791*24fb275aSStefano Zampini PetscCall(MatDestroy(&tr->qnB_pre)); 792*24fb275aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 793*24fb275aSStefano Zampini } 794*24fb275aSStefano Zampini 795d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes) 796d71ae5a4SJacob Faibussowitsch { 7973a40ed3dSBarry Smith PetscFunctionBegin; 798*24fb275aSStefano Zampini PetscCall(SNESReset_NEWTONTR(snes)); 7999566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 8003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8014800dd8cSBarry Smith } 8024800dd8cSBarry Smith 803d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NEWTONTR(SNES snes, PetscOptionItems *PetscOptionsObject) 804d71ae5a4SJacob Faibussowitsch { 80504d7464bSBarry Smith SNES_NEWTONTR *ctx = (SNES_NEWTONTR *)snes->data; 806*24fb275aSStefano Zampini SNESNewtonTRQNType qn; 807*24fb275aSStefano Zampini SNESNewtonTRFallbackType fallback; 808*24fb275aSStefano Zampini NormType norm; 809*24fb275aSStefano Zampini PetscReal deltatol; 810*24fb275aSStefano Zampini PetscBool flg; 8114800dd8cSBarry Smith 8123a40ed3dSBarry Smith PetscFunctionBegin; 813d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations"); 8144a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_eta1", "eta1", "None", ctx->eta1, &ctx->eta1, NULL)); 8154a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_eta2", "eta2", "None", ctx->eta2, &ctx->eta2, NULL)); 8164a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_eta3", "eta3", "None", ctx->eta3, &ctx->eta3, NULL)); 8174a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_t1", "t1", "None", ctx->t1, &ctx->t1, NULL)); 8184a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_t2", "t2", "None", ctx->t2, &ctx->t2, NULL)); 8194a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_deltaM", "deltaM", "None", ctx->deltaM, &ctx->deltaM, NULL)); 8209566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_tr_delta0", "delta0", "None", ctx->delta0, &ctx->delta0, NULL)); 8214b0a5c37SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_kmdc", "sufficient decrease parameter", "None", ctx->kmdc, &ctx->kmdc, NULL)); 822*24fb275aSStefano Zampini 823*24fb275aSStefano Zampini deltatol = snes->deltatol; 824*24fb275aSStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_tol", "Trust region tolerance", "SNESSetTrustRegionTolerance", deltatol, &deltatol, &flg)); 825*24fb275aSStefano Zampini if (flg) PetscCall(SNESSetTrustRegionTolerance(snes, deltatol)); 826*24fb275aSStefano Zampini 827*24fb275aSStefano Zampini fallback = ctx->fallback; 828*24fb275aSStefano Zampini PetscCall(PetscOptionsEnum("-snes_tr_fallback_type", "Type of fallback if subproblem solution is outside of the trust region", "SNESNewtonTRSetFallbackType", SNESNewtonTRFallbackTypes, (PetscEnum)fallback, (PetscEnum *)&fallback, &flg)); 829*24fb275aSStefano Zampini if (flg) PetscCall(SNESNewtonTRSetFallbackType(snes, fallback)); 830*24fb275aSStefano Zampini 831*24fb275aSStefano Zampini qn = ctx->qn; 832*24fb275aSStefano Zampini PetscCall(PetscOptionsEnum("-snes_tr_qn", "Use Quasi-Newton approximations for the model", "SNESNewtonTRSetQNType", SNESNewtonTRQNTypes, (PetscEnum)qn, (PetscEnum *)&qn, &flg)); 833*24fb275aSStefano Zampini if (flg) PetscCall(SNESNewtonTRSetQNType(snes, qn)); 834*24fb275aSStefano Zampini 835*24fb275aSStefano Zampini norm = ctx->norm; 836*24fb275aSStefano Zampini PetscCall(PetscOptionsEnum("-snes_tr_norm_type", "Type of norm for trust region bounds", "SNESNewtonTRSetNormType", NormTypes, (PetscEnum)norm, (PetscEnum *)&norm, &flg)); 837*24fb275aSStefano Zampini if (flg) PetscCall(SNESNewtonTRSetNormType(snes, norm)); 838*24fb275aSStefano Zampini 839d0609cedSBarry Smith PetscOptionsHeadEnd(); 8403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8414800dd8cSBarry Smith } 8424800dd8cSBarry Smith 843d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NEWTONTR(SNES snes, PetscViewer viewer) 844d71ae5a4SJacob Faibussowitsch { 84504d7464bSBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 846ace3abfcSBarry Smith PetscBool iascii; 847a935fc98SLois Curfman McInnes 8483a40ed3dSBarry Smith PetscFunctionBegin; 8499566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 85032077d6dSBarry Smith if (iascii) { 8514a221d59SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " Trust region tolerance %g\n", (double)snes->deltatol)); 8524a221d59SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " eta1=%g, eta2=%g, eta3=%g\n", (double)tr->eta1, (double)tr->eta2, (double)tr->eta3)); 8534a221d59SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " delta0=%g, t1=%g, t2=%g, deltaM=%g\n", (double)tr->delta0, (double)tr->t1, (double)tr->t2, (double)tr->deltaM)); 8544b0a5c37SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " kmdc=%g\n", (double)tr->kmdc)); 8554a221d59SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " fallback=%s\n", SNESNewtonTRFallbackTypes[tr->fallback])); 856*24fb275aSStefano Zampini if (tr->qn) PetscCall(PetscViewerASCIIPrintf(viewer, " qn=%s\n", SNESNewtonTRQNTypes[tr->qn])); 857*24fb275aSStefano Zampini if (tr->norm != NORM_2) PetscCall(PetscViewerASCIIPrintf(viewer, " norm=%s\n", NormTypes[tr->norm])); 85819bcc07fSBarry Smith } 8593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 860a935fc98SLois Curfman McInnes } 861f6dfbefdSBarry Smith 862ebe3b25bSBarry Smith /*MC 8634a221d59SStefano Zampini SNESNEWTONTR - Newton based nonlinear solver that uses trust-region dogleg method with Cauchy direction 864f6dfbefdSBarry Smith 865f6dfbefdSBarry Smith Options Database Keys: 8664a221d59SStefano Zampini + -snes_tr_tol <tol> - trust region tolerance 867*24fb275aSStefano Zampini . -snes_tr_eta1 <eta1> - trust region parameter eta1 <= eta2, rho > eta1 breaks out of the inner iteration (default: eta1=0.001) 868*24fb275aSStefano Zampini . -snes_tr_eta2 <eta2> - trust region parameter, rho <= eta2 shrinks the trust region (default: eta2=0.25) 8694a221d59SStefano Zampini . -snes_tr_eta3 <eta3> - trust region parameter eta3 > eta2, rho >= eta3 expands the trust region (default: eta3=0.75) 8704a221d59SStefano Zampini . -snes_tr_t1 <t1> - trust region parameter, shrinking factor of trust region (default: 0.25) 8714a221d59SStefano Zampini . -snes_tr_t2 <t2> - trust region parameter, expanding factor of trust region (default: 2.0) 8724a221d59SStefano Zampini . -snes_tr_deltaM <deltaM> - trust region parameter, max size of trust region (default: MAX_REAL) 8734a221d59SStefano Zampini . -snes_tr_delta0 <delta0> - trust region parameter, initial size of trust region (default: 0.2) 8744a221d59SStefano Zampini - -snes_tr_fallback_type <newton,cauchy,dogleg> - Solution strategy to test reduction when step is outside of trust region. Can use scaled Newton direction, Cauchy point (Steepest Descent direction) or dogleg method. 875ebe3b25bSBarry Smith 876f6dfbefdSBarry Smith Reference: 8774a221d59SStefano Zampini . * - "Numerical Optimization" by Nocedal and Wright, chapter 4. 878ebe3b25bSBarry Smith 879420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESSetTrustRegionTolerance()`, 8804a221d59SStefano Zampini `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`, 881*24fb275aSStefano Zampini `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRSetFallbackType()`, `SNESNewtonTRSetQNType()` 882ebe3b25bSBarry Smith M*/ 883d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes) 884d71ae5a4SJacob Faibussowitsch { 88504d7464bSBarry Smith SNES_NEWTONTR *neP; 8864800dd8cSBarry Smith 8873a40ed3dSBarry Smith PetscFunctionBegin; 88804d7464bSBarry Smith snes->ops->setup = SNESSetUp_NEWTONTR; 88904d7464bSBarry Smith snes->ops->solve = SNESSolve_NEWTONTR; 890*24fb275aSStefano Zampini snes->ops->reset = SNESReset_NEWTONTR; 89104d7464bSBarry Smith snes->ops->destroy = SNESDestroy_NEWTONTR; 89204d7464bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR; 89304d7464bSBarry Smith snes->ops->view = SNESView_NEWTONTR; 894fbe28522SBarry Smith 8951a58d6d3SStefano Zampini snes->stol = 0.0; 89642f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 8974b0a5c37SStefano Zampini snes->npcside = PC_RIGHT; 8984b0a5c37SStefano Zampini snes->usesnpc = PETSC_TRUE; 89942f4f86dSBarry Smith 9004fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 9014fc747eaSLawrence Mitchell 9024dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&neP)); 903fbe28522SBarry Smith snes->data = (void *)neP; 904fbe28522SBarry Smith neP->delta = 0.0; 905fbe28522SBarry Smith neP->delta0 = 0.2; 9064a221d59SStefano Zampini neP->eta1 = 0.001; 9074a221d59SStefano Zampini neP->eta2 = 0.25; 9084a221d59SStefano Zampini neP->eta3 = 0.75; 9094a221d59SStefano Zampini neP->t1 = 0.25; 9104a221d59SStefano Zampini neP->t2 = 2.0; 9114a221d59SStefano Zampini neP->deltaM = 1.e10; 912*24fb275aSStefano Zampini neP->norm = NORM_2; 9134a221d59SStefano Zampini neP->fallback = SNES_TR_FALLBACK_NEWTON; 9144b0a5c37SStefano Zampini neP->kmdc = 0.0; /* by default do not use sufficient decrease */ 9153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9164800dd8cSBarry Smith } 917