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}; 1124fb275aSStefano Zampini const char *const SNESNewtonTRQNTypes[] = {"NONE", "SAME", "DIFFERENT", "SNESNewtonTRQNType", "SNES_TR_QN_", NULL}; 1224fb275aSStefano Zampini 1324fb275aSStefano Zampini static PetscErrorCode SNESComputeJacobian_MATLMVM(SNES snes, Vec X, Mat J, Mat B, void *dummy) 1424fb275aSStefano Zampini { 1524fb275aSStefano Zampini PetscFunctionBegin; 1624fb275aSStefano Zampini // PetscCall(MatLMVMSymBroydenSetDelta(B, _some_delta)); 1724fb275aSStefano Zampini PetscCall(MatLMVMUpdate(B, X, snes->vec_func)); 1824fb275aSStefano Zampini PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 1924fb275aSStefano Zampini PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 2024fb275aSStefano Zampini if (J != B) { 2124fb275aSStefano Zampini // PetscCall(MatLMVMSymBroydenSetDelta(J, _some_delta)); 2224fb275aSStefano Zampini PetscCall(MatLMVMUpdate(J, X, snes->vec_func)); 2324fb275aSStefano Zampini PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 2424fb275aSStefano Zampini PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 2524fb275aSStefano Zampini } 2624fb275aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 2724fb275aSStefano 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; 38a935fc98SLois Curfman McInnes /* Determine norm of solution */ 399566063dSJacob Faibussowitsch PetscCall(KSPBuildSolution(ksp, NULL, &x)); 4024fb275aSStefano Zampini PetscCall(VecNorm(x, neP->norm, &nrm)); 41064f8208SBarry Smith if (nrm >= neP->delta) { 426d51442aSBarry Smith PetscCall(PetscInfo(snes, "Ending linear iteration early due to exiting trust region, delta=%g, length=%g\n", (double)neP->delta, (double)nrm)); 43329f5518SBarry Smith *reason = KSP_CONVERGED_STEP_LENGTH; 446d51442aSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 45df60cc22SBarry Smith } 466d51442aSBarry Smith PetscCall((*ctx->convtest)(ksp, n, rnorm, reason, ctx->convctx)); 476d51442aSBarry Smith if (*reason) PetscCall(PetscInfo(snes, "Default or user provided convergence test KSP iterations=%" PetscInt_FMT ", rnorm=%g\n", n, (double)rnorm)); 483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49df60cc22SBarry Smith } 5082bf6240SBarry Smith 51d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTR_KSPConverged_Destroy(void *cctx) 52d71ae5a4SJacob Faibussowitsch { 53971273eeSBarry Smith SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx *)cctx; 54971273eeSBarry Smith 55971273eeSBarry Smith PetscFunctionBegin; 569566063dSJacob Faibussowitsch PetscCall((*ctx->convdestroy)(ctx->convctx)); 579566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 59971273eeSBarry Smith } 60971273eeSBarry Smith 61d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTR_Converged_Private(SNES snes, PetscInt it, PetscReal xnorm, PetscReal pnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy) 62d71ae5a4SJacob Faibussowitsch { 6304d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 6485385478SLisandro Dalcin 6585385478SLisandro Dalcin PetscFunctionBegin; 6685385478SLisandro Dalcin *reason = SNES_CONVERGED_ITERATING; 674a221d59SStefano Zampini if (neP->delta < snes->deltatol) { 684a221d59SStefano Zampini PetscCall(PetscInfo(snes, "Diverged due to too small a trust region %g<%g\n", (double)neP->delta, (double)snes->deltatol)); 691c6b2ff8SBarry Smith *reason = SNES_DIVERGED_TR_DELTA; 70e71169deSBarry Smith } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 7163a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n", snes->max_funcs)); 7285385478SLisandro Dalcin *reason = SNES_DIVERGED_FUNCTION_COUNT; 7385385478SLisandro Dalcin } 743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7585385478SLisandro Dalcin } 7685385478SLisandro Dalcin 774a221d59SStefano Zampini /*@ 7824fb275aSStefano Zampini SNESNewtonTRSetNormType - Specify the type of norm to use for the computation of the trust region. 7924fb275aSStefano Zampini 8024fb275aSStefano Zampini Input Parameters: 8124fb275aSStefano Zampini + snes - the nonlinear solver object 8224fb275aSStefano Zampini - norm - the norm type 8324fb275aSStefano Zampini 8424fb275aSStefano Zampini Level: intermediate 8524fb275aSStefano Zampini 8624fb275aSStefano Zampini .seealso: `SNESNEWTONTR`, `NormType` 8724fb275aSStefano Zampini @*/ 8824fb275aSStefano Zampini PetscErrorCode SNESNewtonTRSetNormType(SNES snes, NormType norm) 8924fb275aSStefano Zampini { 9024fb275aSStefano Zampini PetscBool flg; 9124fb275aSStefano Zampini 9224fb275aSStefano Zampini PetscFunctionBegin; 9324fb275aSStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 9424fb275aSStefano Zampini PetscValidLogicalCollectiveEnum(snes, norm, 2); 9524fb275aSStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 9624fb275aSStefano Zampini if (flg) { 9724fb275aSStefano Zampini SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 9824fb275aSStefano Zampini 9924fb275aSStefano Zampini tr->norm = norm; 10024fb275aSStefano Zampini } 10124fb275aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 10224fb275aSStefano Zampini } 10324fb275aSStefano Zampini 10424fb275aSStefano Zampini /*@ 10524fb275aSStefano Zampini SNESNewtonTRSetQNType - Specify to use a quasi-Newton model. 10624fb275aSStefano Zampini 10724fb275aSStefano Zampini Input Parameters: 10824fb275aSStefano Zampini + snes - the nonlinear solver object 10924fb275aSStefano Zampini - use - the type of approximations to be used 11024fb275aSStefano Zampini 11124fb275aSStefano Zampini Level: intermediate 11224fb275aSStefano Zampini 11324fb275aSStefano Zampini Notes: 11424fb275aSStefano Zampini Options for the approximations can be set with the snes_tr_qn_ and snes_tr_qn_pre_ prefixes. 11524fb275aSStefano Zampini 11624fb275aSStefano Zampini .seealso: `SNESNEWTONTR`, `SNESNewtonTRQNType`, `MATLMVM` 11724fb275aSStefano Zampini @*/ 11824fb275aSStefano Zampini PetscErrorCode SNESNewtonTRSetQNType(SNES snes, SNESNewtonTRQNType use) 11924fb275aSStefano Zampini { 12024fb275aSStefano Zampini PetscBool flg; 12124fb275aSStefano Zampini 12224fb275aSStefano Zampini PetscFunctionBegin; 12324fb275aSStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 12424fb275aSStefano Zampini PetscValidLogicalCollectiveEnum(snes, use, 2); 12524fb275aSStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 12624fb275aSStefano Zampini if (flg) { 12724fb275aSStefano Zampini SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 12824fb275aSStefano Zampini 12924fb275aSStefano Zampini tr->qn = use; 13024fb275aSStefano Zampini } 13124fb275aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 13224fb275aSStefano Zampini } 13324fb275aSStefano Zampini 13424fb275aSStefano Zampini /*@ 135420bcc1bSBarry Smith SNESNewtonTRSetFallbackType - Set the type of fallback to use if the solution of the trust region subproblem is outside the radius 1364a221d59SStefano Zampini 1374a221d59SStefano Zampini Input Parameters: 1384a221d59SStefano Zampini + snes - the nonlinear solver object 1394a221d59SStefano Zampini - ftype - the fallback type, see `SNESNewtonTRFallbackType` 1404a221d59SStefano Zampini 1414a221d59SStefano Zampini Level: intermediate 1424a221d59SStefano Zampini 143420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPreCheck()`, 1444a221d59SStefano Zampini `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()` 1454a221d59SStefano Zampini @*/ 1464a221d59SStefano Zampini PetscErrorCode SNESNewtonTRSetFallbackType(SNES snes, SNESNewtonTRFallbackType ftype) 1474a221d59SStefano Zampini { 1484a221d59SStefano Zampini SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 1494a221d59SStefano Zampini PetscBool flg; 1504a221d59SStefano Zampini 1514a221d59SStefano Zampini PetscFunctionBegin; 1524a221d59SStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1534a221d59SStefano Zampini PetscValidLogicalCollectiveEnum(snes, ftype, 2); 1544a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 1554a221d59SStefano Zampini if (flg) tr->fallback = ftype; 1564a221d59SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 1574a221d59SStefano Zampini } 1584a221d59SStefano Zampini 1597cb011f5SBarry Smith /*@C 160c9368356SGlenn Hammond SNESNewtonTRSetPreCheck - Sets a user function that is called before the search step has been determined. 1614a221d59SStefano Zampini Allows the user a chance to change or override the trust region decision. 162f6dfbefdSBarry Smith 163c3339decSBarry Smith Logically Collective 164c9368356SGlenn Hammond 165c9368356SGlenn Hammond Input Parameters: 166c9368356SGlenn Hammond + snes - the nonlinear solver object 16720f4b53cSBarry Smith . func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPreCheck()` 16820f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 169c9368356SGlenn Hammond 170acba1f63SStefano Zampini Level: intermediate 171c9368356SGlenn Hammond 172f6dfbefdSBarry Smith Note: 1734a221d59SStefano Zampini This function is called BEFORE the function evaluation within the solver. 174c9368356SGlenn Hammond 175420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`, 176c9368356SGlenn Hammond @*/ 177d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, PetscBool *, void *), void *ctx) 178d71ae5a4SJacob Faibussowitsch { 179c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 1804a221d59SStefano Zampini PetscBool flg; 181c9368356SGlenn Hammond 182c9368356SGlenn Hammond PetscFunctionBegin; 183c9368356SGlenn Hammond PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1844a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 1854a221d59SStefano Zampini if (flg) { 186c9368356SGlenn Hammond if (func) tr->precheck = func; 187c9368356SGlenn Hammond if (ctx) tr->precheckctx = ctx; 1884a221d59SStefano Zampini } 1893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 190c9368356SGlenn Hammond } 191c9368356SGlenn Hammond 192c9368356SGlenn Hammond /*@C 193c9368356SGlenn Hammond SNESNewtonTRGetPreCheck - Gets the pre-check function 194c9368356SGlenn Hammond 19520f4b53cSBarry Smith Not Collective 196c9368356SGlenn Hammond 197c9368356SGlenn Hammond Input Parameter: 198c9368356SGlenn Hammond . snes - the nonlinear solver context 199c9368356SGlenn Hammond 200c9368356SGlenn Hammond Output Parameters: 20120f4b53cSBarry Smith + func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPreCheck()` 20220f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 203c9368356SGlenn Hammond 204acba1f63SStefano Zampini Level: intermediate 205c9368356SGlenn Hammond 206420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRPreCheck()` 207c9368356SGlenn Hammond @*/ 208d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, PetscBool *, void *), void **ctx) 209d71ae5a4SJacob Faibussowitsch { 210c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 2114a221d59SStefano Zampini PetscBool flg; 212c9368356SGlenn Hammond 213c9368356SGlenn Hammond PetscFunctionBegin; 214c9368356SGlenn Hammond PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 2154a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 2164a221d59SStefano Zampini PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 217c9368356SGlenn Hammond if (func) *func = tr->precheck; 218c9368356SGlenn Hammond if (ctx) *ctx = tr->precheckctx; 2193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 220c9368356SGlenn Hammond } 221c9368356SGlenn Hammond 222c9368356SGlenn Hammond /*@C 2237cb011f5SBarry Smith SNESNewtonTRSetPostCheck - Sets a user function that is called after the search step has been determined but before the next 2244a221d59SStefano Zampini function evaluation. Allows the user a chance to change or override the internal decision of the solver 225f6dfbefdSBarry Smith 226c3339decSBarry Smith Logically Collective 2277cb011f5SBarry Smith 2287cb011f5SBarry Smith Input Parameters: 2297cb011f5SBarry Smith + snes - the nonlinear solver object 23020f4b53cSBarry Smith . func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPostCheck()` 23120f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 2327cb011f5SBarry Smith 233acba1f63SStefano Zampini Level: intermediate 2347cb011f5SBarry Smith 235f6dfbefdSBarry Smith Note: 2364a221d59SStefano Zampini This function is called BEFORE the function evaluation within the solver while the function set in 237f6dfbefdSBarry Smith `SNESLineSearchSetPostCheck()` is called AFTER the function evaluation. 2387cb011f5SBarry Smith 239420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRPostCheck()`, `SNESNewtonTRGetPostCheck()`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()` 2407cb011f5SBarry Smith @*/ 241d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx) 242d71ae5a4SJacob Faibussowitsch { 2437cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 2444a221d59SStefano Zampini PetscBool flg; 2457cb011f5SBarry Smith 2467cb011f5SBarry Smith PetscFunctionBegin; 2477cb011f5SBarry Smith PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 2484a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 2494a221d59SStefano Zampini if (flg) { 2507cb011f5SBarry Smith if (func) tr->postcheck = func; 2517cb011f5SBarry Smith if (ctx) tr->postcheckctx = ctx; 2524a221d59SStefano Zampini } 2533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2547cb011f5SBarry Smith } 2557cb011f5SBarry Smith 2567cb011f5SBarry Smith /*@C 2577cb011f5SBarry Smith SNESNewtonTRGetPostCheck - Gets the post-check function 2587cb011f5SBarry Smith 25920f4b53cSBarry Smith Not Collective 2607cb011f5SBarry Smith 2617cb011f5SBarry Smith Input Parameter: 2627cb011f5SBarry Smith . snes - the nonlinear solver context 2637cb011f5SBarry Smith 2647cb011f5SBarry Smith Output Parameters: 26520f4b53cSBarry Smith + func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPostCheck()` 26620f4b53cSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 2677cb011f5SBarry Smith 2687cb011f5SBarry Smith Level: intermediate 2697cb011f5SBarry Smith 270420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRPostCheck()` 2717cb011f5SBarry Smith @*/ 272d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) 273d71ae5a4SJacob Faibussowitsch { 2747cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 2754a221d59SStefano Zampini PetscBool flg; 2767cb011f5SBarry Smith 2777cb011f5SBarry Smith PetscFunctionBegin; 2787cb011f5SBarry Smith PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 2794a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 2804a221d59SStefano Zampini PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 2817cb011f5SBarry Smith if (func) *func = tr->postcheck; 2827cb011f5SBarry Smith if (ctx) *ctx = tr->postcheckctx; 2833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2847cb011f5SBarry Smith } 2857cb011f5SBarry Smith 2867cb011f5SBarry Smith /*@C 2874a221d59SStefano Zampini SNESNewtonTRPreCheck - Runs the precheck routine 288c9368356SGlenn Hammond 289c3339decSBarry Smith Logically Collective 290c9368356SGlenn Hammond 291c9368356SGlenn Hammond Input Parameters: 292c9368356SGlenn Hammond + snes - the solver 293c9368356SGlenn Hammond . X - The last solution 294c9368356SGlenn Hammond - Y - The step direction 295c9368356SGlenn Hammond 2962fe279fdSBarry Smith Output Parameter: 2972fe279fdSBarry Smith . changed_Y - Indicator that the step direction `Y` has been changed. 298c9368356SGlenn Hammond 2994a221d59SStefano Zampini Level: intermediate 300c9368356SGlenn Hammond 301420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRPostCheck()` 302c9368356SGlenn Hammond @*/ 3034a221d59SStefano Zampini PetscErrorCode SNESNewtonTRPreCheck(SNES snes, Vec X, Vec Y, PetscBool *changed_Y) 304d71ae5a4SJacob Faibussowitsch { 305c9368356SGlenn Hammond SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 3064a221d59SStefano Zampini PetscBool flg; 307c9368356SGlenn Hammond 308c9368356SGlenn Hammond PetscFunctionBegin; 3094a221d59SStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 3104a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 3114a221d59SStefano Zampini PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 312c9368356SGlenn Hammond *changed_Y = PETSC_FALSE; 313c9368356SGlenn Hammond if (tr->precheck) { 3149566063dSJacob Faibussowitsch PetscCall((*tr->precheck)(snes, X, Y, changed_Y, tr->precheckctx)); 315c9368356SGlenn Hammond PetscValidLogicalCollectiveBool(snes, *changed_Y, 4); 316c9368356SGlenn Hammond } 3173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 318c9368356SGlenn Hammond } 319c9368356SGlenn Hammond 320c9368356SGlenn Hammond /*@C 3214a221d59SStefano Zampini SNESNewtonTRPostCheck - Runs the postcheck routine 3227cb011f5SBarry Smith 323c3339decSBarry Smith Logically Collective 3247cb011f5SBarry Smith 3257cb011f5SBarry Smith Input Parameters: 3266b867d5aSJose E. Roman + snes - the solver 3276b867d5aSJose E. Roman . X - The last solution 3287cb011f5SBarry Smith . Y - The full step direction 3293312a946SBarry Smith - W - The updated solution, W = X - Y 3307cb011f5SBarry Smith 3317cb011f5SBarry Smith Output Parameters: 3323312a946SBarry Smith + changed_Y - indicator if step has been changed 3333312a946SBarry Smith - changed_W - Indicator if the new candidate solution W has been changed. 3347cb011f5SBarry Smith 335f6dfbefdSBarry Smith Note: 3363312a946SBarry Smith If Y is changed then W is recomputed as X - Y 3377cb011f5SBarry Smith 3384a221d59SStefano Zampini Level: intermediate 3397cb011f5SBarry Smith 340420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`, `SNESNewtonTRPreCheck()` 3417cb011f5SBarry Smith @*/ 3424a221d59SStefano Zampini PetscErrorCode SNESNewtonTRPostCheck(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) 343d71ae5a4SJacob Faibussowitsch { 3447cb011f5SBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 3454a221d59SStefano Zampini PetscBool flg; 3467cb011f5SBarry Smith 3477cb011f5SBarry Smith PetscFunctionBegin; 3484a221d59SStefano Zampini PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 3494a221d59SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 3504a221d59SStefano Zampini PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 351c9368356SGlenn Hammond *changed_Y = PETSC_FALSE; 3527cb011f5SBarry Smith *changed_W = PETSC_FALSE; 3537cb011f5SBarry Smith if (tr->postcheck) { 3549566063dSJacob Faibussowitsch PetscCall((*tr->postcheck)(snes, X, Y, W, changed_Y, changed_W, tr->postcheckctx)); 355c9368356SGlenn Hammond PetscValidLogicalCollectiveBool(snes, *changed_Y, 5); 3567cb011f5SBarry Smith PetscValidLogicalCollectiveBool(snes, *changed_W, 6); 3577cb011f5SBarry Smith } 3583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3597cb011f5SBarry Smith } 36085385478SLisandro Dalcin 36124fb275aSStefano Zampini /* stable implementation of roots of a*x^2 + b*x + c = 0 */ 3624a221d59SStefano Zampini static inline void PetscQuadraticRoots(PetscReal a, PetscReal b, PetscReal c, PetscReal *xm, PetscReal *xp) 3634a221d59SStefano Zampini { 3644a221d59SStefano Zampini PetscReal temp = -0.5 * (b + PetscCopysignReal(1.0, b) * PetscSqrtReal(b * b - 4 * a * c)); 3654a221d59SStefano Zampini PetscReal x1 = temp / a; 3664a221d59SStefano Zampini PetscReal x2 = c / temp; 3674a221d59SStefano Zampini *xm = PetscMin(x1, x2); 3684a221d59SStefano Zampini *xp = PetscMax(x1, x2); 3694a221d59SStefano Zampini } 3704a221d59SStefano Zampini 3717aa289d8SStefano Zampini /* Computes the quadratic model difference */ 37224fb275aSStefano Zampini static PetscErrorCode SNESNewtonTRQuadraticDelta(SNES snes, Mat J, PetscBool has_objective, Vec Y, Vec GradF, Vec W, PetscReal *yTHy_, PetscReal *gTy_, PetscReal *deltaqm) 3737aa289d8SStefano Zampini { 37424fb275aSStefano Zampini PetscReal yTHy, gTy; 37524fb275aSStefano Zampini 3767aa289d8SStefano Zampini PetscFunctionBegin; 37724fb275aSStefano Zampini PetscCall(MatMult(J, Y, W)); 37824fb275aSStefano Zampini if (has_objective) PetscCall(VecDotRealPart(Y, W, &yTHy)); 37924fb275aSStefano Zampini else PetscCall(VecDotRealPart(W, W, &yTHy)); /* Gauss-Newton approximation J^t * J */ 38024fb275aSStefano Zampini PetscCall(VecDotRealPart(GradF, Y, &gTy)); 38124fb275aSStefano Zampini *deltaqm = -(-(gTy) + 0.5 * (yTHy)); /* difference in quadratic model, -gTy because SNES solves it this way */ 38224fb275aSStefano Zampini if (yTHy_) *yTHy_ = yTHy; 38324fb275aSStefano Zampini if (gTy_) *gTy_ = gTy; 38424fb275aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 38524fb275aSStefano Zampini } 38624fb275aSStefano Zampini 38724fb275aSStefano Zampini /* Computes the new objective given X = Xk, Y = direction 38824fb275aSStefano Zampini W work vector, on output W = X - Y 38924fb275aSStefano Zampini G work vector, on output G = SNESFunction(W) */ 39024fb275aSStefano Zampini static PetscErrorCode SNESNewtonTRObjective(SNES snes, PetscBool has_objective, Vec X, Vec Y, Vec W, Vec G, PetscReal *gnorm, PetscReal *fkp1) 39124fb275aSStefano Zampini { 39224fb275aSStefano Zampini PetscBool changed_y, changed_w; 39324fb275aSStefano Zampini 39424fb275aSStefano Zampini PetscFunctionBegin; 39524fb275aSStefano Zampini /* TODO: we can add a linesearch here */ 39624fb275aSStefano Zampini PetscCall(SNESNewtonTRPreCheck(snes, X, Y, &changed_y)); 39724fb275aSStefano Zampini PetscCall(VecWAXPY(W, -1.0, Y, X)); /* Xkp1 */ 39824fb275aSStefano Zampini PetscCall(SNESNewtonTRPostCheck(snes, X, Y, W, &changed_y, &changed_w)); 39924fb275aSStefano Zampini if (changed_y && !changed_w) PetscCall(VecWAXPY(W, -1.0, Y, X)); 40024fb275aSStefano Zampini 40124fb275aSStefano Zampini PetscCall(SNESComputeFunction(snes, W, G)); /* F(Xkp1) = G */ 40224fb275aSStefano Zampini PetscCall(VecNorm(G, NORM_2, gnorm)); 40324fb275aSStefano Zampini if (has_objective) PetscCall(SNESComputeObjective(snes, W, fkp1)); 40424fb275aSStefano Zampini else *fkp1 = 0.5 * PetscSqr(*gnorm); 40524fb275aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 40624fb275aSStefano Zampini } 40724fb275aSStefano Zampini 40824fb275aSStefano Zampini static PetscErrorCode SNESSetUpQN_NEWTONTR(SNES snes) 40924fb275aSStefano Zampini { 41024fb275aSStefano Zampini SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 41124fb275aSStefano Zampini 41224fb275aSStefano Zampini PetscFunctionBegin; 41324fb275aSStefano Zampini PetscCall(MatDestroy(&tr->qnB)); 41424fb275aSStefano Zampini PetscCall(MatDestroy(&tr->qnB_pre)); 41524fb275aSStefano Zampini if (tr->qn) { 41624fb275aSStefano Zampini PetscInt n, N; 41724fb275aSStefano Zampini const char *optionsprefix; 41824fb275aSStefano Zampini Mat B; 41924fb275aSStefano Zampini 42024fb275aSStefano Zampini PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &B)); 42124fb275aSStefano Zampini PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix)); 42224fb275aSStefano Zampini PetscCall(MatSetOptionsPrefix(B, "snes_tr_qn_")); 42324fb275aSStefano Zampini PetscCall(MatAppendOptionsPrefix(B, optionsprefix)); 42424fb275aSStefano Zampini PetscCall(MatSetType(B, MATLMVMBFGS)); 42524fb275aSStefano Zampini PetscCall(VecGetLocalSize(snes->vec_sol, &n)); 42624fb275aSStefano Zampini PetscCall(VecGetSize(snes->vec_sol, &N)); 42724fb275aSStefano Zampini PetscCall(MatSetSizes(B, n, n, N, N)); 42824fb275aSStefano Zampini PetscCall(MatSetUp(B)); 42924fb275aSStefano Zampini PetscCall(MatSetFromOptions(B)); 43024fb275aSStefano Zampini PetscCall(MatLMVMAllocate(B, snes->vec_sol, snes->vec_func)); 43124fb275aSStefano Zampini tr->qnB = B; 43224fb275aSStefano Zampini if (tr->qn == SNES_TR_QN_DIFFERENT) { 43324fb275aSStefano Zampini PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &B)); 43424fb275aSStefano Zampini PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix)); 43524fb275aSStefano Zampini PetscCall(MatSetOptionsPrefix(B, "snes_tr_qn_pre_")); 43624fb275aSStefano Zampini PetscCall(MatAppendOptionsPrefix(B, optionsprefix)); 43724fb275aSStefano Zampini PetscCall(MatSetType(B, MATLMVMBFGS)); 43824fb275aSStefano Zampini PetscCall(MatSetSizes(B, n, n, N, N)); 43924fb275aSStefano Zampini PetscCall(MatSetUp(B)); 44024fb275aSStefano Zampini PetscCall(MatSetFromOptions(B)); 44124fb275aSStefano Zampini PetscCall(MatLMVMAllocate(B, snes->vec_sol, snes->vec_func)); 44224fb275aSStefano Zampini tr->qnB_pre = B; 44324fb275aSStefano Zampini } else { 44424fb275aSStefano Zampini PetscCall(PetscObjectReference((PetscObject)tr->qnB)); 44524fb275aSStefano Zampini tr->qnB_pre = tr->qnB; 44624fb275aSStefano Zampini } 44724fb275aSStefano Zampini } 4487aa289d8SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 4497aa289d8SStefano Zampini } 4507aa289d8SStefano Zampini 451df60cc22SBarry Smith /* 4524a221d59SStefano Zampini SNESSolve_NEWTONTR - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy 4534a221d59SStefano Zampini (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of 4544a221d59SStefano Zampini nonlinear equations 4554800dd8cSBarry Smith 4564800dd8cSBarry Smith */ 457d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NEWTONTR(SNES snes) 458d71ae5a4SJacob Faibussowitsch { 45904d7464bSBarry Smith SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 46024fb275aSStefano Zampini Vec X, F, Y, G, W, GradF, YU, Yc; 4614a221d59SStefano Zampini PetscInt maxits, lits; 4624b0a5c37SStefano Zampini PetscReal rho, fnorm, gnorm = 0.0, xnorm = 0.0, delta, ynorm; 46324fb275aSStefano Zampini PetscReal deltaM, fk, fkp1, deltaqm = 0.0, gTy = 0.0, yTHy = 0.0; 46424fb275aSStefano Zampini PetscReal auk, tauk, gfnorm, gfnorm_k, ycnorm, gTBg, objmin = 0.0, beta_k = 1.0; 465a0254a93SStefano Zampini PC pc; 46624fb275aSStefano Zampini Mat J, Jp; 46724fb275aSStefano Zampini PetscBool already_done = PETSC_FALSE, on_boundary; 4687aa289d8SStefano Zampini PetscBool clear_converged_test, rho_satisfied, has_objective; 469df8705c3SBarry Smith SNES_TR_KSPConverged_Ctx *ctx; 4705e28bcb6Sprj- void *convctx; 47124fb275aSStefano Zampini 4724a221d59SStefano Zampini PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), (*convdestroy)(void *); 4734a221d59SStefano Zampini PetscErrorCode (*objective)(SNES, Vec, PetscReal *, void *); 4744800dd8cSBarry Smith 4753a40ed3dSBarry Smith PetscFunctionBegin; 4764a221d59SStefano Zampini PetscCall(SNESGetObjective(snes, &objective, NULL)); 4777aa289d8SStefano Zampini has_objective = objective ? PETSC_TRUE : PETSC_FALSE; 478c579b300SPatrick Farrell 479fbe28522SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 480fbe28522SBarry Smith X = snes->vec_sol; /* solution vector */ 48139e2f89bSBarry Smith F = snes->vec_func; /* residual vector */ 4824a221d59SStefano Zampini Y = snes->vec_sol_update; /* update vector */ 4834a221d59SStefano Zampini G = snes->work[0]; /* updated residual */ 4844a221d59SStefano Zampini W = snes->work[1]; /* temporary vector */ 4857aa289d8SStefano Zampini GradF = !has_objective ? snes->work[2] : snes->vec_func; /* grad f = J^T F */ 4864a221d59SStefano Zampini YU = snes->work[3]; /* work vector for dogleg method */ 48724fb275aSStefano Zampini Yc = snes->work[4]; /* Cauchy point */ 4884a221d59SStefano Zampini 4894a221d59SStefano 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); 4904800dd8cSBarry Smith 4919566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 4924c49b128SBarry Smith snes->iter = 0; 4939566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 4944800dd8cSBarry Smith 49524fb275aSStefano Zampini /* setup QN matrices if needed */ 49624fb275aSStefano Zampini PetscCall(SNESSetUpQN_NEWTONTR(snes)); 49724fb275aSStefano Zampini 4984a221d59SStefano Zampini /* Set the linear stopping criteria to use the More' trick if needed */ 4994a221d59SStefano Zampini clear_converged_test = PETSC_FALSE; 500a0254a93SStefano Zampini PetscCall(SNESGetKSP(snes, &snes->ksp)); 501a0254a93SStefano Zampini PetscCall(KSPGetConvergenceTest(snes->ksp, &convtest, &convctx, &convdestroy)); 502fcc61681SStefano Zampini if (convtest != SNESTR_KSPConverged_Private) { 5034a221d59SStefano Zampini clear_converged_test = PETSC_TRUE; 5049566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 505df8705c3SBarry Smith ctx->snes = snes; 506a0254a93SStefano Zampini PetscCall(KSPGetAndClearConvergenceTest(snes->ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 507a0254a93SStefano Zampini PetscCall(KSPSetConvergenceTest(snes->ksp, SNESTR_KSPConverged_Private, ctx, SNESTR_KSPConverged_Destroy)); 5089566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTR_KSPConverged_Private\n")); 509df8705c3SBarry Smith } 510df8705c3SBarry Smith 511e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 5129566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */ 5131aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 514e4ed7901SPeter Brune 5159566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */ 516422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 5179566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */ 5184a221d59SStefano Zampini 5199566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 520fbe28522SBarry Smith snes->norm = fnorm; 5219566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 5224a221d59SStefano Zampini delta = neP->delta0; 5234a221d59SStefano Zampini deltaM = neP->deltaM; 5244800dd8cSBarry Smith neP->delta = delta; 5259566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 526b37302c6SLois Curfman McInnes 52785385478SLisandro Dalcin /* test convergence */ 5284a221d59SStefano Zampini rho_satisfied = PETSC_FALSE; 5292d157150SStefano Zampini PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm)); 5302d157150SStefano Zampini PetscCall(SNESMonitor(snes, 0, fnorm)); 5313ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 5323f149594SLisandro Dalcin 5337aa289d8SStefano Zampini if (has_objective) PetscCall(SNESComputeObjective(snes, X, &fk)); 5344a221d59SStefano Zampini else fk = 0.5 * PetscSqr(fnorm); /* obj(x) = 0.5 * ||F(x)||^2 */ 53599a96b7cSMatthew Knepley 536a0254a93SStefano Zampini /* hook state vector to BFGS preconditioner */ 537a0254a93SStefano Zampini PetscCall(KSPGetPC(snes->ksp, &pc)); 538a0254a93SStefano Zampini PetscCall(PCLMVMSetUpdateVec(pc, X)); 539a0254a93SStefano Zampini 54024fb275aSStefano Zampini if (neP->kmdc) PetscCall(KSPSetComputeEigenvalues(snes->ksp, PETSC_TRUE)); 5416b5873e3SBarry Smith 54224fb275aSStefano Zampini while (snes->iter < maxits) { 54312d0050eSStefano Zampini /* calculating Jacobian and GradF of minimization function only once */ 5444a221d59SStefano Zampini if (!already_done) { 54512d0050eSStefano Zampini /* Call general purpose update function */ 54612d0050eSStefano Zampini PetscTryTypeMethod(snes, update, snes->iter); 54712d0050eSStefano Zampini 5484b0a5c37SStefano Zampini /* apply the nonlinear preconditioner */ 5494b0a5c37SStefano Zampini if (snes->npc && snes->npcside == PC_RIGHT) { 5504b0a5c37SStefano Zampini SNESConvergedReason reason; 5514b0a5c37SStefano Zampini 5524b0a5c37SStefano Zampini PetscCall(SNESSetInitialFunction(snes->npc, F)); 5534b0a5c37SStefano Zampini PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0)); 5544b0a5c37SStefano Zampini PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X)); 5554b0a5c37SStefano Zampini PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0)); 5564b0a5c37SStefano Zampini PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 5574b0a5c37SStefano Zampini if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) { 5584b0a5c37SStefano Zampini snes->reason = SNES_DIVERGED_INNER; 5594b0a5c37SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 5604b0a5c37SStefano Zampini } 5614b0a5c37SStefano Zampini // XXX 5624b0a5c37SStefano Zampini PetscCall(SNESGetNPCFunction(snes, F, &fnorm)); 5634b0a5c37SStefano Zampini } else if (snes->ops->update) { /* if update is present, recompute objective function and function norm */ 56412d0050eSStefano Zampini PetscCall(SNESComputeFunction(snes, X, F)); 56512d0050eSStefano Zampini } 56612d0050eSStefano Zampini 56712d0050eSStefano Zampini /* Jacobian */ 56824fb275aSStefano Zampini J = NULL; 56924fb275aSStefano Zampini Jp = NULL; 57024fb275aSStefano Zampini if (!neP->qnB) { 5714a221d59SStefano Zampini PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 57224fb275aSStefano Zampini J = snes->jacobian; 57324fb275aSStefano Zampini Jp = snes->jacobian_pre; 57424fb275aSStefano Zampini } else { /* QN model */ 57524fb275aSStefano Zampini PetscCall(SNESComputeJacobian_MATLMVM(snes, X, neP->qnB, neP->qnB_pre, NULL)); 57624fb275aSStefano Zampini J = neP->qnB; 57724fb275aSStefano Zampini Jp = neP->qnB_pre; 57824fb275aSStefano Zampini } 5794a221d59SStefano Zampini SNESCheckJacobianDomainerror(snes); 58012d0050eSStefano Zampini 58124fb275aSStefano Zampini /* objective function */ 58224fb275aSStefano Zampini PetscCall(VecNorm(F, NORM_2, &fnorm)); 58324fb275aSStefano Zampini if (has_objective) PetscCall(SNESComputeObjective(snes, X, &fk)); 58424fb275aSStefano Zampini else fk = 0.5 * PetscSqr(fnorm); /* obj(x) = 0.5 * ||F(x)||^2 */ 58524fb275aSStefano Zampini 58612d0050eSStefano Zampini /* GradF */ 5877aa289d8SStefano Zampini if (has_objective) gfnorm = fnorm; 5887aa289d8SStefano Zampini else { 58924fb275aSStefano Zampini PetscCall(MatMultTranspose(J, F, GradF)); /* grad f = J^T F */ 5907aa289d8SStefano Zampini PetscCall(VecNorm(GradF, NORM_2, &gfnorm)); 591fbe28522SBarry Smith } 59224fb275aSStefano Zampini PetscCall(VecNorm(GradF, neP->norm, &gfnorm_k)); 5937aa289d8SStefano Zampini } 5947aa289d8SStefano Zampini already_done = PETSC_TRUE; 5957aa289d8SStefano Zampini 5964b0a5c37SStefano Zampini /* solve trust-region subproblem */ 5974b0a5c37SStefano Zampini 59824fb275aSStefano Zampini /* first compute Cauchy Point */ 59924fb275aSStefano Zampini PetscCall(MatMult(J, GradF, W)); 60024fb275aSStefano Zampini if (has_objective) PetscCall(VecDotRealPart(GradF, W, &gTBg)); 60124fb275aSStefano Zampini else PetscCall(VecDotRealPart(W, W, &gTBg)); /* B = J^t * J */ 60224fb275aSStefano Zampini /* Eqs 4.11 and 4.12 in Nocedal and Wright 2nd Edition (4.7 and 4.8 in 1st Edition) */ 60324fb275aSStefano Zampini auk = delta / gfnorm_k; 60424fb275aSStefano Zampini if (gTBg < 0.0) tauk = 1.0; 60524fb275aSStefano Zampini else tauk = PetscMin(gfnorm * gfnorm * gfnorm_k / (delta * gTBg), 1); 60624fb275aSStefano Zampini auk *= tauk; 60724fb275aSStefano Zampini ycnorm = auk * gfnorm; 60824fb275aSStefano Zampini PetscCall(VecAXPBY(Yc, auk, 0.0, GradF)); 60924fb275aSStefano Zampini 61024fb275aSStefano Zampini on_boundary = PETSC_FALSE; 61124fb275aSStefano Zampini if (tauk != 1.0) { 61224fb275aSStefano Zampini KSPConvergedReason reason; 61324fb275aSStefano Zampini 6144b0a5c37SStefano Zampini /* sufficient decrease (see 6.3.27 in Conn, Gould, Toint "Trust Region Methods") 61524fb275aSStefano Zampini beta_k the largest eigenvalue of the Hessian. Here we use the previous estimated value */ 61624fb275aSStefano Zampini objmin = -neP->kmdc * gnorm * PetscMin(gnorm / beta_k, delta); 617fb01098fSStefano Zampini PetscCall(KSPCGSetObjectiveTarget(snes->ksp, objmin)); 6184b0a5c37SStefano Zampini 61924fb275aSStefano Zampini /* specify radius if looking for Newton step and trust region norm is the l2 norm */ 62024fb275aSStefano Zampini PetscCall(KSPCGSetRadius(snes->ksp, neP->fallback == SNES_TR_FALLBACK_NEWTON && neP->norm == NORM_2 ? delta : 0.0)); 62124fb275aSStefano Zampini PetscCall(KSPSetOperators(snes->ksp, J, Jp)); 6224a221d59SStefano Zampini PetscCall(KSPSolve(snes->ksp, F, Y)); 6234a221d59SStefano Zampini SNESCheckKSPSolve(snes); 6244a221d59SStefano Zampini PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 62524fb275aSStefano Zampini PetscCall(KSPGetConvergedReason(snes->ksp, &reason)); 62624fb275aSStefano Zampini on_boundary = (PetscBool)(reason == KSP_CONVERGED_STEP_LENGTH); 6274b0a5c37SStefano Zampini PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits)); 62824fb275aSStefano Zampini if (neP->kmdc) { /* update estimated Hessian largest eigenvalue */ 62924fb275aSStefano Zampini PetscReal emax, emin; 63024fb275aSStefano Zampini PetscCall(KSPComputeExtremeSingularValues(snes->ksp, &emax, &emin)); 63124fb275aSStefano Zampini if (emax > 0.0) beta_k = emax + 1; 63224fb275aSStefano Zampini } 63324fb275aSStefano Zampini } else { /* Cauchy point is on the boundary, accept it */ 63424fb275aSStefano Zampini on_boundary = PETSC_TRUE; 63524fb275aSStefano Zampini PetscCall(VecCopy(Yc, Y)); 63624fb275aSStefano Zampini PetscCall(PetscInfo(snes, "CP evaluated on boundary. delta: %g, ycnorm: %g, gTBg: %g\n", (double)delta, (double)ycnorm, (double)gTBg)); 63724fb275aSStefano Zampini } 63824fb275aSStefano Zampini PetscCall(VecNorm(Y, neP->norm, &ynorm)); 6394800dd8cSBarry Smith 6404a221d59SStefano Zampini /* decide what to do when the update is outside of trust region */ 641*a02bbafeSStefano Zampini if (tauk != 1.0 && (ynorm > delta || ynorm == 0.0)) { 6425ec2728bSStefano Zampini SNESNewtonTRFallbackType fallback = ynorm > 0.0 ? neP->fallback : SNES_TR_FALLBACK_CAUCHY; 6435ec2728bSStefano Zampini 64424fb275aSStefano Zampini PetscCheck(neP->norm == NORM_2 || fallback != SNES_TR_FALLBACK_DOGLEG, PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "DOGLEG without l2 norm not implemented"); 6455ec2728bSStefano Zampini switch (fallback) { 6464a221d59SStefano Zampini case SNES_TR_FALLBACK_NEWTON: 6474a221d59SStefano Zampini auk = delta / ynorm; 6484a221d59SStefano Zampini PetscCall(VecScale(Y, auk)); 6494b0a5c37SStefano Zampini PetscCall(PetscInfo(snes, "SN evaluated. delta: %g, ynorm: %g\n", (double)delta, (double)ynorm)); 6504a221d59SStefano Zampini break; 6514a221d59SStefano Zampini case SNES_TR_FALLBACK_CAUCHY: 6524a221d59SStefano Zampini case SNES_TR_FALLBACK_DOGLEG: 6535ec2728bSStefano Zampini if (fallback == SNES_TR_FALLBACK_CAUCHY || gTBg <= 0.0) { 65424fb275aSStefano Zampini PetscCall(VecCopy(Yc, Y)); 6554a221d59SStefano Zampini PetscCall(PetscInfo(snes, "CP evaluated. delta: %g, ynorm: %g, ycnorm: %g, gTBg: %g\n", (double)delta, (double)ynorm, (double)ycnorm, (double)gTBg)); 6564a221d59SStefano Zampini } else { /* take linear combination of Cauchy and Newton direction and step */ 65724fb275aSStefano Zampini auk = gfnorm * gfnorm / gTBg; 65824fb275aSStefano Zampini if (gfnorm_k * auk >= delta) { /* first leg: Cauchy point outside of trust region */ 65924fb275aSStefano Zampini PetscCall(VecAXPBY(Y, delta / gfnorm_k, 0.0, GradF)); 66024fb275aSStefano Zampini PetscCall(PetscInfo(snes, "CP evaluated (outside region). delta: %g, ynorm: %g, ycnorm: %g\n", (double)delta, (double)ynorm, (double)ycnorm)); 66124fb275aSStefano Zampini } else { /* second leg */ 6624a221d59SStefano Zampini PetscReal c0, c1, c2, tau = 0.0, tpos, tneg; 6634a221d59SStefano Zampini PetscBool noroots; 664284fb49fSHeeho Park 66524fb275aSStefano Zampini /* Find solutions of (Eq. 4.16 in Nocedal and Wright) 66624fb275aSStefano Zampini ||p_U + lambda * (p_B - p_U)||^2 - delta^2 = 0, 66724fb275aSStefano Zampini where p_U the Cauchy direction, p_B the Newton direction */ 6684a221d59SStefano Zampini PetscCall(VecAXPBY(YU, auk, 0.0, GradF)); 6694a221d59SStefano Zampini PetscCall(VecAXPY(Y, -1.0, YU)); 6704a221d59SStefano Zampini PetscCall(VecNorm(Y, NORM_2, &c0)); 6714a221d59SStefano Zampini PetscCall(VecDotRealPart(YU, Y, &c1)); 6724a221d59SStefano Zampini c0 = PetscSqr(c0); 6734a221d59SStefano Zampini c2 = PetscSqr(ycnorm) - PetscSqr(delta); 67424fb275aSStefano Zampini PetscQuadraticRoots(c0, 2 * c1, c2, &tneg, &tpos); 6754a221d59SStefano Zampini 67624fb275aSStefano Zampini /* In principle the DL strategy as a unique solution in [0,1] 67724fb275aSStefano Zampini here we check that for some reason we numerically failed 67824fb275aSStefano Zampini to compute it. In that case, we use the Cauchy point */ 6794a221d59SStefano Zampini noroots = PetscIsInfOrNanReal(tneg); 68024fb275aSStefano Zampini if (!noroots) { 68124fb275aSStefano Zampini if (tpos > 1) { 68224fb275aSStefano Zampini if (tneg >= 0 && tneg <= 1) { 68324fb275aSStefano Zampini tau = tneg; 68424fb275aSStefano Zampini } else noroots = PETSC_TRUE; 68524fb275aSStefano Zampini } else if (tpos >= 0) { 68624fb275aSStefano Zampini tau = tpos; 68724fb275aSStefano Zampini } else noroots = PETSC_TRUE; 68824fb275aSStefano Zampini } 6894a221d59SStefano Zampini if (noroots) { /* No roots, select Cauchy point */ 69024fb275aSStefano Zampini PetscCall(VecCopy(Yc, Y)); 6914a221d59SStefano Zampini } else { 69224fb275aSStefano Zampini PetscCall(VecAXPBY(Y, 1.0, tau, YU)); 6934a221d59SStefano Zampini } 69424fb275aSStefano 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)); 6954a221d59SStefano Zampini } 6964a221d59SStefano Zampini } 6974a221d59SStefano Zampini break; 6984a221d59SStefano Zampini default: 6994a221d59SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "Unknown fallback mode"); 700454a90a3SBarry Smith break; 70152392280SLois Curfman McInnes } 7024800dd8cSBarry Smith } 7034a221d59SStefano Zampini 7047aa289d8SStefano Zampini /* compute the quadratic model difference */ 70524fb275aSStefano Zampini PetscCall(SNESNewtonTRQuadraticDelta(snes, J, has_objective, Y, GradF, W, &yTHy, &gTy, &deltaqm)); 7064a221d59SStefano Zampini 7074a221d59SStefano Zampini /* Compute new objective function */ 70824fb275aSStefano Zampini PetscCall(SNESNewtonTRObjective(snes, has_objective, X, Y, W, G, &gnorm, &fkp1)); 70924fb275aSStefano Zampini if (PetscIsInfOrNanReal(fkp1)) rho = neP->eta1; 71024fb275aSStefano Zampini else { 7114a221d59SStefano Zampini if (deltaqm > 0.0) rho = (fk - fkp1) / deltaqm; /* actual improvement over predicted improvement */ 71224fb275aSStefano Zampini else rho = neP->eta1; /* no reduction in quadratic model, step must be rejected */ 71324fb275aSStefano Zampini } 7144a221d59SStefano Zampini 71524fb275aSStefano Zampini PetscCall(VecNorm(Y, neP->norm, &ynorm)); 71624fb275aSStefano 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)); 71724fb275aSStefano Zampini 71824fb275aSStefano Zampini /* update the size of the trust region */ 7194a221d59SStefano Zampini if (rho < neP->eta2) delta *= neP->t1; /* shrink the region */ 72024fb275aSStefano Zampini else if (rho > neP->eta3 && on_boundary) delta *= neP->t2; /* expand the region */ 7214a221d59SStefano Zampini delta = PetscMin(delta, deltaM); /* but not greater than deltaM */ 7224a221d59SStefano Zampini 72324fb275aSStefano Zampini /* log 2-norm of update for moniroting routines */ 72424fb275aSStefano Zampini PetscCall(VecNorm(Y, NORM_2, &ynorm)); 72524fb275aSStefano Zampini 72624fb275aSStefano Zampini /* decide on new step */ 7274a221d59SStefano Zampini neP->delta = delta; 72824fb275aSStefano Zampini if (rho > neP->eta1) { 7294a221d59SStefano Zampini rho_satisfied = PETSC_TRUE; 7304a221d59SStefano Zampini } else { 7314a221d59SStefano Zampini rho_satisfied = PETSC_FALSE; 7324a221d59SStefano Zampini PetscCall(PetscInfo(snes, "Trying again in smaller region\n")); 7334a221d59SStefano Zampini /* check to see if progress is hopeless */ 7344a221d59SStefano Zampini PetscCall(SNESTR_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP)); 7352d157150SStefano Zampini if (!snes->reason) PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm)); 7364b0a5c37SStefano Zampini if (snes->reason == SNES_CONVERGED_SNORM_RELATIVE) snes->reason = SNES_DIVERGED_TR_DELTA; 7374a221d59SStefano Zampini snes->numFailures++; 7384a221d59SStefano Zampini /* We're not progressing, so return with the current iterate */ 7394a221d59SStefano Zampini if (snes->reason) break; 7404a221d59SStefano Zampini } 7414a221d59SStefano Zampini if (rho_satisfied) { 7424a221d59SStefano Zampini /* Update function values */ 7434a221d59SStefano Zampini already_done = PETSC_FALSE; 7444800dd8cSBarry Smith fnorm = gnorm; 7454a221d59SStefano Zampini fk = fkp1; 7464a221d59SStefano Zampini 7474a221d59SStefano Zampini /* New residual and linearization point */ 7489566063dSJacob Faibussowitsch PetscCall(VecCopy(G, F)); 7499566063dSJacob Faibussowitsch PetscCall(VecCopy(W, X)); 7504a221d59SStefano Zampini 75185385478SLisandro Dalcin /* Monitor convergence */ 7529566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 7534a221d59SStefano Zampini snes->iter++; 754fbe28522SBarry Smith snes->norm = fnorm; 755c1e67a49SFande Kong snes->xnorm = xnorm; 756c1e67a49SFande Kong snes->ynorm = ynorm; 7579566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 7589566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 7594a221d59SStefano Zampini 76085385478SLisandro Dalcin /* Test for convergence, xnorm = || X || */ 7614a221d59SStefano Zampini PetscCall(VecNorm(X, NORM_2, &xnorm)); 7622d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm)); 7632d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 7644a221d59SStefano Zampini if (snes->reason) break; 7654a221d59SStefano Zampini } 76638442cffSBarry Smith } 767284fb49fSHeeho Park 7684a221d59SStefano Zampini if (clear_converged_test) { 769a0254a93SStefano Zampini PetscCall(KSPGetAndClearConvergenceTest(snes->ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 7709566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 771a0254a93SStefano Zampini PetscCall(KSPSetConvergenceTest(snes->ksp, convtest, convctx, convdestroy)); 7725e28bcb6Sprj- } 7733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7744800dd8cSBarry Smith } 775284fb49fSHeeho Park 776d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes) 777d71ae5a4SJacob Faibussowitsch { 7783a40ed3dSBarry Smith PetscFunctionBegin; 77924fb275aSStefano Zampini PetscCall(SNESSetWorkVecs(snes, 5)); 7809566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 7813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7824800dd8cSBarry Smith } 7836b8b9a38SLisandro Dalcin 78424fb275aSStefano Zampini static PetscErrorCode SNESReset_NEWTONTR(SNES snes) 78524fb275aSStefano Zampini { 78624fb275aSStefano Zampini SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 78724fb275aSStefano Zampini 78824fb275aSStefano Zampini PetscFunctionBegin; 78924fb275aSStefano Zampini PetscCall(MatDestroy(&tr->qnB)); 79024fb275aSStefano Zampini PetscCall(MatDestroy(&tr->qnB_pre)); 79124fb275aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 79224fb275aSStefano Zampini } 79324fb275aSStefano Zampini 794d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes) 795d71ae5a4SJacob Faibussowitsch { 7963a40ed3dSBarry Smith PetscFunctionBegin; 79724fb275aSStefano Zampini PetscCall(SNESReset_NEWTONTR(snes)); 7989566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 7993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8004800dd8cSBarry Smith } 8014800dd8cSBarry Smith 802d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NEWTONTR(SNES snes, PetscOptionItems *PetscOptionsObject) 803d71ae5a4SJacob Faibussowitsch { 80404d7464bSBarry Smith SNES_NEWTONTR *ctx = (SNES_NEWTONTR *)snes->data; 80524fb275aSStefano Zampini SNESNewtonTRQNType qn; 80624fb275aSStefano Zampini SNESNewtonTRFallbackType fallback; 80724fb275aSStefano Zampini NormType norm; 80824fb275aSStefano Zampini PetscReal deltatol; 80924fb275aSStefano Zampini PetscBool flg; 8104800dd8cSBarry Smith 8113a40ed3dSBarry Smith PetscFunctionBegin; 812d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations"); 8134a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_eta1", "eta1", "None", ctx->eta1, &ctx->eta1, NULL)); 8144a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_eta2", "eta2", "None", ctx->eta2, &ctx->eta2, NULL)); 8154a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_eta3", "eta3", "None", ctx->eta3, &ctx->eta3, NULL)); 8164a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_t1", "t1", "None", ctx->t1, &ctx->t1, NULL)); 8174a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_t2", "t2", "None", ctx->t2, &ctx->t2, NULL)); 8184a221d59SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_deltaM", "deltaM", "None", ctx->deltaM, &ctx->deltaM, NULL)); 8199566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_tr_delta0", "delta0", "None", ctx->delta0, &ctx->delta0, NULL)); 8204b0a5c37SStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_kmdc", "sufficient decrease parameter", "None", ctx->kmdc, &ctx->kmdc, NULL)); 82124fb275aSStefano Zampini 82224fb275aSStefano Zampini deltatol = snes->deltatol; 82324fb275aSStefano Zampini PetscCall(PetscOptionsReal("-snes_tr_tol", "Trust region tolerance", "SNESSetTrustRegionTolerance", deltatol, &deltatol, &flg)); 82424fb275aSStefano Zampini if (flg) PetscCall(SNESSetTrustRegionTolerance(snes, deltatol)); 82524fb275aSStefano Zampini 82624fb275aSStefano Zampini fallback = ctx->fallback; 82724fb275aSStefano 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)); 82824fb275aSStefano Zampini if (flg) PetscCall(SNESNewtonTRSetFallbackType(snes, fallback)); 82924fb275aSStefano Zampini 83024fb275aSStefano Zampini qn = ctx->qn; 83124fb275aSStefano Zampini PetscCall(PetscOptionsEnum("-snes_tr_qn", "Use Quasi-Newton approximations for the model", "SNESNewtonTRSetQNType", SNESNewtonTRQNTypes, (PetscEnum)qn, (PetscEnum *)&qn, &flg)); 83224fb275aSStefano Zampini if (flg) PetscCall(SNESNewtonTRSetQNType(snes, qn)); 83324fb275aSStefano Zampini 83424fb275aSStefano Zampini norm = ctx->norm; 83524fb275aSStefano Zampini PetscCall(PetscOptionsEnum("-snes_tr_norm_type", "Type of norm for trust region bounds", "SNESNewtonTRSetNormType", NormTypes, (PetscEnum)norm, (PetscEnum *)&norm, &flg)); 83624fb275aSStefano Zampini if (flg) PetscCall(SNESNewtonTRSetNormType(snes, norm)); 83724fb275aSStefano Zampini 838d0609cedSBarry Smith PetscOptionsHeadEnd(); 8393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8404800dd8cSBarry Smith } 8414800dd8cSBarry Smith 842d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NEWTONTR(SNES snes, PetscViewer viewer) 843d71ae5a4SJacob Faibussowitsch { 84404d7464bSBarry Smith SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 845ace3abfcSBarry Smith PetscBool iascii; 846a935fc98SLois Curfman McInnes 8473a40ed3dSBarry Smith PetscFunctionBegin; 8489566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 84932077d6dSBarry Smith if (iascii) { 8504a221d59SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " Trust region tolerance %g\n", (double)snes->deltatol)); 8514a221d59SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " eta1=%g, eta2=%g, eta3=%g\n", (double)tr->eta1, (double)tr->eta2, (double)tr->eta3)); 8524a221d59SStefano 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)); 8534b0a5c37SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " kmdc=%g\n", (double)tr->kmdc)); 8544a221d59SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " fallback=%s\n", SNESNewtonTRFallbackTypes[tr->fallback])); 85524fb275aSStefano Zampini if (tr->qn) PetscCall(PetscViewerASCIIPrintf(viewer, " qn=%s\n", SNESNewtonTRQNTypes[tr->qn])); 85624fb275aSStefano Zampini if (tr->norm != NORM_2) PetscCall(PetscViewerASCIIPrintf(viewer, " norm=%s\n", NormTypes[tr->norm])); 85719bcc07fSBarry Smith } 8583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 859a935fc98SLois Curfman McInnes } 860f6dfbefdSBarry Smith 861ebe3b25bSBarry Smith /*MC 8621d27aa22SBarry Smith SNESNEWTONTR - Newton based nonlinear solver that uses trust-region dogleg method with Cauchy direction {cite}`nocedal2006numerical` 863f6dfbefdSBarry Smith 864f6dfbefdSBarry Smith Options Database Keys: 8654a221d59SStefano Zampini + -snes_tr_tol <tol> - trust region tolerance 86624fb275aSStefano Zampini . -snes_tr_eta1 <eta1> - trust region parameter eta1 <= eta2, rho > eta1 breaks out of the inner iteration (default: eta1=0.001) 86724fb275aSStefano Zampini . -snes_tr_eta2 <eta2> - trust region parameter, rho <= eta2 shrinks the trust region (default: eta2=0.25) 8684a221d59SStefano Zampini . -snes_tr_eta3 <eta3> - trust region parameter eta3 > eta2, rho >= eta3 expands the trust region (default: eta3=0.75) 8694a221d59SStefano Zampini . -snes_tr_t1 <t1> - trust region parameter, shrinking factor of trust region (default: 0.25) 8704a221d59SStefano Zampini . -snes_tr_t2 <t2> - trust region parameter, expanding factor of trust region (default: 2.0) 8714a221d59SStefano Zampini . -snes_tr_deltaM <deltaM> - trust region parameter, max size of trust region (default: MAX_REAL) 8724a221d59SStefano Zampini . -snes_tr_delta0 <delta0> - trust region parameter, initial size of trust region (default: 0.2) 8734a221d59SStefano 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. 874acba1f63SStefano Zampini 875acba1f63SStefano Zampini Level: beginner 876b3113221SBarry Smith 877420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESSetTrustRegionTolerance()`, 8784a221d59SStefano Zampini `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`, 87924fb275aSStefano Zampini `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRSetFallbackType()`, `SNESNewtonTRSetQNType()` 880ebe3b25bSBarry Smith M*/ 881d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes) 882d71ae5a4SJacob Faibussowitsch { 88304d7464bSBarry Smith SNES_NEWTONTR *neP; 8844800dd8cSBarry Smith 8853a40ed3dSBarry Smith PetscFunctionBegin; 88604d7464bSBarry Smith snes->ops->setup = SNESSetUp_NEWTONTR; 88704d7464bSBarry Smith snes->ops->solve = SNESSolve_NEWTONTR; 88824fb275aSStefano Zampini snes->ops->reset = SNESReset_NEWTONTR; 88904d7464bSBarry Smith snes->ops->destroy = SNESDestroy_NEWTONTR; 89004d7464bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR; 89104d7464bSBarry Smith snes->ops->view = SNESView_NEWTONTR; 892fbe28522SBarry Smith 8931a58d6d3SStefano Zampini snes->stol = 0.0; 89442f4f86dSBarry Smith snes->usesksp = PETSC_TRUE; 8954b0a5c37SStefano Zampini snes->npcside = PC_RIGHT; 8964b0a5c37SStefano Zampini snes->usesnpc = PETSC_TRUE; 89742f4f86dSBarry Smith 8984fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 8994fc747eaSLawrence Mitchell 9004dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&neP)); 901fbe28522SBarry Smith snes->data = (void *)neP; 902fbe28522SBarry Smith neP->delta = 0.0; 903fbe28522SBarry Smith neP->delta0 = 0.2; 9044a221d59SStefano Zampini neP->eta1 = 0.001; 9054a221d59SStefano Zampini neP->eta2 = 0.25; 9064a221d59SStefano Zampini neP->eta3 = 0.75; 9074a221d59SStefano Zampini neP->t1 = 0.25; 9084a221d59SStefano Zampini neP->t2 = 2.0; 9094a221d59SStefano Zampini neP->deltaM = 1.e10; 91024fb275aSStefano Zampini neP->norm = NORM_2; 9114a221d59SStefano Zampini neP->fallback = SNES_TR_FALLBACK_NEWTON; 9124b0a5c37SStefano Zampini neP->kmdc = 0.0; /* by default do not use sufficient decrease */ 9133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9144800dd8cSBarry Smith } 915