142f4f86dSBarry Smith #include <../src/snes/impls/richardson/snesrichardsonimpl.h> 2028b6e76SBarry Smith 3bf7f4e0aSPeter Brune 4028b6e76SBarry Smith #undef __FUNCT__ 5d5c3842bSBarry Smith #define __FUNCT__ "SNESReset_NRichardson" 6d5c3842bSBarry Smith PetscErrorCode SNESReset_NRichardson(SNES snes) 7028b6e76SBarry Smith { 8028b6e76SBarry Smith 9028b6e76SBarry Smith PetscFunctionBegin; 10028b6e76SBarry Smith PetscFunctionReturn(0); 11028b6e76SBarry Smith } 12028b6e76SBarry Smith 13028b6e76SBarry Smith /* 14d5c3842bSBarry Smith SNESDestroy_NRichardson - Destroys the private SNES_NRichardson context that was created with SNESCreate_NRichardson(). 15028b6e76SBarry Smith 16028b6e76SBarry Smith Input Parameter: 17028b6e76SBarry Smith . snes - the SNES context 18028b6e76SBarry Smith 19028b6e76SBarry Smith Application Interface Routine: SNESDestroy() 20028b6e76SBarry Smith */ 21028b6e76SBarry Smith #undef __FUNCT__ 22d5c3842bSBarry Smith #define __FUNCT__ "SNESDestroy_NRichardson" 23d5c3842bSBarry Smith PetscErrorCode SNESDestroy_NRichardson(SNES snes) 24028b6e76SBarry Smith { 25028b6e76SBarry Smith PetscErrorCode ierr; 26028b6e76SBarry Smith 27028b6e76SBarry Smith PetscFunctionBegin; 28d5c3842bSBarry Smith ierr = SNESReset_NRichardson(snes);CHKERRQ(ierr); 29028b6e76SBarry Smith ierr = PetscFree(snes->data);CHKERRQ(ierr); 30028b6e76SBarry Smith PetscFunctionReturn(0); 31028b6e76SBarry Smith } 32028b6e76SBarry Smith 33028b6e76SBarry Smith /* 34d5c3842bSBarry Smith SNESSetUp_NRichardson - Sets up the internal data structures for the later use 35d5c3842bSBarry Smith of the SNESNRICHARDSON nonlinear solver. 36028b6e76SBarry Smith 37028b6e76SBarry Smith Input Parameters: 38028b6e76SBarry Smith + snes - the SNES context 39028b6e76SBarry Smith - x - the solution vector 40028b6e76SBarry Smith 41028b6e76SBarry Smith Application Interface Routine: SNESSetUp() 42028b6e76SBarry Smith */ 43028b6e76SBarry Smith #undef __FUNCT__ 44d5c3842bSBarry Smith #define __FUNCT__ "SNESSetUp_NRichardson" 45d5c3842bSBarry Smith PetscErrorCode SNESSetUp_NRichardson(SNES snes) 46028b6e76SBarry Smith { 47028b6e76SBarry Smith PetscFunctionBegin; 48028b6e76SBarry Smith PetscFunctionReturn(0); 49028b6e76SBarry Smith } 50028b6e76SBarry Smith 51028b6e76SBarry Smith /* 52*04d7464bSBarry Smith SNESSetFromOptions_NRichardson - Sets various parameters for the SNESNEWTONLS method. 53028b6e76SBarry Smith 54028b6e76SBarry Smith Input Parameter: 55028b6e76SBarry Smith . snes - the SNES context 56028b6e76SBarry Smith 57028b6e76SBarry Smith Application Interface Routine: SNESSetFromOptions() 58028b6e76SBarry Smith */ 59028b6e76SBarry Smith #undef __FUNCT__ 60d5c3842bSBarry Smith #define __FUNCT__ "SNESSetFromOptions_NRichardson" 61d5c3842bSBarry Smith static PetscErrorCode SNESSetFromOptions_NRichardson(SNES snes) 62028b6e76SBarry Smith { 63028b6e76SBarry Smith PetscErrorCode ierr; 64f1c6b773SPeter Brune SNESLineSearch linesearch; 65028b6e76SBarry Smith PetscFunctionBegin; 6642f4f86dSBarry Smith ierr = PetscOptionsHead("SNES Richardson options");CHKERRQ(ierr); 67028b6e76SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 689e764e56SPeter Brune if (!snes->linesearch) { 69f1c6b773SPeter Brune ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 701a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 719e764e56SPeter Brune } 72028b6e76SBarry Smith PetscFunctionReturn(0); 73028b6e76SBarry Smith } 74028b6e76SBarry Smith 75028b6e76SBarry Smith /* 76d5c3842bSBarry Smith SNESView_NRichardson - Prints info from the SNESRichardson data structure. 77028b6e76SBarry Smith 78028b6e76SBarry Smith Input Parameters: 79028b6e76SBarry Smith + SNES - the SNES context 80028b6e76SBarry Smith - viewer - visualization context 81028b6e76SBarry Smith 82028b6e76SBarry Smith Application Interface Routine: SNESView() 83028b6e76SBarry Smith */ 84028b6e76SBarry Smith #undef __FUNCT__ 85d5c3842bSBarry Smith #define __FUNCT__ "SNESView_NRichardson" 86d5c3842bSBarry Smith static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer) 87028b6e76SBarry Smith { 88028b6e76SBarry Smith PetscBool iascii; 89028b6e76SBarry Smith PetscErrorCode ierr; 90028b6e76SBarry Smith 91028b6e76SBarry Smith PetscFunctionBegin; 92251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 93028b6e76SBarry Smith if (iascii) { 94028b6e76SBarry Smith } 95028b6e76SBarry Smith PetscFunctionReturn(0); 96028b6e76SBarry Smith } 97028b6e76SBarry Smith 98028b6e76SBarry Smith /* 99d5c3842bSBarry Smith SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method. 100028b6e76SBarry Smith 101028b6e76SBarry Smith Input Parameters: 102028b6e76SBarry Smith . snes - the SNES context 103028b6e76SBarry Smith 104028b6e76SBarry Smith Output Parameter: 105028b6e76SBarry Smith . outits - number of iterations until termination 106028b6e76SBarry Smith 107028b6e76SBarry Smith Application Interface Routine: SNESSolve() 108028b6e76SBarry Smith */ 109028b6e76SBarry Smith #undef __FUNCT__ 110d5c3842bSBarry Smith #define __FUNCT__ "SNESSolve_NRichardson" 111d5c3842bSBarry Smith PetscErrorCode SNESSolve_NRichardson(SNES snes) 112028b6e76SBarry Smith { 113bf7f4e0aSPeter Brune Vec X, Y, F; 114e7058c64SPeter Brune PetscReal xnorm, fnorm, ynorm; 115028b6e76SBarry Smith PetscInt maxits, i; 116028b6e76SBarry Smith PetscErrorCode ierr; 117bf7f4e0aSPeter Brune PetscBool lsSuccess; 118028b6e76SBarry Smith SNESConvergedReason reason; 119028b6e76SBarry Smith 120028b6e76SBarry Smith PetscFunctionBegin; 121028b6e76SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 122028b6e76SBarry Smith 123028b6e76SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 124028b6e76SBarry Smith X = snes->vec_sol; /* X^n */ 125028b6e76SBarry Smith Y = snes->vec_sol_update; /* \tilde X */ 126028b6e76SBarry Smith F = snes->vec_func; /* residual vector */ 127028b6e76SBarry Smith 128028b6e76SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 129028b6e76SBarry Smith snes->iter = 0; 130028b6e76SBarry Smith snes->norm = 0.; 131028b6e76SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 132e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 133028b6e76SBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 134028b6e76SBarry Smith if (snes->domainerror) { 135028b6e76SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 136028b6e76SBarry Smith PetscFunctionReturn(0); 137028b6e76SBarry Smith } 138e4ed7901SPeter Brune } else { 139e4ed7901SPeter Brune snes->vec_func_init_set = PETSC_FALSE; 140e4ed7901SPeter Brune } 141e4ed7901SPeter Brune if (!snes->norm_init_set) { 142028b6e76SBarry Smith ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 143028b6e76SBarry Smith if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 144e4ed7901SPeter Brune } else { 145e4ed7901SPeter Brune fnorm = snes->norm_init; 146e4ed7901SPeter Brune snes->norm_init_set = PETSC_FALSE; 147e4ed7901SPeter Brune } 148e4ed7901SPeter Brune 149028b6e76SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 150028b6e76SBarry Smith snes->norm = fnorm; 151028b6e76SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 152028b6e76SBarry Smith SNESLogConvHistory(snes,fnorm,0); 153028b6e76SBarry Smith ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 154028b6e76SBarry Smith 155028b6e76SBarry Smith /* set parameter for default relative tolerance convergence test */ 156028b6e76SBarry Smith snes->ttol = fnorm*snes->rtol; 157028b6e76SBarry Smith /* test convergence */ 158028b6e76SBarry Smith ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 159028b6e76SBarry Smith if (snes->reason) PetscFunctionReturn(0); 160028b6e76SBarry Smith 161028b6e76SBarry Smith for (i = 0; i < maxits; i++) { 162bf7f4e0aSPeter Brune lsSuccess = PETSC_TRUE; 163028b6e76SBarry Smith /* Call general purpose update function */ 164028b6e76SBarry Smith if (snes->ops->update) { 165028b6e76SBarry Smith ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 166028b6e76SBarry Smith } 167c40d0f55SPeter Brune if (snes->pc && snes->pcside == PC_RIGHT) { 168028b6e76SBarry Smith ierr = VecCopy(X,Y);CHKERRQ(ierr); 169217b9c2eSPeter Brune ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr); 170217b9c2eSPeter Brune ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr); 171c90fad12SPeter Brune ierr = SNESSolve(snes->pc, snes->vec_rhs, Y);CHKERRQ(ierr); 172028b6e76SBarry Smith ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 17334b225ddSBarry Smith if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 174028b6e76SBarry Smith snes->reason = SNES_DIVERGED_INNER; 175028b6e76SBarry Smith PetscFunctionReturn(0); 176028b6e76SBarry Smith } 177e5583369SPeter Brune ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr); 178a3685310SPeter Brune } else { 179a3685310SPeter Brune ierr = VecCopy(F,Y);CHKERRQ(ierr); 180028b6e76SBarry Smith } 181f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 182f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 183f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lsSuccess);CHKERRQ(ierr); 184028b6e76SBarry Smith if (!lsSuccess) { 185028b6e76SBarry Smith if (++snes->numFailures >= snes->maxFailures) { 186028b6e76SBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 187028b6e76SBarry Smith break; 188028b6e76SBarry Smith } 189028b6e76SBarry Smith } 190028b6e76SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 191028b6e76SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 192028b6e76SBarry Smith break; 193028b6e76SBarry Smith } 194028b6e76SBarry Smith if (snes->domainerror) { 195028b6e76SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 196028b6e76SBarry Smith PetscFunctionReturn(0); 197028b6e76SBarry Smith } 198028b6e76SBarry Smith /* Monitor convergence */ 199028b6e76SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 200028b6e76SBarry Smith snes->iter = i+1; 201028b6e76SBarry Smith snes->norm = fnorm; 202028b6e76SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 203028b6e76SBarry Smith SNESLogConvHistory(snes,snes->norm,0); 204028b6e76SBarry Smith ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 205028b6e76SBarry Smith /* Test for convergence */ 206e7058c64SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 207028b6e76SBarry Smith if (snes->reason) break; 208028b6e76SBarry Smith } 209028b6e76SBarry Smith if (i == maxits) { 210028b6e76SBarry Smith ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 211028b6e76SBarry Smith if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 212028b6e76SBarry Smith } 213028b6e76SBarry Smith PetscFunctionReturn(0); 214028b6e76SBarry Smith } 215028b6e76SBarry Smith 216028b6e76SBarry Smith /*MC 217d5c3842bSBarry Smith SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration. 218028b6e76SBarry Smith 219028b6e76SBarry Smith Level: beginner 220028b6e76SBarry Smith 221028b6e76SBarry Smith Options Database: 22272835e02SPeter Brune + -snes_linesearch_type <l2,cp,basic> Line search type. 22372835e02SPeter Brune - -snes_linesearch_damping<1.0> Damping for the line search. 224028b6e76SBarry Smith 225af60355fSPeter Brune Notes: If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda 226af60355fSPeter Brune (F(x^n) - b) where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search. If 227af60355fSPeter Brune an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetPC()) then the inner 228af60355fSPeter Brune solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n} 229af60355fSPeter Brune where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver. 230028b6e76SBarry Smith 23172835e02SPeter Brune The update, especially without inner nonlinear preconditioner, may be ill-scaled. If using the basic 23272835e02SPeter Brune linesearch, one may have to scale the update with -snes_linesearch_damping 23372835e02SPeter Brune 234028b6e76SBarry Smith This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls 235028b6e76SBarry Smith 236*04d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESNCG 237028b6e76SBarry Smith M*/ 238028b6e76SBarry Smith EXTERN_C_BEGIN 239028b6e76SBarry Smith #undef __FUNCT__ 240d5c3842bSBarry Smith #define __FUNCT__ "SNESCreate_NRichardson" 241d5c3842bSBarry Smith PetscErrorCode SNESCreate_NRichardson(SNES snes) 242028b6e76SBarry Smith { 24392c02d66SPeter Brune PetscErrorCode ierr; 244bf7f4e0aSPeter Brune SNES_NRichardson *neP; 245028b6e76SBarry Smith PetscFunctionBegin; 246d5c3842bSBarry Smith snes->ops->destroy = SNESDestroy_NRichardson; 247d5c3842bSBarry Smith snes->ops->setup = SNESSetUp_NRichardson; 248d5c3842bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NRichardson; 249d5c3842bSBarry Smith snes->ops->view = SNESView_NRichardson; 250d5c3842bSBarry Smith snes->ops->solve = SNESSolve_NRichardson; 251d5c3842bSBarry Smith snes->ops->reset = SNESReset_NRichardson; 252028b6e76SBarry Smith 253028b6e76SBarry Smith snes->usesksp = PETSC_FALSE; 25442f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 255028b6e76SBarry Smith 256bf7f4e0aSPeter Brune ierr = PetscNewLog(snes, SNES_NRichardson, &neP);CHKERRQ(ierr); 257bf7f4e0aSPeter Brune snes->data = (void*) neP; 258bf7f4e0aSPeter Brune 25988976e71SPeter Brune if (!snes->tolerancesset) { 2600e444f03SPeter Brune snes->max_funcs = 30000; 2610e444f03SPeter Brune snes->max_its = 10000; 262c522fa08SPeter Brune snes->stol = 1e-20; 26388976e71SPeter Brune } 2640e444f03SPeter Brune 265028b6e76SBarry Smith PetscFunctionReturn(0); 266028b6e76SBarry Smith } 267028b6e76SBarry Smith EXTERN_C_END 268