142f4f86dSBarry Smith #include <../src/snes/impls/richardson/snesrichardsonimpl.h> 2028b6e76SBarry Smith 366976f2fSJacob Faibussowitsch static PetscErrorCode SNESDestroy_NRichardson(SNES snes) 4d71ae5a4SJacob Faibussowitsch { 5028b6e76SBarry Smith PetscFunctionBegin; 69566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 73ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8028b6e76SBarry Smith } 9028b6e76SBarry Smith 1066976f2fSJacob Faibussowitsch static PetscErrorCode SNESSetUp_NRichardson(SNES snes) 11d71ae5a4SJacob Faibussowitsch { 12028b6e76SBarry Smith PetscFunctionBegin; 1308401ef6SPierre Jolivet PetscCheck(snes->npcside != PC_RIGHT, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "NRichardson only supports left preconditioning"); 146c67d002SPeter Brune if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED; 153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16028b6e76SBarry Smith } 17028b6e76SBarry Smith 18ce78bad3SBarry Smith static PetscErrorCode SNESSetFromOptions_NRichardson(SNES snes, PetscOptionItems PetscOptionsObject) 19d71ae5a4SJacob Faibussowitsch { 20028b6e76SBarry Smith PetscFunctionBegin; 21d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES Richardson options"); 22d0609cedSBarry Smith PetscOptionsHeadEnd(); 233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 24028b6e76SBarry Smith } 25028b6e76SBarry Smith 26d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer) 27d71ae5a4SJacob Faibussowitsch { 28028b6e76SBarry Smith PetscBool iascii; 29028b6e76SBarry Smith 30028b6e76SBarry Smith PetscFunctionBegin; 319566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 329371c9d4SSatish Balay if (iascii) { } 333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 34028b6e76SBarry Smith } 35028b6e76SBarry Smith 3666976f2fSJacob Faibussowitsch static PetscErrorCode SNESSolve_NRichardson(SNES snes) 37d71ae5a4SJacob Faibussowitsch { 38bf7f4e0aSPeter Brune Vec X, Y, F; 39e7058c64SPeter Brune PetscReal xnorm, fnorm, ynorm; 40028b6e76SBarry Smith PetscInt maxits, i; 41422a814eSBarry Smith SNESLineSearchReason lsresult; 42028b6e76SBarry Smith SNESConvergedReason reason; 43028b6e76SBarry Smith 44028b6e76SBarry Smith PetscFunctionBegin; 450b121fc5SBarry 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); 46c579b300SPatrick Farrell 47028b6e76SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 48028b6e76SBarry Smith 49028b6e76SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 50028b6e76SBarry Smith X = snes->vec_sol; /* X^n */ 51028b6e76SBarry Smith Y = snes->vec_sol_update; /* \tilde X */ 52028b6e76SBarry Smith F = snes->vec_func; /* residual vector */ 53028b6e76SBarry Smith 549566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 55028b6e76SBarry Smith snes->iter = 0; 56028b6e76SBarry Smith snes->norm = 0.; 579566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 58b7281c8aSPeter Brune 59efd4aadfSBarry Smith if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 609566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, NULL, F)); 619566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 62b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 63b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 65b7281c8aSPeter Brune } 669566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); 67b7281c8aSPeter Brune } else { 68e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 699566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 701aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 71c1c75074SPeter Brune 729566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); 73422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 74b7281c8aSPeter Brune } 75efd4aadfSBarry Smith if (snes->npc && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 769566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, F, Y)); 779566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 78b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 79b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 81b7281c8aSPeter Brune } 82b7281c8aSPeter Brune } else { 839566063dSJacob Faibussowitsch PetscCall(VecCopy(F, Y)); 84b7281c8aSPeter Brune } 85e4ed7901SPeter Brune 869566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 87028b6e76SBarry Smith snes->norm = fnorm; 889566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 899566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 90028b6e76SBarry Smith 91028b6e76SBarry Smith /* test convergence */ 922d157150SStefano Zampini PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm)); 932d157150SStefano Zampini PetscCall(SNESMonitor(snes, 0, fnorm)); 943ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 95028b6e76SBarry Smith 96028b6e76SBarry Smith /* Call general purpose update function */ 97dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 98b7281c8aSPeter Brune 99b7281c8aSPeter Brune for (i = 1; i < maxits + 1; i++) { 1009566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y)); 1019566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetReason(snes->linesearch, &lsresult)); 1029566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm)); 103422a814eSBarry Smith if (lsresult) { 104028b6e76SBarry Smith if (++snes->numFailures >= snes->maxFailures) { 105028b6e76SBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 106028b6e76SBarry Smith break; 107028b6e76SBarry Smith } 108028b6e76SBarry Smith } 109e71169deSBarry Smith if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 110028b6e76SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 111028b6e76SBarry Smith break; 112028b6e76SBarry Smith } 113b7281c8aSPeter Brune 114028b6e76SBarry Smith /* Monitor convergence */ 1159566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 116b7281c8aSPeter Brune snes->iter = i; 117028b6e76SBarry Smith snes->norm = fnorm; 118c1e67a49SFande Kong snes->xnorm = xnorm; 119c1e67a49SFande Kong snes->ynorm = ynorm; 1209566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 1219566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 122028b6e76SBarry Smith /* Test for convergence */ 1232d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm)); 1242d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 125028b6e76SBarry Smith if (snes->reason) break; 126b7281c8aSPeter Brune 127b7281c8aSPeter Brune /* Call general purpose update function */ 128dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 129b7281c8aSPeter Brune 130efd4aadfSBarry Smith if (snes->npc) { 131b7281c8aSPeter Brune if (snes->functype == SNES_FUNCTION_PRECONDITIONED) { 1329566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, NULL, Y)); 1339566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); 1349566063dSJacob Faibussowitsch PetscCall(VecCopy(Y, F)); 135b7281c8aSPeter Brune } else { 1369566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, F, Y)); 137b7281c8aSPeter Brune } 1389566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 139b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 140b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 1413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 142b7281c8aSPeter Brune } 143b7281c8aSPeter Brune } else { 1449566063dSJacob Faibussowitsch PetscCall(VecCopy(F, Y)); 145b7281c8aSPeter Brune } 146b7281c8aSPeter Brune } 147b7281c8aSPeter Brune if (i == maxits + 1) { 14863a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits)); 149028b6e76SBarry Smith if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 150028b6e76SBarry Smith } 1513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 152028b6e76SBarry Smith } 153028b6e76SBarry Smith 154028b6e76SBarry Smith /*MC 155d5c3842bSBarry Smith SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration. 156028b6e76SBarry Smith 157f6dfbefdSBarry Smith Options Database Keys: 15867b8a455SSatish Balay + -snes_linesearch_type <l2,cp,basic> - Line search type. 15967b8a455SSatish Balay - -snes_linesearch_damping <1.0> - Damping for the line search. 160028b6e76SBarry Smith 161f6dfbefdSBarry Smith Level: beginner 162f6dfbefdSBarry Smith 16395452b02SPatrick Sanan Notes: 1640b4b7b1cSBarry Smith If no inner nonlinear preconditioner is provided then solves $F(x) - b = 0 $ using $x^{n+1} = x^{n} - \lambda 1650b4b7b1cSBarry Smith (F(x^n) - b) $ where $ \lambda$ is obtained with either `SNESLineSearchSetDamping()`, `-snes_damping` or a line search. If 1660b4b7b1cSBarry Smith an inner nonlinear preconditioner is provided (either with `-npc_snes_typ`e or `SNESSetNPC()`) then the inner 1670b4b7b1cSBarry Smith solver is called on the initial solution $x^n$ and the nonlinear Richardson uses $ x^{n+1} = x^{n} + \lambda d^{n}$ 1680b4b7b1cSBarry Smith where $d^{n} = \hat{x}^{n} - x^{n} $ where $\hat{x}^{n} $ is the solution returned from the inner solver. 169028b6e76SBarry Smith 17072835e02SPeter Brune The update, especially without inner nonlinear preconditioner, may be ill-scaled. If using the basic 1710b4b7b1cSBarry Smith linesearch, one may have to scale the update with `-snes_linesearch_damping` 17272835e02SPeter Brune 173*e87b5d96SPierre Jolivet This uses no derivative information provided with `SNESSetJacobian()` thus it will be much slower than Newton's method obtained with `-snes_type ls` 174028b6e76SBarry Smith 1752d547940SBarry Smith Only supports left non-linear preconditioning. 1762d547940SBarry Smith 177420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESNGMRES`, `SNESQN`, `SNESNCG`, 178420bcc1bSBarry Smith `SNESLineSearchSetDamping()` 179028b6e76SBarry Smith M*/ 180d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes) 181d71ae5a4SJacob Faibussowitsch { 182bf7f4e0aSPeter Brune SNES_NRichardson *neP; 183d8d34be6SBarry Smith SNESLineSearch linesearch; 1845fd66863SKarl Rupp 185028b6e76SBarry Smith PetscFunctionBegin; 186d5c3842bSBarry Smith snes->ops->destroy = SNESDestroy_NRichardson; 187d5c3842bSBarry Smith snes->ops->setup = SNESSetUp_NRichardson; 188d5c3842bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NRichardson; 189d5c3842bSBarry Smith snes->ops->view = SNESView_NRichardson; 190d5c3842bSBarry Smith snes->ops->solve = SNESSolve_NRichardson; 191028b6e76SBarry Smith 192028b6e76SBarry Smith snes->usesksp = PETSC_FALSE; 193efd4aadfSBarry Smith snes->usesnpc = PETSC_TRUE; 194028b6e76SBarry Smith 195efd4aadfSBarry Smith snes->npcside = PC_LEFT; 196c6b63b32SPeter Brune 1979566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 19848a46eb9SPierre Jolivet if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2)); 199d8d34be6SBarry Smith 2004fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 2014fc747eaSLawrence Mitchell 20277e5a1f9SBarry Smith PetscCall(SNESParametersInitialize(snes)); 20377e5a1f9SBarry Smith PetscObjectParameterSetDefault(snes, max_funcs, 30000); 20477e5a1f9SBarry Smith PetscObjectParameterSetDefault(snes, max_its, 10000); 20577e5a1f9SBarry Smith PetscObjectParameterSetDefault(snes, stol, 1e-20); 20677e5a1f9SBarry Smith 2074dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&neP)); 208bf7f4e0aSPeter Brune snes->data = (void *)neP; 2093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 210028b6e76SBarry Smith } 211