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 PetscErrorCode ierr; 8fef7b6d8SPeter Brune 9fef7b6d8SPeter Brune PetscFunctionBegin; 10fef7b6d8SPeter Brune if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);} 11fef7b6d8SPeter Brune PetscFunctionReturn(0); 12fef7b6d8SPeter Brune } 13fef7b6d8SPeter Brune 14fef7b6d8SPeter Brune /* 15fef7b6d8SPeter Brune SNESDestroy_NCG - Destroys the private SNES_NCG context that was created with SNESCreate_NCG(). 16fef7b6d8SPeter Brune 17fef7b6d8SPeter Brune Input Parameter: 18fef7b6d8SPeter Brune . snes - the SNES context 19fef7b6d8SPeter Brune 20fef7b6d8SPeter Brune Application Interface Routine: SNESDestroy() 21fef7b6d8SPeter Brune */ 22fef7b6d8SPeter Brune #undef __FUNCT__ 23fef7b6d8SPeter Brune #define __FUNCT__ "SNESDestroy_NCG" 24fef7b6d8SPeter Brune PetscErrorCode SNESDestroy_NCG(SNES snes) 25fef7b6d8SPeter Brune { 26fef7b6d8SPeter Brune PetscErrorCode ierr; 27fef7b6d8SPeter Brune 28fef7b6d8SPeter Brune PetscFunctionBegin; 29fef7b6d8SPeter Brune ierr = SNESReset_NCG(snes);CHKERRQ(ierr); 30fef7b6d8SPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 31fef7b6d8SPeter Brune PetscFunctionReturn(0); 32fef7b6d8SPeter Brune } 33fef7b6d8SPeter Brune 34fef7b6d8SPeter Brune /* 35fef7b6d8SPeter Brune SNESSetUp_NCG - Sets up the internal data structures for the later use 36fef7b6d8SPeter Brune of the SNESNCG nonlinear solver. 37fef7b6d8SPeter Brune 38fef7b6d8SPeter Brune Input Parameters: 39fef7b6d8SPeter Brune + snes - the SNES context 40fef7b6d8SPeter Brune - x - the solution vector 41fef7b6d8SPeter Brune 42fef7b6d8SPeter Brune Application Interface Routine: SNESSetUp() 43fef7b6d8SPeter Brune */ 44fef7b6d8SPeter Brune #undef __FUNCT__ 45fef7b6d8SPeter Brune #define __FUNCT__ "SNESSetUp_NCG" 46fef7b6d8SPeter Brune PetscErrorCode SNESSetUp_NCG(SNES snes) 47fef7b6d8SPeter Brune { 48fef7b6d8SPeter Brune PetscErrorCode ierr; 49fef7b6d8SPeter Brune 50fef7b6d8SPeter Brune PetscFunctionBegin; 51*a5892487SPeter Brune ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr); 52fef7b6d8SPeter Brune PetscFunctionReturn(0); 53fef7b6d8SPeter Brune } 54fef7b6d8SPeter Brune /* 55fef7b6d8SPeter Brune SNESSetFromOptions_NCG - Sets various parameters for the SNESLS method. 56fef7b6d8SPeter Brune 57fef7b6d8SPeter Brune Input Parameter: 58fef7b6d8SPeter Brune . snes - the SNES context 59fef7b6d8SPeter Brune 60fef7b6d8SPeter Brune Application Interface Routine: SNESSetFromOptions() 61fef7b6d8SPeter Brune */ 62fef7b6d8SPeter Brune #undef __FUNCT__ 63fef7b6d8SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NCG" 64fef7b6d8SPeter Brune static PetscErrorCode SNESSetFromOptions_NCG(SNES snes) 65fef7b6d8SPeter Brune { 66fef7b6d8SPeter Brune PetscErrorCode ierr; 67fef7b6d8SPeter Brune 68fef7b6d8SPeter Brune PetscFunctionBegin; 69fef7b6d8SPeter Brune ierr = PetscOptionsHead("SNES NCG options");CHKERRQ(ierr); 70fef7b6d8SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 71fef7b6d8SPeter Brune PetscFunctionReturn(0); 72fef7b6d8SPeter Brune } 73fef7b6d8SPeter Brune 74fef7b6d8SPeter Brune /* 75fef7b6d8SPeter Brune SNESView_NCG - Prints info from the SNESNCG data structure. 76fef7b6d8SPeter Brune 77fef7b6d8SPeter Brune Input Parameters: 78fef7b6d8SPeter Brune + SNES - the SNES context 79fef7b6d8SPeter Brune - viewer - visualization context 80fef7b6d8SPeter Brune 81fef7b6d8SPeter Brune Application Interface Routine: SNESView() 82fef7b6d8SPeter Brune */ 83fef7b6d8SPeter Brune #undef __FUNCT__ 84fef7b6d8SPeter Brune #define __FUNCT__ "SNESView_NCG" 85fef7b6d8SPeter Brune static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer) 86fef7b6d8SPeter Brune { 87fef7b6d8SPeter Brune PetscBool iascii; 88fef7b6d8SPeter Brune PetscErrorCode ierr; 89fef7b6d8SPeter Brune 90fef7b6d8SPeter Brune PetscFunctionBegin; 91fef7b6d8SPeter Brune ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 92fef7b6d8SPeter Brune if (iascii) { 93ea630c6eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer," line search type variant: %s\n", SNESLineSearchTypeName(snes->ls_type));CHKERRQ(ierr); 94fef7b6d8SPeter Brune } 95fef7b6d8SPeter Brune PetscFunctionReturn(0); 96fef7b6d8SPeter Brune } 97fef7b6d8SPeter Brune 98fef7b6d8SPeter Brune #undef __FUNCT__ 99fef7b6d8SPeter Brune #define __FUNCT__ "SNESLineSearchQuadratic_NCG" 100fef7b6d8SPeter Brune PetscErrorCode SNESLineSearchQuadratic_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) 101fef7b6d8SPeter Brune { 102fef7b6d8SPeter Brune PetscInt i; 103fef7b6d8SPeter Brune PetscReal alphas[3] = {0.0, 0.5, 1.0}; 104fef7b6d8SPeter Brune PetscReal norms[3]; 105fef7b6d8SPeter Brune PetscReal alpha,a,b; 106fef7b6d8SPeter Brune PetscErrorCode ierr; 107fef7b6d8SPeter Brune 108fef7b6d8SPeter Brune PetscFunctionBegin; 109fef7b6d8SPeter Brune norms[0] = fnorm; 110fef7b6d8SPeter Brune for(i=1; i < 3; ++i) { 111fef7b6d8SPeter Brune ierr = VecWAXPY(W, alphas[i], Y, X);CHKERRQ(ierr); /* W = X^n - \alpha Y */ 112*a5892487SPeter Brune ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr); 113*a5892487SPeter Brune ierr = VecNorm(G, NORM_2, &norms[i]);CHKERRQ(ierr); 114fef7b6d8SPeter Brune } 115fef7b6d8SPeter Brune for(i = 0; i < 3; ++i) { 116fef7b6d8SPeter Brune norms[i] = PetscSqr(norms[i]); 117fef7b6d8SPeter Brune } 118fef7b6d8SPeter Brune /* Fit a quadratic: 119fef7b6d8SPeter Brune If we have x_{0,1,2} = 0, x_1, x_2 which generate norms y_{0,1,2} 120fef7b6d8SPeter Brune 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) 121fef7b6d8SPeter Brune 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) 122fef7b6d8SPeter Brune c = y_0 123fef7b6d8SPeter Brune x_min = -b/2a 124fef7b6d8SPeter Brune 125fef7b6d8SPeter Brune If we let x_{0,1,2} = 0, 0.5, 1.0 126fef7b6d8SPeter Brune a = 2 y_2 - 4 y_1 + 2 y_0 127fef7b6d8SPeter Brune b = -y_2 + 4 y_1 - 3 y_0 128fef7b6d8SPeter Brune c = y_0 129fef7b6d8SPeter Brune */ 130fef7b6d8SPeter Brune 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])); 131fef7b6d8SPeter Brune 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]); 132fef7b6d8SPeter Brune /* Check for positive a (concave up) */ 133fef7b6d8SPeter Brune if (a >= 0.0) { 134fef7b6d8SPeter Brune alpha = -b/(2.0*a); 135fef7b6d8SPeter Brune alpha = PetscMin(alpha, alphas[2]); 136fef7b6d8SPeter Brune alpha = PetscMax(alpha, alphas[0]); 137fef7b6d8SPeter Brune } else { 138fef7b6d8SPeter Brune alpha = 1.0; 139fef7b6d8SPeter Brune } 140ea630c6eSPeter Brune if (snes->ls_monitor) { 141ea630c6eSPeter Brune ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 142ea630c6eSPeter Brune ierr = PetscViewerASCIIPrintf(snes->ls_monitor," Line search: norms[0] = %g, norms[1] = %g, norms[2] = %g alpha %g\n", sqrt(norms[0]),sqrt(norms[1]),sqrt(norms[2]),alpha);CHKERRQ(ierr); 143ea630c6eSPeter Brune ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 144fef7b6d8SPeter Brune } 145*a5892487SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 146*a5892487SPeter Brune ierr = VecAXPY(W, alpha, Y);CHKERRQ(ierr); 147fef7b6d8SPeter Brune if (alpha != 1.0) { 148*a5892487SPeter Brune ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr); 149*a5892487SPeter Brune ierr = VecNorm(G, NORM_2, gnorm);CHKERRQ(ierr); 150fef7b6d8SPeter Brune } else { 151fef7b6d8SPeter Brune *gnorm = PetscSqrtReal(norms[2]); 152fef7b6d8SPeter Brune } 153fef7b6d8SPeter Brune if (alpha == 0.0) *flag = PETSC_FALSE; 154fef7b6d8SPeter Brune else *flag = PETSC_TRUE; 155fef7b6d8SPeter Brune PetscFunctionReturn(0); 156fef7b6d8SPeter Brune } 157fef7b6d8SPeter Brune 158ea630c6eSPeter Brune 159ea630c6eSPeter Brune 160ea630c6eSPeter Brune #undef __FUNCT__ 161ea630c6eSPeter Brune #define __FUNCT__ "SNESLineSearchExactLinear_NCG" 162ea630c6eSPeter 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) 163ea630c6eSPeter Brune { 164ea630c6eSPeter Brune PetscScalar alpha, ptAp; 165ea630c6eSPeter Brune PetscErrorCode ierr; 166ea630c6eSPeter Brune /* SNES_NCG *ncg = (SNES_NCG *) snes->data; */ 167ea630c6eSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 168ea630c6eSPeter Brune 169ea630c6eSPeter Brune PetscFunctionBegin; 170ea630c6eSPeter Brune /* 171ea630c6eSPeter Brune 172ea630c6eSPeter Brune The exact step size for linear CG is just: 173ea630c6eSPeter Brune 174ea630c6eSPeter Brune alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy) 175ea630c6eSPeter Brune 176ea630c6eSPeter Brune */ 177ea630c6eSPeter Brune 178*a5892487SPeter Brune ierr = SNESComputeJacobian(snes, X, &snes->jacobian, &snes->jacobian_pre, &flg);CHKERRQ(ierr); 179ea630c6eSPeter Brune ierr = VecDot(F, F, &alpha);CHKERRQ(ierr); 180ea630c6eSPeter Brune ierr = MatMult(snes->jacobian, Y, W);CHKERRQ(ierr); 181ea630c6eSPeter Brune ierr = VecDot(Y, W, &ptAp);CHKERRQ(ierr); 182ea630c6eSPeter Brune alpha = alpha / ptAp; 183ea630c6eSPeter Brune ierr = PetscPrintf(((PetscObject)snes)->comm, "alpha: %g\n", alpha);CHKERRQ(ierr); 184*a5892487SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 185*a5892487SPeter Brune ierr = VecAXPY(W, alpha, Y);CHKERRQ(ierr); 186*a5892487SPeter Brune ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr); 187*a5892487SPeter Brune ierr = VecNorm(G, NORM_2, gnorm);CHKERRQ(ierr); 188ea630c6eSPeter Brune PetscFunctionReturn(0); 189ea630c6eSPeter Brune } 190ea630c6eSPeter Brune 191fef7b6d8SPeter Brune /* 192fef7b6d8SPeter Brune SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method. 193fef7b6d8SPeter Brune 194fef7b6d8SPeter Brune Input Parameters: 195fef7b6d8SPeter Brune . snes - the SNES context 196fef7b6d8SPeter Brune 197fef7b6d8SPeter Brune Output Parameter: 198fef7b6d8SPeter Brune . outits - number of iterations until termination 199fef7b6d8SPeter Brune 200fef7b6d8SPeter Brune Application Interface Routine: SNESSolve() 201fef7b6d8SPeter Brune */ 202fef7b6d8SPeter Brune #undef __FUNCT__ 203fef7b6d8SPeter Brune #define __FUNCT__ "SNESSolve_NCG" 204fef7b6d8SPeter Brune PetscErrorCode SNESSolve_NCG(SNES snes) 205fef7b6d8SPeter Brune { 206*a5892487SPeter Brune Vec X, dX, lX, F, W, B, G; 207*a5892487SPeter Brune PetscReal fnorm, gnorm, dummyNorm; 208fef7b6d8SPeter Brune PetscScalar dXdot, dXdot_old; 209fef7b6d8SPeter Brune PetscInt maxits, i; 210fef7b6d8SPeter Brune PetscErrorCode ierr; 211fef7b6d8SPeter Brune SNESConvergedReason reason; 212fef7b6d8SPeter Brune PetscBool lsSuccess = PETSC_TRUE; 213fef7b6d8SPeter Brune PetscFunctionBegin; 214fef7b6d8SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 215fef7b6d8SPeter Brune 216fef7b6d8SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 217fef7b6d8SPeter Brune X = snes->vec_sol; /* X^n */ 218169e2be8SPeter Brune dX = snes->work[1]; /* the steepest direction */ 219169e2be8SPeter Brune lX = snes->vec_sol_update; /* the conjugate direction */ 220fef7b6d8SPeter Brune F = snes->vec_func; /* residual vector */ 221*a5892487SPeter Brune W = snes->work[0]; /* work vector and solution output for the line search */ 222302dbbaeSPeter Brune B = snes->vec_rhs; /* the right hand side */ 223*a5892487SPeter Brune G = snes->work[2]; /* the line search result function evaluation */ 224fef7b6d8SPeter Brune 225fef7b6d8SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 226fef7b6d8SPeter Brune snes->iter = 0; 227fef7b6d8SPeter Brune snes->norm = 0.; 228fef7b6d8SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 229fef7b6d8SPeter Brune 230fef7b6d8SPeter Brune /* compute the initial function and preconditioned update delX */ 231fef7b6d8SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 232fef7b6d8SPeter Brune if (snes->domainerror) { 233fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 234fef7b6d8SPeter Brune PetscFunctionReturn(0); 235fef7b6d8SPeter Brune } 236fef7b6d8SPeter Brune /* convergence test */ 237fef7b6d8SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 238fef7b6d8SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 239fef7b6d8SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 240fef7b6d8SPeter Brune snes->norm = fnorm; 241fef7b6d8SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 242fef7b6d8SPeter Brune SNESLogConvHistory(snes,fnorm,0); 243fef7b6d8SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 244fef7b6d8SPeter Brune 245fef7b6d8SPeter Brune /* set parameter for default relative tolerance convergence test */ 246fef7b6d8SPeter Brune snes->ttol = fnorm*snes->rtol; 247fef7b6d8SPeter Brune /* test convergence */ 248fef7b6d8SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 249fef7b6d8SPeter Brune if (snes->reason) PetscFunctionReturn(0); 250fef7b6d8SPeter Brune 251fef7b6d8SPeter Brune /* Call general purpose update function */ 252fef7b6d8SPeter Brune if (snes->ops->update) { 253fef7b6d8SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 254fef7b6d8SPeter Brune } 255fef7b6d8SPeter Brune 256fef7b6d8SPeter Brune /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */ 257fef7b6d8SPeter Brune if (!snes->pc) { 258fef7b6d8SPeter Brune ierr = VecCopy(F, lX);CHKERRQ(ierr); 259fef7b6d8SPeter Brune ierr = VecScale(lX, -1.0);CHKERRQ(ierr); 260fef7b6d8SPeter Brune } else { 261fef7b6d8SPeter Brune ierr = VecCopy(X, lX);CHKERRQ(ierr); 26251e86f29SPeter Brune ierr = SNESSolve(snes->pc, snes->vec_rhs, lX);CHKERRQ(ierr); 263fef7b6d8SPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 264302dbbaeSPeter Brune ierr = VecAXPY(lX, -1.0, X);CHKERRQ(ierr); 26551e86f29SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 266fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_INNER; 267fef7b6d8SPeter Brune PetscFunctionReturn(0); 268fef7b6d8SPeter Brune } 269fef7b6d8SPeter Brune } 270fef7b6d8SPeter Brune ierr = VecDot(lX, lX, &dXdot);CHKERRQ(ierr); 271fef7b6d8SPeter Brune for(i = 1; i < maxits; i++) { 272fef7b6d8SPeter Brune lsSuccess = PETSC_TRUE; 273*a5892487SPeter Brune ierr = (*snes->ops->linesearch)(snes, snes->lsP, X, F, lX, fnorm, 0.0, G, W, &dummyNorm, &gnorm, &lsSuccess);CHKERRQ(ierr); 274fef7b6d8SPeter Brune if (!lsSuccess) { 275fef7b6d8SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 276fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 277fef7b6d8SPeter Brune break; 278fef7b6d8SPeter Brune } 279fef7b6d8SPeter Brune } 280fef7b6d8SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 281fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 282fef7b6d8SPeter Brune break; 283fef7b6d8SPeter Brune } 284fef7b6d8SPeter Brune if (snes->domainerror) { 285fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 286fef7b6d8SPeter Brune PetscFunctionReturn(0); 287fef7b6d8SPeter Brune } 288302dbbaeSPeter Brune 289*a5892487SPeter Brune /* copy over the solution */ 290*a5892487SPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 291*a5892487SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 292*a5892487SPeter Brune fnorm = gnorm; 293*a5892487SPeter Brune 294fef7b6d8SPeter Brune /* Monitor convergence */ 295fef7b6d8SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 296169e2be8SPeter Brune snes->iter = i; 297fef7b6d8SPeter Brune snes->norm = fnorm; 298fef7b6d8SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 299fef7b6d8SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 300fef7b6d8SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 301302dbbaeSPeter Brune 302fef7b6d8SPeter Brune /* Test for convergence */ 303fef7b6d8SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 304fef7b6d8SPeter Brune if (snes->reason) break; 305302dbbaeSPeter Brune 306302dbbaeSPeter Brune /* Call general purpose update function */ 307302dbbaeSPeter Brune if (snes->ops->update) { 308302dbbaeSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 309302dbbaeSPeter Brune } 310ee78dd50SPeter Brune if (!snes->pc) { 311302dbbaeSPeter Brune ierr = VecCopy(F,dX);CHKERRQ(ierr); 312302dbbaeSPeter Brune ierr = VecScale(dX,-1.0);CHKERRQ(ierr); 313302dbbaeSPeter Brune } else { 314302dbbaeSPeter Brune ierr = VecCopy(X,dX);CHKERRQ(ierr); 315ee78dd50SPeter Brune ierr = SNESSolve(snes->pc, snes->vec_rhs, dX);CHKERRQ(ierr); 316302dbbaeSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 31751e86f29SPeter Brune if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 318302dbbaeSPeter Brune snes->reason = SNES_DIVERGED_INNER; 319302dbbaeSPeter Brune PetscFunctionReturn(0); 320302dbbaeSPeter Brune } 321302dbbaeSPeter Brune ierr = VecAXPY(dX,-1.0,X);CHKERRQ(ierr); 322302dbbaeSPeter Brune } 323302dbbaeSPeter Brune 324169e2be8SPeter Brune /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/ 325302dbbaeSPeter Brune dXdot_old = dXdot; 326302dbbaeSPeter Brune ierr = VecDot(dX, dX, &dXdot);CHKERRQ(ierr); 327ee78dd50SPeter Brune #if 0 328421d9b32SPeter Brune if (1) 329169e2be8SPeter Brune ierr = PetscPrintf(PETSC_COMM_WORLD, "beta %e = %e / %e \n", dXdot / dXdot_old, dXdot, dXdot_old);CHKERRQ(ierr); 330ee78dd50SPeter Brune #endif 331302dbbaeSPeter Brune ierr = VecAYPX(lX, dXdot / dXdot_old, dX);CHKERRQ(ierr); 332302dbbaeSPeter Brune /* line search for the proper contribution of lX to the solution */ 333fef7b6d8SPeter Brune } 334fef7b6d8SPeter Brune if (i == maxits) { 335fef7b6d8SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 336fef7b6d8SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 337fef7b6d8SPeter Brune } 338fef7b6d8SPeter Brune PetscFunctionReturn(0); 339fef7b6d8SPeter Brune } 340fef7b6d8SPeter Brune 341fef7b6d8SPeter Brune /*MC 342fef7b6d8SPeter Brune SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems. 343fef7b6d8SPeter Brune 344fef7b6d8SPeter Brune Level: beginner 345fef7b6d8SPeter Brune 346fef7b6d8SPeter Brune Options Database: 347fef7b6d8SPeter Brune + -snes_ls_damping - damping factor to apply to F(x) (used only if -snes_ls is basic or basicnonorms) 348fef7b6d8SPeter Brune - -snes_ls <basic,basicnormnorms,quadratic> 349fef7b6d8SPeter Brune 350fef7b6d8SPeter Brune Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate 351fef7b6d8SPeter Brune gradient method. This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise 352fef7b6d8SPeter Brune chooses the initial search direction as F(x) for the initial guess x. 353fef7b6d8SPeter Brune 354fef7b6d8SPeter Brune 355fef7b6d8SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNGMRES, SNESNQN 356fef7b6d8SPeter Brune M*/ 357fef7b6d8SPeter Brune EXTERN_C_BEGIN 358fef7b6d8SPeter Brune #undef __FUNCT__ 359fef7b6d8SPeter Brune #define __FUNCT__ "SNESCreate_NCG" 360fef7b6d8SPeter Brune PetscErrorCode SNESCreate_NCG(SNES snes) 361fef7b6d8SPeter Brune { 362fef7b6d8SPeter Brune PetscErrorCode ierr; 363ea630c6eSPeter Brune SNES_NCG * neP; 364fef7b6d8SPeter Brune 365fef7b6d8SPeter Brune PetscFunctionBegin; 366fef7b6d8SPeter Brune snes->ops->destroy = SNESDestroy_NCG; 367fef7b6d8SPeter Brune snes->ops->setup = SNESSetUp_NCG; 368fef7b6d8SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NCG; 369fef7b6d8SPeter Brune snes->ops->view = SNESView_NCG; 370fef7b6d8SPeter Brune snes->ops->solve = SNESSolve_NCG; 371fef7b6d8SPeter Brune snes->ops->reset = SNESReset_NCG; 372fef7b6d8SPeter Brune 373fef7b6d8SPeter Brune snes->usesksp = PETSC_FALSE; 374fef7b6d8SPeter Brune snes->usespc = PETSC_TRUE; 375fef7b6d8SPeter Brune 376fef7b6d8SPeter Brune ierr = PetscNewLog(snes, SNES_NCG, &neP);CHKERRQ(ierr); 377fef7b6d8SPeter Brune snes->data = (void*) neP; 378ea630c6eSPeter Brune snes->ops->linesearchquadratic = SNESLineSearchQuadratic_NCG; 379*a5892487SPeter Brune snes->ops->linesearchcubic = SNESLineSearchCubic; 380*a5892487SPeter Brune snes->ops->linesearchno = SNESLineSearchNo; 381*a5892487SPeter Brune snes->ops->linesearchnonorms = SNESLineSearchNoNorms; 382ea630c6eSPeter Brune snes->ops->linesearchtest = SNESLineSearchExactLinear_NCG; 383ea630c6eSPeter Brune snes->lsP = PETSC_NULL; 384ea630c6eSPeter Brune snes->ops->postcheckstep = PETSC_NULL; 385ea630c6eSPeter Brune snes->postcheck = PETSC_NULL; 386ea630c6eSPeter Brune snes->ops->precheckstep = PETSC_NULL; 387ea630c6eSPeter Brune snes->precheck = PETSC_NULL; 388*a5892487SPeter Brune 389*a5892487SPeter Brune ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr); 390*a5892487SPeter Brune 391fef7b6d8SPeter Brune PetscFunctionReturn(0); 392fef7b6d8SPeter Brune } 393fef7b6d8SPeter Brune EXTERN_C_END 394