1028b6e76SBarry Smith 242f4f86dSBarry Smith #include <../src/snes/impls/richardson/snesrichardsonimpl.h> 3028b6e76SBarry Smith 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 PetscErrorCode ierr; 48028b6e76SBarry Smith 49028b6e76SBarry Smith PetscFunctionBegin; 50f9a8e2e6SPeter Brune ierr = SNESDefaultGetWork(snes,2);CHKERRQ(ierr); 51028b6e76SBarry Smith PetscFunctionReturn(0); 52028b6e76SBarry Smith } 53028b6e76SBarry Smith 54028b6e76SBarry Smith /* 55d5c3842bSBarry Smith SNESSetFromOptions_NRichardson - Sets various parameters for the SNESLS method. 56028b6e76SBarry Smith 57028b6e76SBarry Smith Input Parameter: 58028b6e76SBarry Smith . snes - the SNES context 59028b6e76SBarry Smith 60028b6e76SBarry Smith Application Interface Routine: SNESSetFromOptions() 61028b6e76SBarry Smith */ 62028b6e76SBarry Smith #undef __FUNCT__ 63d5c3842bSBarry Smith #define __FUNCT__ "SNESSetFromOptions_NRichardson" 64d5c3842bSBarry Smith static PetscErrorCode SNESSetFromOptions_NRichardson(SNES snes) 65028b6e76SBarry Smith { 66028b6e76SBarry Smith PetscErrorCode ierr; 67028b6e76SBarry Smith PetscFunctionBegin; 6842f4f86dSBarry Smith ierr = PetscOptionsHead("SNES Richardson options");CHKERRQ(ierr); 69028b6e76SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 70028b6e76SBarry Smith PetscFunctionReturn(0); 71028b6e76SBarry Smith } 72028b6e76SBarry Smith 73028b6e76SBarry Smith /* 74d5c3842bSBarry Smith SNESView_NRichardson - Prints info from the SNESRichardson data structure. 75028b6e76SBarry Smith 76028b6e76SBarry Smith Input Parameters: 77028b6e76SBarry Smith + SNES - the SNES context 78028b6e76SBarry Smith - viewer - visualization context 79028b6e76SBarry Smith 80028b6e76SBarry Smith Application Interface Routine: SNESView() 81028b6e76SBarry Smith */ 82028b6e76SBarry Smith #undef __FUNCT__ 83d5c3842bSBarry Smith #define __FUNCT__ "SNESView_NRichardson" 84d5c3842bSBarry Smith static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer) 85028b6e76SBarry Smith { 86028b6e76SBarry Smith PetscBool iascii; 87028b6e76SBarry Smith PetscErrorCode ierr; 88028b6e76SBarry Smith 89028b6e76SBarry Smith PetscFunctionBegin; 90028b6e76SBarry Smith ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 91028b6e76SBarry Smith if (iascii) { 92ea630c6eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer," richardson variant: %s\n", SNESLineSearchTypeName(snes->ls_type));CHKERRQ(ierr); 93028b6e76SBarry Smith } 94028b6e76SBarry Smith PetscFunctionReturn(0); 95028b6e76SBarry Smith } 96028b6e76SBarry Smith 97028b6e76SBarry Smith #undef __FUNCT__ 98d5c3842bSBarry Smith #define __FUNCT__ "SNESLineSearchQuadratic_NRichardson" 99d5c3842bSBarry Smith PetscErrorCode SNESLineSearchQuadratic_NRichardson(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec G,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag) 100028b6e76SBarry Smith { 101028b6e76SBarry Smith PetscInt i; 102028b6e76SBarry Smith PetscReal alphas[3] = {0.0, 0.5, 1.0}; 103028b6e76SBarry Smith PetscReal norms[3]; 104028b6e76SBarry Smith PetscReal alpha,a,b; 105028b6e76SBarry Smith PetscErrorCode ierr; 106028b6e76SBarry Smith 107028b6e76SBarry Smith PetscFunctionBegin; 108028b6e76SBarry Smith norms[0] = fnorm; 109028b6e76SBarry Smith for(i=1; i < 3; ++i) { 110eb1825c3SPeter Brune ierr = VecWAXPY(W, -alphas[i], Y, X);CHKERRQ(ierr); /* W = X^n - \alpha Y */ 111f9a8e2e6SPeter Brune ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr); 112f9a8e2e6SPeter Brune ierr = VecNorm(G, NORM_2, &norms[i]);CHKERRQ(ierr); 113028b6e76SBarry Smith } 114028b6e76SBarry Smith for(i = 0; i < 3; ++i) { 115028b6e76SBarry Smith norms[i] = PetscSqr(norms[i]); 116028b6e76SBarry Smith } 117028b6e76SBarry Smith /* Fit a quadratic: 118028b6e76SBarry Smith If we have x_{0,1,2} = 0, x_1, x_2 which generate norms y_{0,1,2} 119028b6e76SBarry Smith a = (x_1 y_2 - x_2 y_1 + (x_2 - x_1) y_0)/(x^2_2 x_1 - x_2 x^2_1) 120028b6e76SBarry Smith b = (x^2_1 y_2 - x^2_2 y_1 + (x^2_2 - x^2_1) y_0)/(x_2 x^2_1 - x^2_2 x_1) 121028b6e76SBarry Smith c = y_0 122028b6e76SBarry Smith x_min = -b/2a 123028b6e76SBarry Smith 124028b6e76SBarry Smith If we let x_{0,1,2} = 0, 0.5, 1.0 125028b6e76SBarry Smith a = 2 y_2 - 4 y_1 + 2 y_0 126028b6e76SBarry Smith b = -y_2 + 4 y_1 - 3 y_0 127028b6e76SBarry Smith c = y_0 128028b6e76SBarry Smith */ 129028b6e76SBarry Smith a = (alphas[1]*norms[2] - alphas[2]*norms[1] + (alphas[2] - alphas[1])*norms[0])/(PetscSqr(alphas[2])*alphas[1] - alphas[2]*PetscSqr(alphas[1])); 130028b6e76SBarry Smith b = (PetscSqr(alphas[1])*norms[2] - PetscSqr(alphas[2])*norms[1] + (PetscSqr(alphas[2]) - PetscSqr(alphas[1]))*norms[0])/(alphas[2]*PetscSqr(alphas[1]) - PetscSqr(alphas[2])*alphas[1]); 131028b6e76SBarry Smith /* Check for positive a (concave up) */ 132028b6e76SBarry Smith if (a >= 0.0) { 133028b6e76SBarry Smith alpha = -b/(2.0*a); 134028b6e76SBarry Smith alpha = PetscMin(alpha, alphas[2]); 135028b6e76SBarry Smith alpha = PetscMax(alpha, alphas[0]); 136028b6e76SBarry Smith } else { 137028b6e76SBarry Smith alpha = 1.0; 138028b6e76SBarry Smith } 139ea630c6eSPeter Brune if (snes->ls_monitor) { 140ea630c6eSPeter Brune ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 141cf22de31SPeter Brune ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line search: norms[0] = %g, norms[1] = %g, norms[2] = %g alpha %g\n", 142cf22de31SPeter Brune PetscSqrtReal(norms[0]),PetscSqrtReal(norms[1]),PetscSqrtReal(norms[2]),alpha);CHKERRQ(ierr); 143ea630c6eSPeter Brune ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 144028b6e76SBarry Smith } 145f9a8e2e6SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 146eb1825c3SPeter Brune ierr = VecAXPY(W, -alpha, Y);CHKERRQ(ierr); 147028b6e76SBarry Smith if (alpha != 1.0) { 148f9a8e2e6SPeter Brune ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr); 149f9a8e2e6SPeter Brune ierr = VecNorm(G, NORM_2, gnorm);CHKERRQ(ierr); 150028b6e76SBarry Smith } else { 151028b6e76SBarry Smith *gnorm = PetscSqrtReal(norms[2]); 152028b6e76SBarry Smith } 153028b6e76SBarry Smith if (alpha == 0.0) *flag = PETSC_FALSE; 154028b6e76SBarry Smith else *flag = PETSC_TRUE; 155028b6e76SBarry Smith PetscFunctionReturn(0); 156028b6e76SBarry Smith } 157028b6e76SBarry Smith 158028b6e76SBarry Smith /* 159d5c3842bSBarry Smith SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method. 160028b6e76SBarry Smith 161028b6e76SBarry Smith Input Parameters: 162028b6e76SBarry Smith . snes - the SNES context 163028b6e76SBarry Smith 164028b6e76SBarry Smith Output Parameter: 165028b6e76SBarry Smith . outits - number of iterations until termination 166028b6e76SBarry Smith 167028b6e76SBarry Smith Application Interface Routine: SNESSolve() 168028b6e76SBarry Smith */ 169028b6e76SBarry Smith #undef __FUNCT__ 170d5c3842bSBarry Smith #define __FUNCT__ "SNESSolve_NRichardson" 171d5c3842bSBarry Smith PetscErrorCode SNESSolve_NRichardson(SNES snes) 172028b6e76SBarry Smith { 173f9a8e2e6SPeter Brune Vec X, Y, F, W, G; 174028b6e76SBarry Smith PetscReal fnorm; 175028b6e76SBarry Smith PetscInt maxits, i; 176028b6e76SBarry Smith PetscErrorCode ierr; 177028b6e76SBarry Smith SNESConvergedReason reason; 178028b6e76SBarry Smith 179028b6e76SBarry Smith PetscFunctionBegin; 180028b6e76SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 181028b6e76SBarry Smith 182028b6e76SBarry Smith maxits = snes->max_its; /* maximum number of iterations */ 183028b6e76SBarry Smith X = snes->vec_sol; /* X^n */ 184028b6e76SBarry Smith Y = snes->vec_sol_update; /* \tilde X */ 185028b6e76SBarry Smith F = snes->vec_func; /* residual vector */ 186028b6e76SBarry Smith W = snes->work[0]; /* work vector */ 187f9a8e2e6SPeter Brune G = snes->work[1]; /* line search function vector */ 188028b6e76SBarry Smith 189028b6e76SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 190028b6e76SBarry Smith snes->iter = 0; 191028b6e76SBarry Smith snes->norm = 0.; 192028b6e76SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 193028b6e76SBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 194028b6e76SBarry Smith if (snes->domainerror) { 195028b6e76SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 196028b6e76SBarry Smith PetscFunctionReturn(0); 197028b6e76SBarry Smith } 198028b6e76SBarry Smith ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 199028b6e76SBarry Smith if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 200028b6e76SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 201028b6e76SBarry Smith snes->norm = fnorm; 202028b6e76SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 203028b6e76SBarry Smith SNESLogConvHistory(snes,fnorm,0); 204028b6e76SBarry Smith ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 205028b6e76SBarry Smith 206028b6e76SBarry Smith /* set parameter for default relative tolerance convergence test */ 207028b6e76SBarry Smith snes->ttol = fnorm*snes->rtol; 208028b6e76SBarry Smith /* test convergence */ 209028b6e76SBarry Smith ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 210028b6e76SBarry Smith if (snes->reason) PetscFunctionReturn(0); 211028b6e76SBarry Smith 212028b6e76SBarry Smith for(i = 0; i < maxits; i++) { 213028b6e76SBarry Smith PetscBool lsSuccess = PETSC_TRUE; 214028b6e76SBarry Smith PetscReal dummyNorm; 215028b6e76SBarry Smith 216028b6e76SBarry Smith /* Call general purpose update function */ 217028b6e76SBarry Smith if (snes->ops->update) { 218028b6e76SBarry Smith ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 219028b6e76SBarry Smith } 22083947a81SPeter Brune else if (snes->pc) { 221028b6e76SBarry Smith ierr = VecCopy(X,Y);CHKERRQ(ierr); 222c90fad12SPeter Brune ierr = SNESSolve(snes->pc, snes->vec_rhs, Y);CHKERRQ(ierr); 223028b6e76SBarry Smith ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 22434b225ddSBarry Smith if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 225028b6e76SBarry Smith snes->reason = SNES_DIVERGED_INNER; 226028b6e76SBarry Smith PetscFunctionReturn(0); 227028b6e76SBarry Smith } 228e5583369SPeter Brune ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr); 229a3685310SPeter Brune } else { 230a3685310SPeter Brune ierr = VecCopy(F,Y);CHKERRQ(ierr); 231028b6e76SBarry Smith } 23205b53524SPeter Brune ierr = SNESLineSearchPreCheckApply(snes, X, Y, PETSC_NULL);CHKERRQ(ierr); 23305b53524SPeter Brune ierr = SNESLineSearchApply(snes, X, F, Y, fnorm, 0.0, W, G, &dummyNorm, &fnorm, &lsSuccess);CHKERRQ(ierr); 234028b6e76SBarry Smith if (!lsSuccess) { 235028b6e76SBarry Smith if (++snes->numFailures >= snes->maxFailures) { 236028b6e76SBarry Smith snes->reason = SNES_DIVERGED_LINE_SEARCH; 237028b6e76SBarry Smith break; 238028b6e76SBarry Smith } 239028b6e76SBarry Smith } 240028b6e76SBarry Smith if (snes->nfuncs >= snes->max_funcs) { 241028b6e76SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 242028b6e76SBarry Smith break; 243028b6e76SBarry Smith } 244028b6e76SBarry Smith if (snes->domainerror) { 245028b6e76SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 246028b6e76SBarry Smith PetscFunctionReturn(0); 247028b6e76SBarry Smith } 248f9a8e2e6SPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 249f9a8e2e6SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 250028b6e76SBarry Smith /* Monitor convergence */ 251028b6e76SBarry Smith ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 252028b6e76SBarry Smith snes->iter = i+1; 253028b6e76SBarry Smith snes->norm = fnorm; 254028b6e76SBarry Smith ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 255028b6e76SBarry Smith SNESLogConvHistory(snes,snes->norm,0); 256028b6e76SBarry Smith ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 257028b6e76SBarry Smith /* Test for convergence */ 258028b6e76SBarry Smith ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 259028b6e76SBarry Smith if (snes->reason) break; 260028b6e76SBarry Smith } 261028b6e76SBarry Smith if (i == maxits) { 262028b6e76SBarry Smith ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 263028b6e76SBarry Smith if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 264028b6e76SBarry Smith } 265028b6e76SBarry Smith PetscFunctionReturn(0); 266028b6e76SBarry Smith } 267028b6e76SBarry Smith 268028b6e76SBarry Smith 26992c02d66SPeter Brune EXTERN_C_BEGIN 27092c02d66SPeter Brune #undef __FUNCT__ 27192c02d66SPeter Brune #define __FUNCT__ "SNESLineSearchSetType_NRichardson" 27292c02d66SPeter Brune PetscErrorCode SNESLineSearchSetType_NRichardson(SNES snes, SNESLineSearchType type) 27392c02d66SPeter Brune { 27492c02d66SPeter Brune PetscErrorCode ierr; 27592c02d66SPeter Brune PetscFunctionBegin; 27692c02d66SPeter Brune 27792c02d66SPeter Brune switch (type) { 27892c02d66SPeter Brune case SNES_LS_BASIC: 279f9a8e2e6SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 28092c02d66SPeter Brune break; 28192c02d66SPeter Brune case SNES_LS_BASIC_NONORMS: 282f9a8e2e6SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 28392c02d66SPeter Brune break; 28492c02d66SPeter Brune case SNES_LS_QUADRATIC: 28588d374b2SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr); 28692c02d66SPeter Brune break; 287*f3aaf736SPeter Brune case SNES_LS_CRITICAL: 288*f3aaf736SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchCriticalSecant,PETSC_NULL);CHKERRQ(ierr); 289f9a8e2e6SPeter Brune break; 29092c02d66SPeter Brune default: 29192c02d66SPeter Brune SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type"); 29292c02d66SPeter Brune break; 29392c02d66SPeter Brune } 29492c02d66SPeter Brune snes->ls_type = type; 29592c02d66SPeter Brune PetscFunctionReturn(0); 29692c02d66SPeter Brune } 29792c02d66SPeter Brune EXTERN_C_END 29892c02d66SPeter Brune 299028b6e76SBarry Smith /*MC 300d5c3842bSBarry Smith SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration. 301028b6e76SBarry Smith 302028b6e76SBarry Smith Level: beginner 303028b6e76SBarry Smith 304028b6e76SBarry Smith Options Database: 305*f3aaf736SPeter Brune . -snes_ls <basic,basicnormnorms,quadratic,critical> Line search type. 306028b6e76SBarry Smith 307af60355fSPeter Brune Notes: If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda 308af60355fSPeter Brune (F(x^n) - b) where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search. If 309af60355fSPeter Brune an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetPC()) then the inner 310af60355fSPeter Brune solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n} 311af60355fSPeter Brune where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver. 312028b6e76SBarry Smith 313028b6e76SBarry Smith This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls 314028b6e76SBarry Smith 3151867fe5bSPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNGMRES, SNESQN, SNESNCG 316028b6e76SBarry Smith M*/ 317028b6e76SBarry Smith EXTERN_C_BEGIN 318028b6e76SBarry Smith #undef __FUNCT__ 319d5c3842bSBarry Smith #define __FUNCT__ "SNESCreate_NRichardson" 320d5c3842bSBarry Smith PetscErrorCode SNESCreate_NRichardson(SNES snes) 321028b6e76SBarry Smith { 32292c02d66SPeter Brune PetscErrorCode ierr; 323028b6e76SBarry Smith PetscFunctionBegin; 324d5c3842bSBarry Smith snes->ops->destroy = SNESDestroy_NRichardson; 325d5c3842bSBarry Smith snes->ops->setup = SNESSetUp_NRichardson; 326d5c3842bSBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NRichardson; 327d5c3842bSBarry Smith snes->ops->view = SNESView_NRichardson; 328d5c3842bSBarry Smith snes->ops->solve = SNESSolve_NRichardson; 329d5c3842bSBarry Smith snes->ops->reset = SNESReset_NRichardson; 330028b6e76SBarry Smith 331028b6e76SBarry Smith snes->usesksp = PETSC_FALSE; 33242f4f86dSBarry Smith snes->usespc = PETSC_TRUE; 333028b6e76SBarry Smith 3340e444f03SPeter Brune snes->max_funcs = 30000; 3350e444f03SPeter Brune snes->max_its = 10000; 3360e444f03SPeter Brune 33792c02d66SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_NRichardson",SNESLineSearchSetType_NRichardson);CHKERRQ(ierr); 338*f3aaf736SPeter Brune ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr); 339028b6e76SBarry Smith 340028b6e76SBarry Smith PetscFunctionReturn(0); 341028b6e76SBarry Smith } 342028b6e76SBarry Smith EXTERN_C_END 343