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 PetscFunctionBegin; 9028b6e76SBarry Smith PetscFunctionReturn(0); 10028b6e76SBarry Smith } 11028b6e76SBarry Smith 12028b6e76SBarry Smith /* 13d5c3842bSBarry Smith SNESDestroy_NRichardson - Destroys the private SNES_NRichardson context that was created with SNESCreate_NRichardson(). 14028b6e76SBarry Smith 15028b6e76SBarry Smith Input Parameter: 16028b6e76SBarry Smith . snes - the SNES context 17028b6e76SBarry Smith 18028b6e76SBarry Smith Application Interface Routine: SNESDestroy() 19028b6e76SBarry Smith */ 20028b6e76SBarry Smith #undef __FUNCT__ 21d5c3842bSBarry Smith #define __FUNCT__ "SNESDestroy_NRichardson" 22d5c3842bSBarry Smith PetscErrorCode SNESDestroy_NRichardson(SNES snes) 23028b6e76SBarry Smith { 24028b6e76SBarry Smith PetscErrorCode ierr; 25028b6e76SBarry Smith 26028b6e76SBarry Smith PetscFunctionBegin; 27d5c3842bSBarry Smith ierr = SNESReset_NRichardson(snes);CHKERRQ(ierr); 28028b6e76SBarry Smith ierr = PetscFree(snes->data);CHKERRQ(ierr); 29028b6e76SBarry Smith PetscFunctionReturn(0); 30028b6e76SBarry Smith } 31028b6e76SBarry Smith 32028b6e76SBarry Smith /* 33d5c3842bSBarry Smith SNESSetUp_NRichardson - Sets up the internal data structures for the later use 34d5c3842bSBarry Smith of the SNESNRICHARDSON nonlinear solver. 35028b6e76SBarry Smith 36028b6e76SBarry Smith Input Parameters: 37028b6e76SBarry Smith + snes - the SNES context 38028b6e76SBarry Smith - x - the solution vector 39028b6e76SBarry Smith 40028b6e76SBarry Smith Application Interface Routine: SNESSetUp() 41028b6e76SBarry Smith */ 42028b6e76SBarry Smith #undef __FUNCT__ 43d5c3842bSBarry Smith #define __FUNCT__ "SNESSetUp_NRichardson" 44d5c3842bSBarry Smith PetscErrorCode SNESSetUp_NRichardson(SNES snes) 45028b6e76SBarry Smith { 46028b6e76SBarry Smith PetscFunctionBegin; 47c6b63b32SPeter Brune if (snes->pcside == PC_RIGHT) {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"NRichardson only supports left preconditioning");} 486c67d002SPeter Brune if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED; 49028b6e76SBarry Smith PetscFunctionReturn(0); 50028b6e76SBarry Smith } 51028b6e76SBarry Smith 52028b6e76SBarry Smith /* 5304d7464bSBarry Smith SNESSetFromOptions_NRichardson - Sets various parameters for the SNESNEWTONLS method. 54028b6e76SBarry Smith 55028b6e76SBarry Smith Input Parameter: 56028b6e76SBarry Smith . snes - the SNES context 57028b6e76SBarry Smith 58028b6e76SBarry Smith Application Interface Routine: SNESSetFromOptions() 59028b6e76SBarry Smith */ 60028b6e76SBarry Smith #undef __FUNCT__ 61d5c3842bSBarry Smith #define __FUNCT__ "SNESSetFromOptions_NRichardson" 62d5c3842bSBarry Smith static PetscErrorCode SNESSetFromOptions_NRichardson(SNES snes) 63028b6e76SBarry Smith { 64028b6e76SBarry Smith PetscErrorCode ierr; 65f1c6b773SPeter Brune SNESLineSearch linesearch; 665fd66863SKarl Rupp 67028b6e76SBarry Smith PetscFunctionBegin; 6842f4f86dSBarry Smith ierr = PetscOptionsHead("SNES Richardson options");CHKERRQ(ierr); 69028b6e76SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 709e764e56SPeter Brune if (!snes->linesearch) { 717601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 721a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 739e764e56SPeter Brune } 74028b6e76SBarry Smith PetscFunctionReturn(0); 75028b6e76SBarry Smith } 76028b6e76SBarry Smith 77028b6e76SBarry Smith /* 78d5c3842bSBarry Smith SNESView_NRichardson - Prints info from the SNESRichardson data structure. 79028b6e76SBarry Smith 80028b6e76SBarry Smith Input Parameters: 81028b6e76SBarry Smith + SNES - the SNES context 82028b6e76SBarry Smith - viewer - visualization context 83028b6e76SBarry Smith 84028b6e76SBarry Smith Application Interface Routine: SNESView() 85028b6e76SBarry Smith */ 86028b6e76SBarry Smith #undef __FUNCT__ 87d5c3842bSBarry Smith #define __FUNCT__ "SNESView_NRichardson" 88d5c3842bSBarry Smith static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer) 89028b6e76SBarry Smith { 90028b6e76SBarry Smith PetscBool iascii; 91028b6e76SBarry Smith PetscErrorCode ierr; 92028b6e76SBarry Smith 93028b6e76SBarry Smith PetscFunctionBegin; 94251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 95028b6e76SBarry Smith if (iascii) { 96028b6e76SBarry Smith } 97028b6e76SBarry Smith PetscFunctionReturn(0); 98028b6e76SBarry Smith } 99028b6e76SBarry Smith 100028b6e76SBarry Smith /* 101d5c3842bSBarry Smith SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method. 102028b6e76SBarry Smith 103028b6e76SBarry Smith Input Parameters: 104028b6e76SBarry Smith . snes - the SNES context 105028b6e76SBarry Smith 106028b6e76SBarry Smith Output Parameter: 107028b6e76SBarry Smith . outits - number of iterations until termination 108028b6e76SBarry Smith 109028b6e76SBarry Smith Application Interface Routine: SNESSolve() 110028b6e76SBarry Smith */ 111028b6e76SBarry Smith #undef __FUNCT__ 112d5c3842bSBarry Smith #define __FUNCT__ "SNESSolve_NRichardson" 113d5c3842bSBarry Smith PetscErrorCode SNESSolve_NRichardson(SNES snes) 114028b6e76SBarry Smith { 115bf7f4e0aSPeter Brune Vec X, Y, F; 116e7058c64SPeter Brune PetscReal xnorm, fnorm, ynorm; 117028b6e76SBarry Smith PetscInt maxits, i; 118028b6e76SBarry Smith PetscErrorCode ierr; 119bf7f4e0aSPeter Brune PetscBool lsSuccess; 120028b6e76SBarry Smith SNESConvergedReason reason; 121028b6e76SBarry Smith 122028b6e76SBarry Smith PetscFunctionBegin; 123028b6e76SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 124028b6e76SBarry Smith 125028b6e76SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 126028b6e76SBarry Smith X = snes->vec_sol; /* X^n */ 127028b6e76SBarry Smith Y = snes->vec_sol_update; /* \tilde X */ 128028b6e76SBarry Smith F = snes->vec_func; /* residual vector */ 129028b6e76SBarry Smith 130*e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 131028b6e76SBarry Smith snes->iter = 0; 132028b6e76SBarry Smith snes->norm = 0.; 133*e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 134b7281c8aSPeter Brune 135b7281c8aSPeter Brune if (snes->pc && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 136b7281c8aSPeter Brune ierr = SNESApplyPC(snes,X,NULL,NULL,F);CHKERRQ(ierr); 137b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 138b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 139b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 140b7281c8aSPeter Brune PetscFunctionReturn(0); 141b7281c8aSPeter Brune } 142b7281c8aSPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 143b7281c8aSPeter Brune } else { 144e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 145028b6e76SBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 146028b6e76SBarry Smith if (snes->domainerror) { 147028b6e76SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 148028b6e76SBarry Smith PetscFunctionReturn(0); 149028b6e76SBarry Smith } 1501aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 151e4ed7901SPeter Brune if (!snes->norm_init_set) { 152b7281c8aSPeter Brune /* convergence test */ 153b7281c8aSPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 154189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 155189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 156189a9710SBarry Smith PetscFunctionReturn(0); 157189a9710SBarry Smith } 158e4ed7901SPeter Brune } else { 159e4ed7901SPeter Brune fnorm = snes->norm_init; 160e4ed7901SPeter Brune snes->norm_init_set = PETSC_FALSE; 161e4ed7901SPeter Brune } 162b7281c8aSPeter Brune } 163b7281c8aSPeter Brune if (snes->pc && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 164b7281c8aSPeter Brune ierr = SNESApplyPC(snes,X,F,&fnorm,Y);CHKERRQ(ierr); 165b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 166b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 167b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 168b7281c8aSPeter Brune PetscFunctionReturn(0); 169b7281c8aSPeter Brune } 170b7281c8aSPeter Brune } else { 171b7281c8aSPeter Brune ierr = VecCopy(F,Y);CHKERRQ(ierr); 172b7281c8aSPeter Brune } 173e4ed7901SPeter Brune 174*e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 175028b6e76SBarry Smith snes->norm = fnorm; 176*e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 177a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 178028b6e76SBarry Smith ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 179028b6e76SBarry Smith 180028b6e76SBarry Smith /* set parameter for default relative tolerance convergence test */ 181028b6e76SBarry Smith snes->ttol = fnorm*snes->rtol; 182028b6e76SBarry Smith /* test convergence */ 183028b6e76SBarry Smith ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 184028b6e76SBarry Smith if (snes->reason) PetscFunctionReturn(0); 185028b6e76SBarry Smith 186028b6e76SBarry Smith /* Call general purpose update function */ 187028b6e76SBarry Smith if (snes->ops->update) { 188028b6e76SBarry Smith ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 189028b6e76SBarry Smith } 190b7281c8aSPeter Brune 191b7281c8aSPeter Brune /* set parameter for default relative tolerance convergence test */ 192b7281c8aSPeter Brune snes->ttol = fnorm*snes->rtol; 193b7281c8aSPeter Brune /* test convergence */ 194b7281c8aSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 195b7281c8aSPeter Brune if (snes->reason) PetscFunctionReturn(0); 196b7281c8aSPeter Brune 197b7281c8aSPeter Brune for (i = 1; i < maxits+1; i++) { 198b7281c8aSPeter Brune lsSuccess = PETSC_TRUE; 199b7281c8aSPeter Brune 200f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 201f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 202f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lsSuccess);CHKERRQ(ierr); 203028b6e76SBarry Smith if (!lsSuccess) { 204028b6e76SBarry Smith if (++snes->numFailures >= snes->maxFailures) { 205028b6e76SBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 206028b6e76SBarry Smith break; 207028b6e76SBarry Smith } 208028b6e76SBarry Smith } 209028b6e76SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 210028b6e76SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 211028b6e76SBarry Smith break; 212028b6e76SBarry Smith } 213028b6e76SBarry Smith if (snes->domainerror) { 214028b6e76SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 215028b6e76SBarry Smith PetscFunctionReturn(0); 216028b6e76SBarry Smith } 217b7281c8aSPeter Brune 218028b6e76SBarry Smith /* Monitor convergence */ 219*e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 220b7281c8aSPeter Brune snes->iter = i; 221028b6e76SBarry Smith snes->norm = fnorm; 222*e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 223a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 224028b6e76SBarry Smith ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 225028b6e76SBarry Smith /* Test for convergence */ 226e7058c64SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 227028b6e76SBarry Smith if (snes->reason) break; 228b7281c8aSPeter Brune 229b7281c8aSPeter Brune /* Call general purpose update function */ 230b7281c8aSPeter Brune if (snes->ops->update) { 231b7281c8aSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 232028b6e76SBarry Smith } 233b7281c8aSPeter Brune 234b7281c8aSPeter Brune if (snes->pc) { 235b7281c8aSPeter Brune if (snes->functype == SNES_FUNCTION_PRECONDITIONED) { 236b7281c8aSPeter Brune ierr = SNESApplyPC(snes,X,NULL,NULL,Y);CHKERRQ(ierr); 237b7281c8aSPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 238b7281c8aSPeter Brune ierr = VecCopy(Y,F);CHKERRQ(ierr); 239b7281c8aSPeter Brune } else { 240b7281c8aSPeter Brune ierr = SNESApplyPC(snes,X,F,&fnorm,Y);CHKERRQ(ierr); 241b7281c8aSPeter Brune } 242b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 243b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 244b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 245b7281c8aSPeter Brune PetscFunctionReturn(0); 246b7281c8aSPeter Brune } 247b7281c8aSPeter Brune } else { 248b7281c8aSPeter Brune ierr = VecCopy(F,Y);CHKERRQ(ierr); 249b7281c8aSPeter Brune } 250b7281c8aSPeter Brune } 251b7281c8aSPeter Brune if (i == maxits+1) { 252028b6e76SBarry Smith ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 253028b6e76SBarry Smith if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 254028b6e76SBarry Smith } 255028b6e76SBarry Smith PetscFunctionReturn(0); 256028b6e76SBarry Smith } 257028b6e76SBarry Smith 258028b6e76SBarry Smith /*MC 259d5c3842bSBarry Smith SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration. 260028b6e76SBarry Smith 261028b6e76SBarry Smith Level: beginner 262028b6e76SBarry Smith 263028b6e76SBarry Smith Options Database: 26472835e02SPeter Brune + -snes_linesearch_type <l2,cp,basic> Line search type. 26572835e02SPeter Brune - -snes_linesearch_damping<1.0> Damping for the line search. 266028b6e76SBarry Smith 267af60355fSPeter Brune Notes: If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda 268af60355fSPeter Brune (F(x^n) - b) where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search. If 269af60355fSPeter Brune an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetPC()) then the inner 270af60355fSPeter Brune solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n} 271af60355fSPeter Brune where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver. 272028b6e76SBarry Smith 27372835e02SPeter Brune The update, especially without inner nonlinear preconditioner, may be ill-scaled. If using the basic 27472835e02SPeter Brune linesearch, one may have to scale the update with -snes_linesearch_damping 27572835e02SPeter Brune 276028b6e76SBarry Smith This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls 277028b6e76SBarry Smith 27804d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESNCG 279028b6e76SBarry Smith M*/ 280028b6e76SBarry Smith #undef __FUNCT__ 281d5c3842bSBarry Smith #define __FUNCT__ "SNESCreate_NRichardson" 2828cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes) 283028b6e76SBarry Smith { 28492c02d66SPeter Brune PetscErrorCode ierr; 285bf7f4e0aSPeter Brune SNES_NRichardson *neP; 2865fd66863SKarl Rupp 287028b6e76SBarry Smith PetscFunctionBegin; 288d5c3842bSBarry Smith snes->ops->destroy = SNESDestroy_NRichardson; 289d5c3842bSBarry Smith snes->ops->setup = SNESSetUp_NRichardson; 290d5c3842bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NRichardson; 291d5c3842bSBarry Smith snes->ops->view = SNESView_NRichardson; 292d5c3842bSBarry Smith snes->ops->solve = SNESSolve_NRichardson; 293d5c3842bSBarry Smith snes->ops->reset = SNESReset_NRichardson; 294028b6e76SBarry Smith 295028b6e76SBarry Smith snes->usesksp = PETSC_FALSE; 29642f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 297028b6e76SBarry Smith 298c6b63b32SPeter Brune snes->pcside = PC_LEFT; 299c6b63b32SPeter Brune 300bf7f4e0aSPeter Brune ierr = PetscNewLog(snes, SNES_NRichardson, &neP);CHKERRQ(ierr); 301bf7f4e0aSPeter Brune snes->data = (void*) neP; 302bf7f4e0aSPeter Brune 30388976e71SPeter Brune if (!snes->tolerancesset) { 3040e444f03SPeter Brune snes->max_funcs = 30000; 3050e444f03SPeter Brune snes->max_its = 10000; 306c522fa08SPeter Brune snes->stol = 1e-20; 30788976e71SPeter Brune } 308028b6e76SBarry Smith PetscFunctionReturn(0); 309028b6e76SBarry Smith } 310