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