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");} 48028b6e76SBarry Smith PetscFunctionReturn(0); 49028b6e76SBarry Smith } 50028b6e76SBarry Smith 51028b6e76SBarry Smith /* 5204d7464bSBarry 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; 655fd66863SKarl Rupp 66028b6e76SBarry Smith PetscFunctionBegin; 6742f4f86dSBarry Smith ierr = PetscOptionsHead("SNES Richardson options");CHKERRQ(ierr); 68028b6e76SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 699e764e56SPeter Brune if (!snes->linesearch) { 707601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 711a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 729e764e56SPeter Brune } 73028b6e76SBarry Smith PetscFunctionReturn(0); 74028b6e76SBarry Smith } 75028b6e76SBarry Smith 76028b6e76SBarry Smith /* 77d5c3842bSBarry Smith SNESView_NRichardson - Prints info from the SNESRichardson data structure. 78028b6e76SBarry Smith 79028b6e76SBarry Smith Input Parameters: 80028b6e76SBarry Smith + SNES - the SNES context 81028b6e76SBarry Smith - viewer - visualization context 82028b6e76SBarry Smith 83028b6e76SBarry Smith Application Interface Routine: SNESView() 84028b6e76SBarry Smith */ 85028b6e76SBarry Smith #undef __FUNCT__ 86d5c3842bSBarry Smith #define __FUNCT__ "SNESView_NRichardson" 87d5c3842bSBarry Smith static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer) 88028b6e76SBarry Smith { 89028b6e76SBarry Smith PetscBool iascii; 90028b6e76SBarry Smith PetscErrorCode ierr; 91028b6e76SBarry Smith 92028b6e76SBarry Smith PetscFunctionBegin; 93251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 94028b6e76SBarry Smith if (iascii) { 95028b6e76SBarry Smith } 96028b6e76SBarry Smith PetscFunctionReturn(0); 97028b6e76SBarry Smith } 98028b6e76SBarry Smith 99028b6e76SBarry Smith /* 100d5c3842bSBarry Smith SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method. 101028b6e76SBarry Smith 102028b6e76SBarry Smith Input Parameters: 103028b6e76SBarry Smith . snes - the SNES context 104028b6e76SBarry Smith 105028b6e76SBarry Smith Output Parameter: 106028b6e76SBarry Smith . outits - number of iterations until termination 107028b6e76SBarry Smith 108028b6e76SBarry Smith Application Interface Routine: SNESSolve() 109028b6e76SBarry Smith */ 110028b6e76SBarry Smith #undef __FUNCT__ 111d5c3842bSBarry Smith #define __FUNCT__ "SNESSolve_NRichardson" 112d5c3842bSBarry Smith PetscErrorCode SNESSolve_NRichardson(SNES snes) 113028b6e76SBarry Smith { 114bf7f4e0aSPeter Brune Vec X, Y, F; 115e7058c64SPeter Brune PetscReal xnorm, fnorm, ynorm; 116028b6e76SBarry Smith PetscInt maxits, i; 117028b6e76SBarry Smith PetscErrorCode ierr; 118bf7f4e0aSPeter Brune PetscBool lsSuccess; 119028b6e76SBarry Smith SNESConvergedReason reason; 120028b6e76SBarry Smith 121028b6e76SBarry Smith PetscFunctionBegin; 122028b6e76SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 123028b6e76SBarry Smith 124028b6e76SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 125028b6e76SBarry Smith X = snes->vec_sol; /* X^n */ 126028b6e76SBarry Smith Y = snes->vec_sol_update; /* \tilde X */ 127028b6e76SBarry Smith F = snes->vec_func; /* residual vector */ 128028b6e76SBarry Smith 129ce8c27fbSBarry Smith ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 130028b6e76SBarry Smith snes->iter = 0; 131028b6e76SBarry Smith snes->norm = 0.; 132ce8c27fbSBarry Smith ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 133*b7281c8aSPeter Brune 134*b7281c8aSPeter Brune if (snes->pc && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 135*b7281c8aSPeter Brune ierr = SNESApplyPC(snes,X,NULL,NULL,F);CHKERRQ(ierr); 136*b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 137*b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 138*b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 139*b7281c8aSPeter Brune PetscFunctionReturn(0); 140*b7281c8aSPeter Brune } 141*b7281c8aSPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 142*b7281c8aSPeter Brune } else { 143e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 144028b6e76SBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 145028b6e76SBarry Smith if (snes->domainerror) { 146028b6e76SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 147028b6e76SBarry Smith PetscFunctionReturn(0); 148028b6e76SBarry Smith } 1491aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 150e4ed7901SPeter Brune if (!snes->norm_init_set) { 151*b7281c8aSPeter Brune /* convergence test */ 152*b7281c8aSPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 153189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 154189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 155189a9710SBarry Smith PetscFunctionReturn(0); 156189a9710SBarry Smith } 157e4ed7901SPeter Brune } else { 158e4ed7901SPeter Brune fnorm = snes->norm_init; 159e4ed7901SPeter Brune snes->norm_init_set = PETSC_FALSE; 160e4ed7901SPeter Brune } 161*b7281c8aSPeter Brune } 162*b7281c8aSPeter Brune if (snes->pc && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 163*b7281c8aSPeter Brune ierr = SNESApplyPC(snes,X,F,&fnorm,Y);CHKERRQ(ierr); 164*b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 165*b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 166*b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 167*b7281c8aSPeter Brune PetscFunctionReturn(0); 168*b7281c8aSPeter Brune } 169*b7281c8aSPeter Brune } else { 170*b7281c8aSPeter Brune ierr = VecCopy(F,Y);CHKERRQ(ierr); 171*b7281c8aSPeter Brune } 172e4ed7901SPeter Brune 173ce8c27fbSBarry Smith ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 174028b6e76SBarry Smith snes->norm = fnorm; 175ce8c27fbSBarry Smith ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 176a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 177028b6e76SBarry Smith ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 178028b6e76SBarry Smith 179028b6e76SBarry Smith /* set parameter for default relative tolerance convergence test */ 180028b6e76SBarry Smith snes->ttol = fnorm*snes->rtol; 181028b6e76SBarry Smith /* test convergence */ 182028b6e76SBarry Smith ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 183028b6e76SBarry Smith if (snes->reason) PetscFunctionReturn(0); 184028b6e76SBarry Smith 185028b6e76SBarry Smith /* Call general purpose update function */ 186028b6e76SBarry Smith if (snes->ops->update) { 187028b6e76SBarry Smith ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 188028b6e76SBarry Smith } 189*b7281c8aSPeter Brune 190*b7281c8aSPeter Brune /* set parameter for default relative tolerance convergence test */ 191*b7281c8aSPeter Brune snes->ttol = fnorm*snes->rtol; 192*b7281c8aSPeter Brune /* test convergence */ 193*b7281c8aSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 194*b7281c8aSPeter Brune if (snes->reason) PetscFunctionReturn(0); 195*b7281c8aSPeter Brune 196*b7281c8aSPeter Brune for (i = 1; i < maxits+1; i++) { 197*b7281c8aSPeter Brune lsSuccess = PETSC_TRUE; 198*b7281c8aSPeter Brune 199f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr); 200f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 201f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lsSuccess);CHKERRQ(ierr); 202028b6e76SBarry Smith if (!lsSuccess) { 203028b6e76SBarry Smith if (++snes->numFailures >= snes->maxFailures) { 204028b6e76SBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 205028b6e76SBarry Smith break; 206028b6e76SBarry Smith } 207028b6e76SBarry Smith } 208028b6e76SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 209028b6e76SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 210028b6e76SBarry Smith break; 211028b6e76SBarry Smith } 212028b6e76SBarry Smith if (snes->domainerror) { 213028b6e76SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 214028b6e76SBarry Smith PetscFunctionReturn(0); 215028b6e76SBarry Smith } 216*b7281c8aSPeter Brune 217028b6e76SBarry Smith /* Monitor convergence */ 218ce8c27fbSBarry Smith ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 219*b7281c8aSPeter Brune snes->iter = i; 220028b6e76SBarry Smith snes->norm = fnorm; 221ce8c27fbSBarry Smith ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 222a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 223028b6e76SBarry Smith ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 224028b6e76SBarry Smith /* Test for convergence */ 225e7058c64SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 226028b6e76SBarry Smith if (snes->reason) break; 227*b7281c8aSPeter Brune 228*b7281c8aSPeter Brune /* Call general purpose update function */ 229*b7281c8aSPeter Brune if (snes->ops->update) { 230*b7281c8aSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 231028b6e76SBarry Smith } 232*b7281c8aSPeter Brune 233*b7281c8aSPeter Brune if (snes->pc) { 234*b7281c8aSPeter Brune if (snes->functype == SNES_FUNCTION_PRECONDITIONED) { 235*b7281c8aSPeter Brune ierr = SNESApplyPC(snes,X,NULL,NULL,Y);CHKERRQ(ierr); 236*b7281c8aSPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 237*b7281c8aSPeter Brune ierr = VecCopy(Y,F);CHKERRQ(ierr); 238*b7281c8aSPeter Brune } else { 239*b7281c8aSPeter Brune ierr = SNESApplyPC(snes,X,F,&fnorm,Y);CHKERRQ(ierr); 240*b7281c8aSPeter Brune } 241*b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 242*b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 243*b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 244*b7281c8aSPeter Brune PetscFunctionReturn(0); 245*b7281c8aSPeter Brune } 246*b7281c8aSPeter Brune } else { 247*b7281c8aSPeter Brune ierr = VecCopy(F,Y);CHKERRQ(ierr); 248*b7281c8aSPeter Brune } 249*b7281c8aSPeter Brune } 250*b7281c8aSPeter Brune if (i == maxits+1) { 251028b6e76SBarry Smith ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 252028b6e76SBarry Smith if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 253028b6e76SBarry Smith } 254028b6e76SBarry Smith PetscFunctionReturn(0); 255028b6e76SBarry Smith } 256028b6e76SBarry Smith 257028b6e76SBarry Smith /*MC 258d5c3842bSBarry Smith SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration. 259028b6e76SBarry Smith 260028b6e76SBarry Smith Level: beginner 261028b6e76SBarry Smith 262028b6e76SBarry Smith Options Database: 26372835e02SPeter Brune + -snes_linesearch_type <l2,cp,basic> Line search type. 26472835e02SPeter Brune - -snes_linesearch_damping<1.0> Damping for the line search. 265028b6e76SBarry Smith 266af60355fSPeter Brune Notes: If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda 267af60355fSPeter Brune (F(x^n) - b) where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search. If 268af60355fSPeter Brune an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetPC()) then the inner 269af60355fSPeter Brune solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n} 270af60355fSPeter Brune where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver. 271028b6e76SBarry Smith 27272835e02SPeter Brune The update, especially without inner nonlinear preconditioner, may be ill-scaled. If using the basic 27372835e02SPeter Brune linesearch, one may have to scale the update with -snes_linesearch_damping 27472835e02SPeter Brune 275028b6e76SBarry Smith This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls 276028b6e76SBarry Smith 27704d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESNCG 278028b6e76SBarry Smith M*/ 279028b6e76SBarry Smith #undef __FUNCT__ 280d5c3842bSBarry Smith #define __FUNCT__ "SNESCreate_NRichardson" 2818cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes) 282028b6e76SBarry Smith { 28392c02d66SPeter Brune PetscErrorCode ierr; 284bf7f4e0aSPeter Brune SNES_NRichardson *neP; 2855fd66863SKarl Rupp 286028b6e76SBarry Smith PetscFunctionBegin; 287d5c3842bSBarry Smith snes->ops->destroy = SNESDestroy_NRichardson; 288d5c3842bSBarry Smith snes->ops->setup = SNESSetUp_NRichardson; 289d5c3842bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NRichardson; 290d5c3842bSBarry Smith snes->ops->view = SNESView_NRichardson; 291d5c3842bSBarry Smith snes->ops->solve = SNESSolve_NRichardson; 292d5c3842bSBarry Smith snes->ops->reset = SNESReset_NRichardson; 293028b6e76SBarry Smith 294028b6e76SBarry Smith snes->usesksp = PETSC_FALSE; 29542f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 296028b6e76SBarry Smith 297c6b63b32SPeter Brune snes->functype = SNES_FUNCTION_UNPRECONDITIONED; 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