1fef7b6d8SPeter Brune #include <../src/snes/impls/ncg/snesncgimpl.h> 2fef7b6d8SPeter Brune 3fef7b6d8SPeter Brune #undef __FUNCT__ 4fef7b6d8SPeter Brune #define __FUNCT__ "SNESReset_NCG" 5fef7b6d8SPeter Brune PetscErrorCode SNESReset_NCG(SNES snes) 6fef7b6d8SPeter Brune { 7fef7b6d8SPeter Brune 8fef7b6d8SPeter Brune PetscFunctionBegin; 9fef7b6d8SPeter Brune PetscFunctionReturn(0); 10fef7b6d8SPeter Brune } 11fef7b6d8SPeter Brune 12fef7b6d8SPeter Brune /* 13fef7b6d8SPeter Brune SNESDestroy_NCG - Destroys the private SNES_NCG context that was created with SNESCreate_NCG(). 14fef7b6d8SPeter Brune 15fef7b6d8SPeter Brune Input Parameter: 16fef7b6d8SPeter Brune . snes - the SNES context 17fef7b6d8SPeter Brune 18fef7b6d8SPeter Brune Application Interface Routine: SNESDestroy() 19fef7b6d8SPeter Brune */ 20fef7b6d8SPeter Brune #undef __FUNCT__ 21fef7b6d8SPeter Brune #define __FUNCT__ "SNESDestroy_NCG" 22fef7b6d8SPeter Brune PetscErrorCode SNESDestroy_NCG(SNES snes) 23fef7b6d8SPeter Brune { 24fef7b6d8SPeter Brune PetscErrorCode ierr; 25fef7b6d8SPeter Brune 26fef7b6d8SPeter Brune PetscFunctionBegin; 27fef7b6d8SPeter Brune ierr = SNESReset_NCG(snes);CHKERRQ(ierr); 28fef7b6d8SPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 29fef7b6d8SPeter Brune PetscFunctionReturn(0); 30fef7b6d8SPeter Brune } 31fef7b6d8SPeter Brune 32fef7b6d8SPeter Brune /* 33fef7b6d8SPeter Brune SNESSetUp_NCG - Sets up the internal data structures for the later use 34fef7b6d8SPeter Brune of the SNESNCG nonlinear solver. 35fef7b6d8SPeter Brune 36fef7b6d8SPeter Brune Input Parameters: 37fef7b6d8SPeter Brune + snes - the SNES context 38fef7b6d8SPeter Brune - x - the solution vector 39fef7b6d8SPeter Brune 40fef7b6d8SPeter Brune Application Interface Routine: SNESSetUp() 41fef7b6d8SPeter Brune */ 42fef7b6d8SPeter Brune #undef __FUNCT__ 43fef7b6d8SPeter Brune #define __FUNCT__ "SNESSetUp_NCG" 44fef7b6d8SPeter Brune PetscErrorCode SNESSetUp_NCG(SNES snes) 45fef7b6d8SPeter Brune { 46fef7b6d8SPeter Brune PetscErrorCode ierr; 47fef7b6d8SPeter Brune 48fef7b6d8SPeter Brune PetscFunctionBegin; 49dfb256c7SPeter Brune ierr = SNESDefaultGetWork(snes,4);CHKERRQ(ierr); 50fef7b6d8SPeter Brune PetscFunctionReturn(0); 51fef7b6d8SPeter Brune } 52fef7b6d8SPeter Brune /* 53fef7b6d8SPeter Brune SNESSetFromOptions_NCG - Sets various parameters for the SNESLS method. 54fef7b6d8SPeter Brune 55fef7b6d8SPeter Brune Input Parameter: 56fef7b6d8SPeter Brune . snes - the SNES context 57fef7b6d8SPeter Brune 58fef7b6d8SPeter Brune Application Interface Routine: SNESSetFromOptions() 59fef7b6d8SPeter Brune */ 60fef7b6d8SPeter Brune #undef __FUNCT__ 61fef7b6d8SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NCG" 62fef7b6d8SPeter Brune static PetscErrorCode SNESSetFromOptions_NCG(SNES snes) 63fef7b6d8SPeter Brune { 64dfb256c7SPeter Brune SNES_NCG *ncg = (SNES_NCG *)snes->data; 65dfb256c7SPeter Brune const char *betas[] = {"FR","PRP","HS", "DY", "CD"}; 66fef7b6d8SPeter Brune PetscErrorCode ierr; 67dfb256c7SPeter Brune PetscBool debug, flg; 68dfb256c7SPeter Brune PetscInt indx; 69fef7b6d8SPeter Brune 70fef7b6d8SPeter Brune PetscFunctionBegin; 71fef7b6d8SPeter Brune ierr = PetscOptionsHead("SNES NCG options");CHKERRQ(ierr); 72dfb256c7SPeter Brune ierr = PetscOptionsBool("-snes_ncg_monitor", "Monitor NCG iterations", "SNES", ncg->monitor ? PETSC_TRUE: PETSC_FALSE, &debug, PETSC_NULL);CHKERRQ(ierr); 73dfb256c7SPeter Brune ierr = PetscOptionsEList("-snes_ncg_type","Nonlinear CG update used","", betas, 5, "FR",&indx, &flg);CHKERRQ(ierr); 74dfb256c7SPeter Brune if (flg) { 75dfb256c7SPeter Brune ncg->betatype = indx; 76dfb256c7SPeter Brune } 77dfb256c7SPeter Brune if (debug) { 78dfb256c7SPeter Brune ncg->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 79dfb256c7SPeter Brune } 80fef7b6d8SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 81fef7b6d8SPeter Brune PetscFunctionReturn(0); 82fef7b6d8SPeter Brune } 83fef7b6d8SPeter Brune 84fef7b6d8SPeter Brune /* 85fef7b6d8SPeter Brune SNESView_NCG - Prints info from the SNESNCG data structure. 86fef7b6d8SPeter Brune 87fef7b6d8SPeter Brune Input Parameters: 88fef7b6d8SPeter Brune + SNES - the SNES context 89fef7b6d8SPeter Brune - viewer - visualization context 90fef7b6d8SPeter Brune 91fef7b6d8SPeter Brune Application Interface Routine: SNESView() 92fef7b6d8SPeter Brune */ 93fef7b6d8SPeter Brune #undef __FUNCT__ 94fef7b6d8SPeter Brune #define __FUNCT__ "SNESView_NCG" 95fef7b6d8SPeter Brune static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer) 96fef7b6d8SPeter Brune { 97fef7b6d8SPeter Brune PetscBool iascii; 98fef7b6d8SPeter Brune PetscErrorCode ierr; 99fef7b6d8SPeter Brune 100fef7b6d8SPeter Brune PetscFunctionBegin; 101fef7b6d8SPeter Brune ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 102fef7b6d8SPeter Brune if (iascii) { 103ea630c6eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer," line search type variant: %s\n", SNESLineSearchTypeName(snes->ls_type));CHKERRQ(ierr); 104fef7b6d8SPeter Brune } 105fef7b6d8SPeter Brune PetscFunctionReturn(0); 106fef7b6d8SPeter Brune } 107fef7b6d8SPeter Brune 108fef7b6d8SPeter Brune #undef __FUNCT__ 109ea630c6eSPeter Brune #define __FUNCT__ "SNESLineSearchExactLinear_NCG" 110ea630c6eSPeter Brune PetscErrorCode SNESLineSearchExactLinear_NCG(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec G,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag) 111ea630c6eSPeter Brune { 112ea630c6eSPeter Brune PetscScalar alpha, ptAp; 113ea630c6eSPeter Brune PetscErrorCode ierr; 114ea630c6eSPeter Brune /* SNES_NCG *ncg = (SNES_NCG *) snes->data; */ 115ea630c6eSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 116ea630c6eSPeter Brune 117ea630c6eSPeter Brune PetscFunctionBegin; 118ea630c6eSPeter Brune /* 119ea630c6eSPeter Brune 120d2140ca3SPeter Brune The exact step size for unpreconditioned linear CG is just: 121ea630c6eSPeter Brune 122ea630c6eSPeter Brune alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy) 123ea630c6eSPeter Brune 124ea630c6eSPeter Brune */ 125ea630c6eSPeter Brune 126a5892487SPeter Brune ierr = SNESComputeJacobian(snes, X, &snes->jacobian, &snes->jacobian_pre, &flg);CHKERRQ(ierr); 127ea630c6eSPeter Brune ierr = VecDot(F, F, &alpha);CHKERRQ(ierr); 128ea630c6eSPeter Brune ierr = MatMult(snes->jacobian, Y, W);CHKERRQ(ierr); 129ea630c6eSPeter Brune ierr = VecDot(Y, W, &ptAp);CHKERRQ(ierr); 130ea630c6eSPeter Brune alpha = alpha / ptAp; 131d9fe6452SJed Brown ierr = PetscPrintf(((PetscObject)snes)->comm, "alpha: %G\n", PetscRealPart(alpha));CHKERRQ(ierr); 132a5892487SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 133a5892487SPeter Brune ierr = VecAXPY(W, alpha, Y);CHKERRQ(ierr); 134a5892487SPeter Brune ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr); 135a5892487SPeter Brune ierr = VecNorm(G, NORM_2, gnorm);CHKERRQ(ierr); 136ea630c6eSPeter Brune PetscFunctionReturn(0); 137ea630c6eSPeter Brune } 138ea630c6eSPeter Brune 13970d3b23bSPeter Brune EXTERN_C_BEGIN 14070d3b23bSPeter Brune #undef __FUNCT__ 14170d3b23bSPeter Brune #define __FUNCT__ "SNESLineSearchSetType_NCG" 14270d3b23bSPeter Brune PetscErrorCode SNESLineSearchSetType_NCG(SNES snes, SNESLineSearchType type) 14370d3b23bSPeter Brune { 14470d3b23bSPeter Brune PetscErrorCode ierr; 14570d3b23bSPeter Brune PetscFunctionBegin; 14670d3b23bSPeter Brune 14770d3b23bSPeter Brune switch (type) { 14870d3b23bSPeter Brune case SNES_LS_BASIC: 14970d3b23bSPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 15070d3b23bSPeter Brune break; 15170d3b23bSPeter Brune case SNES_LS_BASIC_NONORMS: 15270d3b23bSPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 15370d3b23bSPeter Brune break; 15470d3b23bSPeter Brune case SNES_LS_QUADRATIC: 155d2140ca3SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr); 15670d3b23bSPeter Brune break; 15770d3b23bSPeter Brune case SNES_LS_TEST: 15870d3b23bSPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchExactLinear_NCG,PETSC_NULL);CHKERRQ(ierr); 15970d3b23bSPeter Brune break; 160cf9edfecSPeter Brune case SNES_LS_SECANT: 161cf9edfecSPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchSecant,PETSC_NULL);CHKERRQ(ierr); 162cf9edfecSPeter Brune break; 16370d3b23bSPeter Brune default: 16470d3b23bSPeter Brune SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type"); 16570d3b23bSPeter Brune break; 16670d3b23bSPeter Brune } 16770d3b23bSPeter Brune snes->ls_type = type; 16870d3b23bSPeter Brune PetscFunctionReturn(0); 16970d3b23bSPeter Brune } 17070d3b23bSPeter Brune EXTERN_C_END 17170d3b23bSPeter Brune 172fef7b6d8SPeter Brune /* 173fef7b6d8SPeter Brune SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method. 174fef7b6d8SPeter Brune 175fef7b6d8SPeter Brune Input Parameters: 176fef7b6d8SPeter Brune . snes - the SNES context 177fef7b6d8SPeter Brune 178fef7b6d8SPeter Brune Output Parameter: 179fef7b6d8SPeter Brune . outits - number of iterations until termination 180fef7b6d8SPeter Brune 181fef7b6d8SPeter Brune Application Interface Routine: SNESSolve() 182fef7b6d8SPeter Brune */ 183fef7b6d8SPeter Brune #undef __FUNCT__ 184fef7b6d8SPeter Brune #define __FUNCT__ "SNESSolve_NCG" 185fef7b6d8SPeter Brune PetscErrorCode SNESSolve_NCG(SNES snes) 186fef7b6d8SPeter Brune { 187dfb256c7SPeter Brune SNES_NCG *ncg = (SNES_NCG *)snes->data; 18829c94e34SPeter Brune Vec X, dX, lX, F, W, B, G, dXold; 189dfb256c7SPeter Brune PetscReal fnorm, gnorm, dummyNorm, beta = 0.0; 19029c94e34SPeter Brune PetscScalar dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold; 191fef7b6d8SPeter Brune PetscInt maxits, i; 192fef7b6d8SPeter Brune PetscErrorCode ierr; 193fef7b6d8SPeter Brune SNESConvergedReason reason; 194fef7b6d8SPeter Brune PetscBool lsSuccess = PETSC_TRUE; 195fef7b6d8SPeter Brune PetscFunctionBegin; 196fef7b6d8SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 197fef7b6d8SPeter Brune 198fef7b6d8SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 199fef7b6d8SPeter Brune X = snes->vec_sol; /* X^n */ 20029c94e34SPeter Brune dXold = snes->work[0]; /* The previous iterate of X */ 201169e2be8SPeter Brune dX = snes->work[1]; /* the steepest direction */ 202169e2be8SPeter Brune lX = snes->vec_sol_update; /* the conjugate direction */ 203fef7b6d8SPeter Brune F = snes->vec_func; /* residual vector */ 204dfb256c7SPeter Brune W = snes->work[2]; /* work vector and solution output for the line search */ 205302dbbaeSPeter Brune B = snes->vec_rhs; /* the right hand side */ 206dfb256c7SPeter Brune G = snes->work[3]; /* the line search result function evaluation */ 207fef7b6d8SPeter Brune 208fef7b6d8SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 209fef7b6d8SPeter Brune snes->iter = 0; 210fef7b6d8SPeter Brune snes->norm = 0.; 211fef7b6d8SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 212fef7b6d8SPeter Brune 213fef7b6d8SPeter Brune /* compute the initial function and preconditioned update delX */ 214fef7b6d8SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 215fef7b6d8SPeter Brune if (snes->domainerror) { 216fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 217fef7b6d8SPeter Brune PetscFunctionReturn(0); 218fef7b6d8SPeter Brune } 219fef7b6d8SPeter Brune /* convergence test */ 220fef7b6d8SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 221fef7b6d8SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 222fef7b6d8SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 223fef7b6d8SPeter Brune snes->norm = fnorm; 224fef7b6d8SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 225fef7b6d8SPeter Brune SNESLogConvHistory(snes,fnorm,0); 226fef7b6d8SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 227fef7b6d8SPeter Brune 228fef7b6d8SPeter Brune /* set parameter for default relative tolerance convergence test */ 229fef7b6d8SPeter Brune snes->ttol = fnorm*snes->rtol; 230fef7b6d8SPeter Brune /* test convergence */ 231fef7b6d8SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 232fef7b6d8SPeter Brune if (snes->reason) PetscFunctionReturn(0); 233fef7b6d8SPeter Brune 234fef7b6d8SPeter Brune /* Call general purpose update function */ 235fef7b6d8SPeter Brune if (snes->ops->update) { 236fef7b6d8SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 237fef7b6d8SPeter Brune } 238fef7b6d8SPeter Brune 239fef7b6d8SPeter Brune /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */ 24083947a81SPeter Brune if (snes->pc) { 24129c94e34SPeter Brune ierr = VecCopy(X, dX);CHKERRQ(ierr); 24229c94e34SPeter Brune ierr = SNESSolve(snes->pc, B, dX);CHKERRQ(ierr); 243fef7b6d8SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 24451e86f29SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 245fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_INNER; 246fef7b6d8SPeter Brune PetscFunctionReturn(0); 247fef7b6d8SPeter Brune } 24829c94e34SPeter Brune ierr = VecAYPX(dX,-1.0,X);CHKERRQ(ierr); 2495d10c001SPeter Brune } else { 25029c94e34SPeter Brune ierr = VecCopy(F, dX);CHKERRQ(ierr); 251fef7b6d8SPeter Brune } 25229c94e34SPeter Brune ierr = VecCopy(dX, lX);CHKERRQ(ierr); 25329c94e34SPeter Brune ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr); 25409c08436SPeter Brune for(i = 1; i < maxits + 1; i++) { 255fef7b6d8SPeter Brune lsSuccess = PETSC_TRUE; 25629c94e34SPeter Brune /* some update types require the old update direction or conjugate direction */ 25729c94e34SPeter Brune if (ncg->betatype == 1 || ncg->betatype == 2 || ncg->betatype == 3 || ncg->betatype == 4) { /* prp, hs, dy, cd need dXold*/ 25829c94e34SPeter Brune ierr = VecCopy(dX, dXold);CHKERRQ(ierr); 25929c94e34SPeter Brune } 260a5892487SPeter Brune ierr = (*snes->ops->linesearch)(snes, snes->lsP, X, F, lX, fnorm, 0.0, G, W, &dummyNorm, &gnorm, &lsSuccess);CHKERRQ(ierr); 261fef7b6d8SPeter Brune if (!lsSuccess) { 262fef7b6d8SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 263fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 2645d10c001SPeter Brune PetscFunctionReturn(0); 265fef7b6d8SPeter Brune } 266fef7b6d8SPeter Brune } 267fef7b6d8SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 268fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 2695d10c001SPeter Brune PetscFunctionReturn(0); 270fef7b6d8SPeter Brune } 271fef7b6d8SPeter Brune if (snes->domainerror) { 272fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 273fef7b6d8SPeter Brune PetscFunctionReturn(0); 274fef7b6d8SPeter Brune } 275302dbbaeSPeter Brune 276a5892487SPeter Brune /* copy over the solution */ 277a5892487SPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 278a5892487SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 279a5892487SPeter Brune fnorm = gnorm; 280a5892487SPeter Brune 281fef7b6d8SPeter Brune /* Monitor convergence */ 282fef7b6d8SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 283169e2be8SPeter Brune snes->iter = i; 284fef7b6d8SPeter Brune snes->norm = fnorm; 285fef7b6d8SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 286fef7b6d8SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 287fef7b6d8SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 288302dbbaeSPeter Brune 289fef7b6d8SPeter Brune /* Test for convergence */ 290fef7b6d8SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2915d10c001SPeter Brune if (snes->reason) PetscFunctionReturn(0); 292302dbbaeSPeter Brune 293302dbbaeSPeter Brune /* Call general purpose update function */ 294302dbbaeSPeter Brune if (snes->ops->update) { 295302dbbaeSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 296302dbbaeSPeter Brune } 29783947a81SPeter Brune if (snes->pc) { 298302dbbaeSPeter Brune ierr = VecCopy(X,dX);CHKERRQ(ierr); 299cf949ffaSPeter Brune ierr = SNESSolve(snes->pc, B, dX);CHKERRQ(ierr); 300302dbbaeSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 30151e86f29SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 302302dbbaeSPeter Brune snes->reason = SNES_DIVERGED_INNER; 303302dbbaeSPeter Brune PetscFunctionReturn(0); 304302dbbaeSPeter Brune } 305eb1825c3SPeter Brune ierr = VecAYPX(dX,-1.0,X);CHKERRQ(ierr); 306a3685310SPeter Brune } else { 307a3685310SPeter Brune ierr = VecCopy(F,dX);CHKERRQ(ierr); 308302dbbaeSPeter Brune } 309302dbbaeSPeter Brune 31029c94e34SPeter Brune /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/ 311dfb256c7SPeter Brune switch(ncg->betatype) { 312dfb256c7SPeter Brune case 0: /* Fletcher-Reeves */ 31329c94e34SPeter Brune dXolddotdXold = dXdotdX; 31429c94e34SPeter Brune ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr); 31529c94e34SPeter Brune beta = PetscRealPart(dXdotdX / dXolddotdXold); 316dfb256c7SPeter Brune break; 317dfb256c7SPeter Brune case 1: /* Polak-Ribiere-Poylak */ 31829c94e34SPeter Brune dXolddotdXold = dXdotdX; 31929c94e34SPeter Brune ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr); 32029c94e34SPeter Brune ierr = VecDot(dXold, dX, &dXdotdXold);CHKERRQ(ierr); 32125f8a51bSPeter Brune beta = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold)); 322dfb256c7SPeter Brune if (beta < 0.0) beta = 0.0; /* restart */ 323dfb256c7SPeter Brune break; 324dfb256c7SPeter Brune case 2: /* Hestenes-Stiefel */ 32529c94e34SPeter Brune ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr); 32629c94e34SPeter Brune ierr = VecDot(dX, dXold, &dXdotdXold);CHKERRQ(ierr); 32729c94e34SPeter Brune ierr = VecDot(lX, dX, &lXdotdX);CHKERRQ(ierr); 32829c94e34SPeter Brune ierr = VecDot(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 32929c94e34SPeter Brune beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold)); 330dfb256c7SPeter Brune break; 331dfb256c7SPeter Brune case 3: /* Dai-Yuan */ 33229c94e34SPeter Brune ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr); 33329c94e34SPeter Brune ierr = VecDot(lX, dX, &lXdotdXold);CHKERRQ(ierr); 33429c94e34SPeter Brune ierr = VecDot(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 33529c94e34SPeter Brune beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));CHKERRQ(ierr); 336dfb256c7SPeter Brune break; 337dfb256c7SPeter Brune case 4: /* Conjugate Descent */ 33829c94e34SPeter Brune ierr = VecDot(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 33929c94e34SPeter Brune ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr); 34029c94e34SPeter Brune beta = PetscRealPart(dXdotdX / lXdotdXold);CHKERRQ(ierr); 341dfb256c7SPeter Brune break; 342dfb256c7SPeter Brune } 343dfb256c7SPeter Brune if (ncg->monitor) { 344646217ecSPeter Brune ierr = PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", beta);CHKERRQ(ierr); 345dfb256c7SPeter Brune } 346dfb256c7SPeter Brune ierr = VecAYPX(lX, beta, dX);CHKERRQ(ierr); 347fef7b6d8SPeter Brune } 348fef7b6d8SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 349fef7b6d8SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 350fef7b6d8SPeter Brune PetscFunctionReturn(0); 351fef7b6d8SPeter Brune } 352fef7b6d8SPeter Brune 353fef7b6d8SPeter Brune /*MC 354fef7b6d8SPeter Brune SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems. 355fef7b6d8SPeter Brune 356fef7b6d8SPeter Brune Level: beginner 357fef7b6d8SPeter Brune 358fef7b6d8SPeter Brune Options Database: 359fef7b6d8SPeter Brune + -snes_ls_damping - damping factor to apply to F(x) (used only if -snes_ls is basic or basicnonorms) 360fef7b6d8SPeter Brune - -snes_ls <basic,basicnormnorms,quadratic> 361fef7b6d8SPeter Brune 362fef7b6d8SPeter Brune Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate 363fef7b6d8SPeter Brune gradient method. This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise 364fef7b6d8SPeter Brune chooses the initial search direction as F(x) for the initial guess x. 365fef7b6d8SPeter Brune 366fef7b6d8SPeter Brune 367fef7b6d8SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNGMRES, SNESNQN 368fef7b6d8SPeter Brune M*/ 369fef7b6d8SPeter Brune EXTERN_C_BEGIN 370fef7b6d8SPeter Brune #undef __FUNCT__ 371fef7b6d8SPeter Brune #define __FUNCT__ "SNESCreate_NCG" 372fef7b6d8SPeter Brune PetscErrorCode SNESCreate_NCG(SNES snes) 373fef7b6d8SPeter Brune { 374fef7b6d8SPeter Brune PetscErrorCode ierr; 375ea630c6eSPeter Brune SNES_NCG * neP; 376fef7b6d8SPeter Brune 377fef7b6d8SPeter Brune PetscFunctionBegin; 378fef7b6d8SPeter Brune snes->ops->destroy = SNESDestroy_NCG; 379fef7b6d8SPeter Brune snes->ops->setup = SNESSetUp_NCG; 380fef7b6d8SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NCG; 381fef7b6d8SPeter Brune snes->ops->view = SNESView_NCG; 382fef7b6d8SPeter Brune snes->ops->solve = SNESSolve_NCG; 383fef7b6d8SPeter Brune snes->ops->reset = SNESReset_NCG; 384fef7b6d8SPeter Brune 385fef7b6d8SPeter Brune snes->usesksp = PETSC_FALSE; 386fef7b6d8SPeter Brune snes->usespc = PETSC_TRUE; 387fef7b6d8SPeter Brune 388*0e444f03SPeter Brune snes->max_funcs = 30000; 389*0e444f03SPeter Brune snes->max_its = 10000; 390*0e444f03SPeter Brune 391fef7b6d8SPeter Brune ierr = PetscNewLog(snes, SNES_NCG, &neP);CHKERRQ(ierr); 392fef7b6d8SPeter Brune snes->data = (void*) neP; 393dfb256c7SPeter Brune neP->monitor = PETSC_NULL; 39425f8a51bSPeter Brune neP->betatype = 4; 39570d3b23bSPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_NCG",SNESLineSearchSetType_NCG);CHKERRQ(ierr); 396d2140ca3SPeter Brune ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr); 397fef7b6d8SPeter Brune PetscFunctionReturn(0); 398fef7b6d8SPeter Brune } 399fef7b6d8SPeter Brune EXTERN_C_END 400