10a844d1aSPeter Brune #include <../src/snes/impls/ncg/snesncgimpl.h> /*I "petscsnes.h" I*/ 29e5d0892SLisandro Dalcin const char *const SNESNCGTypes[] = {"FR", "PRP", "HS", "DY", "CD", "SNESNCGType", "SNES_NCG_", NULL}; 3fef7b6d8SPeter Brune 4d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESReset_NCG(SNES snes) 5d71ae5a4SJacob Faibussowitsch { 6fef7b6d8SPeter Brune PetscFunctionBegin; 73ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8fef7b6d8SPeter Brune } 9fef7b6d8SPeter Brune 10d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_NCG(SNES snes) 11d71ae5a4SJacob Faibussowitsch { 12fef7b6d8SPeter Brune PetscFunctionBegin; 132e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNCGSetType_C", NULL)); 149566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16fef7b6d8SPeter Brune } 17fef7b6d8SPeter Brune 18d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_NCG(SNES snes) 19d71ae5a4SJacob Faibussowitsch { 20fef7b6d8SPeter Brune PetscFunctionBegin; 219566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes, 2)); 2208401ef6SPierre Jolivet PetscCheck(snes->npcside != PC_RIGHT, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning"); 236c67d002SPeter Brune if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED; 243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 25fef7b6d8SPeter Brune } 26fef7b6d8SPeter Brune 27d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch) 28d71ae5a4SJacob Faibussowitsch { 29ea630c6eSPeter Brune PetscScalar alpha, ptAp; 30bbd5d0b3SPeter Brune Vec X, Y, F, W; 31bbd5d0b3SPeter Brune SNES snes; 32bbd5d0b3SPeter Brune PetscReal *fnorm, *xnorm, *ynorm; 33ea630c6eSPeter Brune 34ea630c6eSPeter Brune PetscFunctionBegin; 359566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 36bbd5d0b3SPeter Brune X = linesearch->vec_sol; 37bbd5d0b3SPeter Brune W = linesearch->vec_sol_new; 38bbd5d0b3SPeter Brune F = linesearch->vec_func; 39bbd5d0b3SPeter Brune Y = linesearch->vec_update; 40bbd5d0b3SPeter Brune fnorm = &linesearch->fnorm; 41bbd5d0b3SPeter Brune xnorm = &linesearch->xnorm; 42bbd5d0b3SPeter Brune ynorm = &linesearch->ynorm; 43bbd5d0b3SPeter Brune 4448a46eb9SPierre Jolivet if (!snes->jacobian) PetscCall(SNESSetUpMatrices(snes)); 459105210eSPeter Brune 46ea630c6eSPeter Brune /* 47d2140ca3SPeter Brune The exact step size for unpreconditioned linear CG is just: 48ea630c6eSPeter Brune alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy) 49ea630c6eSPeter Brune */ 509566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 519566063dSJacob Faibussowitsch PetscCall(VecDot(F, F, &alpha)); 529566063dSJacob Faibussowitsch PetscCall(MatMult(snes->jacobian, Y, W)); 539566063dSJacob Faibussowitsch PetscCall(VecDot(Y, W, &ptAp)); 54ea630c6eSPeter Brune alpha = alpha / ptAp; 559566063dSJacob Faibussowitsch PetscCall(VecAXPY(X, -alpha, Y)); 569566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 57bbd5d0b3SPeter Brune 589566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, fnorm)); 599566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, xnorm)); 609566063dSJacob Faibussowitsch PetscCall(VecNorm(Y, NORM_2, ynorm)); 613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 62ea630c6eSPeter Brune } 63ea630c6eSPeter Brune 64b5badacbSBarry Smith /*MC 65f6dfbefdSBarry Smith SNESLINESEARCHNCGLINEAR - Special line search only for the nonlinear CG solver `SNESNCG` 66b5badacbSBarry Smith 67b5badacbSBarry Smith This line search uses the length "as if" the problem is linear (that is what is computed by the linear CG method) using the Jacobian of the function. 68b5badacbSBarry Smith alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy) where r (f) is the current residual (function value), p (y) is the current search direction. 69b5badacbSBarry Smith 70420bcc1bSBarry Smith Level: advanced 71420bcc1bSBarry Smith 72f6dfbefdSBarry Smith Notes: 73f6dfbefdSBarry Smith This requires a Jacobian-vector product but does not require the solution of a linear system with the Jacobian 74b5badacbSBarry Smith 75b5badacbSBarry Smith This is a "odd-ball" line search, we don't know if it is in the literature or used in practice by anyone. 76b5badacbSBarry Smith 77420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNCG` `SNESLineSearchCreate()`, `SNESLineSearchSetType()` 78b5badacbSBarry Smith M*/ 79b5badacbSBarry Smith 80d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch) 81d71ae5a4SJacob Faibussowitsch { 82bbd5d0b3SPeter Brune PetscFunctionBegin; 83f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_NCGLinear; 840298fd71SBarry Smith linesearch->ops->destroy = NULL; 850298fd71SBarry Smith linesearch->ops->setfromoptions = NULL; 860298fd71SBarry Smith linesearch->ops->reset = NULL; 870298fd71SBarry Smith linesearch->ops->view = NULL; 880298fd71SBarry Smith linesearch->ops->setup = NULL; 893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 90bbd5d0b3SPeter Brune } 91bbd5d0b3SPeter Brune 92d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NCG(SNES snes, PetscOptionItems *PetscOptionsObject) 93d71ae5a4SJacob Faibussowitsch { 94b5badacbSBarry Smith SNES_NCG *ncg = (SNES_NCG *)snes->data; 95b5badacbSBarry Smith PetscBool debug = PETSC_FALSE; 96b5badacbSBarry Smith SNESNCGType ncgtype = ncg->type; 97b5badacbSBarry Smith SNESLineSearch linesearch; 98b5badacbSBarry Smith 99b5badacbSBarry Smith PetscFunctionBegin; 100d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES NCG options"); 1019566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ncg_monitor", "Monitor the beta values used in the NCG iterations", "SNES", ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL)); 102ad540459SPierre Jolivet if (debug) ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)); 1039566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-snes_ncg_type", "NCG Beta type used", "SNESNCGSetType", SNESNCGTypes, (PetscEnum)ncg->type, (PetscEnum *)&ncgtype, NULL)); 1049566063dSJacob Faibussowitsch PetscCall(SNESNCGSetType(snes, ncgtype)); 105d0609cedSBarry Smith PetscOptionsHeadEnd(); 106b5badacbSBarry Smith if (!snes->linesearch) { 1079566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 108ec786807SJed Brown if (!((PetscObject)linesearch)->type_name) { 109b5badacbSBarry Smith if (!snes->npc) { 1109566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHCP)); 111b5badacbSBarry Smith } else { 1129566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2)); 113b5badacbSBarry Smith } 114b5badacbSBarry Smith } 115ec786807SJed Brown } 1163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 117b5badacbSBarry Smith } 118b5badacbSBarry Smith 119d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer) 120d71ae5a4SJacob Faibussowitsch { 121b5badacbSBarry Smith SNES_NCG *ncg = (SNES_NCG *)snes->data; 122b5badacbSBarry Smith PetscBool iascii; 123b5badacbSBarry Smith 124b5badacbSBarry Smith PetscFunctionBegin; 1259566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 12648a46eb9SPierre Jolivet if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, " type: %s\n", SNESNCGTypes[ncg->type])); 1273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 128b5badacbSBarry Smith } 129b5badacbSBarry Smith 1300a844d1aSPeter Brune /*@ 131f6dfbefdSBarry Smith SNESNCGSetType - Sets the conjugate update type for nonlinear CG `SNESNCG`. 1320a844d1aSPeter Brune 133c3339decSBarry Smith Logically Collective 1340a844d1aSPeter Brune 1350a844d1aSPeter Brune Input Parameters: 1360a844d1aSPeter Brune + snes - the iterative context 137ceaaa498SBarry Smith - btype - update type, see `SNESNCGType` 1380a844d1aSPeter Brune 1393c7db156SBarry Smith Options Database Key: 140b5badacbSBarry Smith . -snes_ncg_type <prp,fr,hs,dy,cd> - strategy for selecting algorithm for computing beta 1410a844d1aSPeter Brune 1420a844d1aSPeter Brune Level: intermediate 1430a844d1aSPeter Brune 1440a844d1aSPeter Brune Notes: 145f6dfbefdSBarry Smith `SNES_NCG_PRP` is the default, and the only one that tolerates generalized search directions. 146b5badacbSBarry Smith 147b5badacbSBarry Smith It is not clear what "generalized search directions" means, does it mean use with a nonlinear preconditioner, 148f6dfbefdSBarry Smith that is using -npc_snes_type <type>, `SNESSetNPC()`, or `SNESGetNPC()`? 1490a844d1aSPeter Brune 150420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESNCGType`, `SNES_NCG_FR`, `SNES_NCG_PRP`, `SNES_NCG_HS`, `SNES_NCG_DY`, `SNES_NCG_CD` 1510a844d1aSPeter Brune @*/ 152d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype) 153d71ae5a4SJacob Faibussowitsch { 1540a844d1aSPeter Brune PetscFunctionBegin; 1550a844d1aSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 156cac4c232SBarry Smith PetscTryMethod(snes, "SNESNCGSetType_C", (SNES, SNESNCGType), (snes, btype)); 1573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1580a844d1aSPeter Brune } 1590a844d1aSPeter Brune 160d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype) 161d71ae5a4SJacob Faibussowitsch { 1620a844d1aSPeter Brune SNES_NCG *ncg = (SNES_NCG *)snes->data; 1636e111a19SKarl Rupp 1640a844d1aSPeter Brune PetscFunctionBegin; 1650a844d1aSPeter Brune ncg->type = btype; 1663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1670a844d1aSPeter Brune } 1680a844d1aSPeter Brune 169fef7b6d8SPeter Brune /* 170fef7b6d8SPeter Brune SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method. 171fef7b6d8SPeter Brune 1722fe279fdSBarry Smith Input Parameter: 1732fe279fdSBarry Smith . snes - the `SNES` context 174fef7b6d8SPeter Brune 175fef7b6d8SPeter Brune Output Parameter: 176fef7b6d8SPeter Brune . outits - number of iterations until termination 177fef7b6d8SPeter Brune 178fef7b6d8SPeter Brune Application Interface Routine: SNESSolve() 179fef7b6d8SPeter Brune */ 180d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NCG(SNES snes) 181d71ae5a4SJacob Faibussowitsch { 182dfb256c7SPeter Brune SNES_NCG *ncg = (SNES_NCG *)snes->data; 18332cc618eSPeter Brune Vec X, dX, lX, F, dXold; 184bbd5d0b3SPeter Brune PetscReal fnorm, ynorm, xnorm, beta = 0.0; 18532cc618eSPeter Brune PetscScalar dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold; 186fef7b6d8SPeter Brune PetscInt maxits, i; 187422a814eSBarry Smith SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED; 188f1c6b773SPeter Brune SNESLineSearch linesearch; 189b7281c8aSPeter Brune SNESConvergedReason reason; 19088d374b2SPeter Brune 191fef7b6d8SPeter Brune PetscFunctionBegin; 1920b121fc5SBarry Smith PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 193c579b300SPatrick Farrell 1949566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite)); 195fef7b6d8SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 196fef7b6d8SPeter Brune 197fef7b6d8SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 198fef7b6d8SPeter Brune X = snes->vec_sol; /* X^n */ 19932cc618eSPeter Brune dXold = snes->work[0]; /* The previous iterate of X */ 20088d374b2SPeter Brune dX = snes->work[1]; /* the preconditioned direction */ 201169e2be8SPeter Brune lX = snes->vec_sol_update; /* the conjugate direction */ 202fef7b6d8SPeter Brune F = snes->vec_func; /* residual vector */ 203fef7b6d8SPeter Brune 2049566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 2059e764e56SPeter Brune 2069566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 207fef7b6d8SPeter Brune snes->iter = 0; 208fef7b6d8SPeter Brune snes->norm = 0.; 2099566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 210fef7b6d8SPeter Brune 211bbd5d0b3SPeter Brune /* compute the initial function and preconditioned update dX */ 212a71552e2SPeter Brune 213efd4aadfSBarry Smith if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 2149566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, NULL, dX)); 2159566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 216b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 217b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 2183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 219b7281c8aSPeter Brune } 2209566063dSJacob Faibussowitsch PetscCall(VecCopy(dX, F)); 2219566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); 222a71552e2SPeter Brune } else { 223e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 2249566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 2251aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 226c1c75074SPeter Brune 227fef7b6d8SPeter Brune /* convergence test */ 2289566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); 229422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 2309566063dSJacob Faibussowitsch PetscCall(VecCopy(F, dX)); 231a71552e2SPeter Brune } 232efd4aadfSBarry Smith if (snes->npc) { 233a71552e2SPeter Brune if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 2349566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, F, dX)); 2359566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 236b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 237b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 2383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 239b7281c8aSPeter Brune } 240a71552e2SPeter Brune } 241a71552e2SPeter Brune } 2429566063dSJacob Faibussowitsch PetscCall(VecCopy(dX, lX)); 2439566063dSJacob Faibussowitsch PetscCall(VecDot(dX, dX, &dXdotdX)); 244a71552e2SPeter Brune 2459566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 246fef7b6d8SPeter Brune snes->norm = fnorm; 2479566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2489566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 249fef7b6d8SPeter Brune 250fef7b6d8SPeter Brune /* test convergence */ 2512d157150SStefano Zampini PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm)); 2522d157150SStefano Zampini PetscCall(SNESMonitor(snes, 0, fnorm)); 2533ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 254fef7b6d8SPeter Brune 255fef7b6d8SPeter Brune /* Call general purpose update function */ 256dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 257fef7b6d8SPeter Brune 258fef7b6d8SPeter Brune /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */ 259bbd5d0b3SPeter Brune 26009c08436SPeter Brune for (i = 1; i < maxits + 1; i++) { 26129c94e34SPeter Brune /* some update types require the old update direction or conjugate direction */ 26248a46eb9SPierre Jolivet if (ncg->type != SNES_NCG_FR) PetscCall(VecCopy(dX, dXold)); 2639566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, lX)); 2649566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetReason(linesearch, &lsresult)); 2659566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm)); 266422a814eSBarry Smith if (lsresult) { 267fef7b6d8SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 268fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 2693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 270fef7b6d8SPeter Brune } 271fef7b6d8SPeter Brune } 272e71169deSBarry Smith if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 273fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 2743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 275fef7b6d8SPeter Brune } 276fef7b6d8SPeter Brune /* Monitor convergence */ 2779566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 278169e2be8SPeter Brune snes->iter = i; 279fef7b6d8SPeter Brune snes->norm = fnorm; 280c1e67a49SFande Kong snes->xnorm = xnorm; 281c1e67a49SFande Kong snes->ynorm = ynorm; 2829566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2839566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 284302dbbaeSPeter Brune 285fef7b6d8SPeter Brune /* Test for convergence */ 2862d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm)); 2872d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 2883ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 289302dbbaeSPeter Brune 290302dbbaeSPeter Brune /* Call general purpose update function */ 291dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 292efd4aadfSBarry Smith if (snes->npc) { 293a71552e2SPeter Brune if (snes->functype == SNES_FUNCTION_PRECONDITIONED) { 2949566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, NULL, dX)); 2959566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 296b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 297b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 2983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 299b7281c8aSPeter Brune } 3009566063dSJacob Faibussowitsch PetscCall(VecCopy(dX, F)); 301a71552e2SPeter Brune } else { 3029566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, F, dX)); 3039566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 304b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 305b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 3063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 307b7281c8aSPeter Brune } 308302dbbaeSPeter Brune } 309a3685310SPeter Brune } else { 3109566063dSJacob Faibussowitsch PetscCall(VecCopy(F, dX)); 311302dbbaeSPeter Brune } 312302dbbaeSPeter Brune 31329c94e34SPeter Brune /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/ 3140a844d1aSPeter Brune switch (ncg->type) { 3150a844d1aSPeter Brune case SNES_NCG_FR: /* Fletcher-Reeves */ 31632cc618eSPeter Brune dXolddotdXold = dXdotdX; 3179566063dSJacob Faibussowitsch PetscCall(VecDot(dX, dX, &dXdotdX)); 31832cc618eSPeter Brune beta = PetscRealPart(dXdotdX / dXolddotdXold); 319dfb256c7SPeter Brune break; 3200a844d1aSPeter Brune case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */ 32132cc618eSPeter Brune dXolddotdXold = dXdotdX; 3229566063dSJacob Faibussowitsch PetscCall(VecDotBegin(dX, dX, &dXdotdXold)); 3239566063dSJacob Faibussowitsch PetscCall(VecDotBegin(dXold, dX, &dXdotdXold)); 3249566063dSJacob Faibussowitsch PetscCall(VecDotEnd(dX, dX, &dXdotdX)); 3259566063dSJacob Faibussowitsch PetscCall(VecDotEnd(dXold, dX, &dXdotdXold)); 32632cc618eSPeter Brune beta = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold)); 327dfb256c7SPeter Brune if (beta < 0.0) beta = 0.0; /* restart */ 328dfb256c7SPeter Brune break; 3290a844d1aSPeter Brune case SNES_NCG_HS: /* Hestenes-Stiefel */ 3309566063dSJacob Faibussowitsch PetscCall(VecDotBegin(dX, dX, &dXdotdX)); 3319566063dSJacob Faibussowitsch PetscCall(VecDotBegin(dX, dXold, &dXdotdXold)); 3329566063dSJacob Faibussowitsch PetscCall(VecDotBegin(lX, dX, &lXdotdX)); 3339566063dSJacob Faibussowitsch PetscCall(VecDotBegin(lX, dXold, &lXdotdXold)); 3349566063dSJacob Faibussowitsch PetscCall(VecDotEnd(dX, dX, &dXdotdX)); 3359566063dSJacob Faibussowitsch PetscCall(VecDotEnd(dX, dXold, &dXdotdXold)); 3369566063dSJacob Faibussowitsch PetscCall(VecDotEnd(lX, dX, &lXdotdX)); 3379566063dSJacob Faibussowitsch PetscCall(VecDotEnd(lX, dXold, &lXdotdXold)); 33832cc618eSPeter Brune beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold)); 339dfb256c7SPeter Brune break; 3400a844d1aSPeter Brune case SNES_NCG_DY: /* Dai-Yuan */ 3419566063dSJacob Faibussowitsch PetscCall(VecDotBegin(dX, dX, &dXdotdX)); 3429566063dSJacob Faibussowitsch PetscCall(VecDotBegin(lX, dX, &lXdotdX)); 3439566063dSJacob Faibussowitsch PetscCall(VecDotBegin(lX, dXold, &lXdotdXold)); 3449566063dSJacob Faibussowitsch PetscCall(VecDotEnd(dX, dX, &dXdotdX)); 3459566063dSJacob Faibussowitsch PetscCall(VecDotEnd(lX, dX, &lXdotdX)); 3469566063dSJacob Faibussowitsch PetscCall(VecDotEnd(lX, dXold, &lXdotdXold)); 3472f613bf5SBarry Smith beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX)); 348dfb256c7SPeter Brune break; 3490a844d1aSPeter Brune case SNES_NCG_CD: /* Conjugate Descent */ 3509566063dSJacob Faibussowitsch PetscCall(VecDotBegin(dX, dX, &dXdotdX)); 3519566063dSJacob Faibussowitsch PetscCall(VecDotBegin(lX, dXold, &lXdotdXold)); 3529566063dSJacob Faibussowitsch PetscCall(VecDotEnd(dX, dX, &dXdotdX)); 3539566063dSJacob Faibussowitsch PetscCall(VecDotEnd(lX, dXold, &lXdotdXold)); 3542f613bf5SBarry Smith beta = PetscRealPart(dXdotdX / lXdotdXold); 355dfb256c7SPeter Brune break; 356dfb256c7SPeter Brune } 35748a46eb9SPierre Jolivet if (ncg->monitor) PetscCall(PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta)); 3589566063dSJacob Faibussowitsch PetscCall(VecAYPX(lX, beta, dX)); 359fef7b6d8SPeter Brune } 3603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 361fef7b6d8SPeter Brune } 362fef7b6d8SPeter Brune 363fef7b6d8SPeter Brune /*MC 364*1d27aa22SBarry Smith SNESNCG - Nonlinear Conjugate-Gradient method for the solution of nonlinear systems {cite}`bruneknepleysmithtu15`. 365fef7b6d8SPeter Brune 366fef7b6d8SPeter Brune Level: beginner 367fef7b6d8SPeter Brune 368f6dfbefdSBarry Smith Options Database Keys: 3694f02bc6aSBarry Smith + -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is prp. 37041a8cb50SPeter Brune . -snes_linesearch_type <cp,l2,basic> - Line search type. 371*1d27aa22SBarry Smith - -snes_ncg_monitor - Print the beta values nonlinear Conjugate-Gradient used in the iteration, . 372fef7b6d8SPeter Brune 37395452b02SPatrick Sanan Notes: 374*1d27aa22SBarry Smith This solves the nonlinear system of equations $ F(x) = 0 $ using the nonlinear generalization of the conjugate 375fef7b6d8SPeter Brune gradient method. This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise 376*1d27aa22SBarry Smith chooses the initial search direction as $ F(x) $ for the initial guess $x$. 377fef7b6d8SPeter Brune 3782d547940SBarry Smith Only supports left non-linear preconditioning. 3792d547940SBarry Smith 380f6dfbefdSBarry Smith Default line search is `SNESLINESEARCHCP`, unless a nonlinear preconditioner is used with -npc_snes_type <type>, `SNESSetNPC()`, or `SNESGetNPC()` then 381f6dfbefdSBarry Smith `SNESLINESEARCHL2` is used. Also supports the special purpose line search `SNESLINESEARCHNCGLINEAR` 382b5badacbSBarry Smith 383420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESNGMRES`, `SNESQN`, `SNESLINESEARCHNCGLINEAR`, `SNESNCGSetType()`, `SNESLineSearchSetType()` 384fef7b6d8SPeter Brune M*/ 385d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes) 386d71ae5a4SJacob Faibussowitsch { 387ea630c6eSPeter Brune SNES_NCG *neP; 388fef7b6d8SPeter Brune 389fef7b6d8SPeter Brune PetscFunctionBegin; 390fef7b6d8SPeter Brune snes->ops->destroy = SNESDestroy_NCG; 391fef7b6d8SPeter Brune snes->ops->setup = SNESSetUp_NCG; 392fef7b6d8SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NCG; 393fef7b6d8SPeter Brune snes->ops->view = SNESView_NCG; 394fef7b6d8SPeter Brune snes->ops->solve = SNESSolve_NCG; 395fef7b6d8SPeter Brune snes->ops->reset = SNESReset_NCG; 396fef7b6d8SPeter Brune 397fef7b6d8SPeter Brune snes->usesksp = PETSC_FALSE; 398efd4aadfSBarry Smith snes->usesnpc = PETSC_TRUE; 399efd4aadfSBarry Smith snes->npcside = PC_LEFT; 400fef7b6d8SPeter Brune 4014fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 4024fc747eaSLawrence Mitchell 40388976e71SPeter Brune if (!snes->tolerancesset) { 4040e444f03SPeter Brune snes->max_funcs = 30000; 4050e444f03SPeter Brune snes->max_its = 10000; 406c522fa08SPeter Brune snes->stol = 1e-20; 40788976e71SPeter Brune } 4080e444f03SPeter Brune 4094dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&neP)); 410fef7b6d8SPeter Brune snes->data = (void *)neP; 4110298fd71SBarry Smith neP->monitor = NULL; 4120a844d1aSPeter Brune neP->type = SNES_NCG_PRP; 4139566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNCGSetType_C", SNESNCGSetType_NCG)); 4143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 415fef7b6d8SPeter Brune } 416