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 SNES_NCG *ncg = (SNES_NCG *)snes->data; 9 10 PetscFunctionBegin; 11 ierr = PetscLineSearchDestroy(&ncg->linesearch);CHKERRQ(ierr); 12 PetscFunctionReturn(0); 13 } 14 15 #define PETSCLINESEARCHNCGLINEAR "linear" 16 17 /* 18 SNESDestroy_NCG - Destroys the private SNES_NCG context that was created with SNESCreate_NCG(). 19 20 Input Parameter: 21 . snes - the SNES context 22 23 Application Interface Routine: SNESDestroy() 24 */ 25 #undef __FUNCT__ 26 #define __FUNCT__ "SNESDestroy_NCG" 27 PetscErrorCode SNESDestroy_NCG(SNES snes) 28 { 29 PetscErrorCode ierr; 30 31 PetscFunctionBegin; 32 ierr = PetscFree(snes->data);CHKERRQ(ierr); 33 PetscFunctionReturn(0); 34 } 35 36 /* 37 SNESSetUp_NCG - Sets up the internal data structures for the later use 38 of the SNESNCG nonlinear solver. 39 40 Input Parameters: 41 + snes - the SNES context 42 - x - the solution vector 43 44 Application Interface Routine: SNESSetUp() 45 */ 46 47 EXTERN_C_BEGIN; 48 extern PetscErrorCode PetscLineSearchCreate_NCGLinear(PetscLineSearch); 49 EXTERN_C_END; 50 51 #undef __FUNCT__ 52 #define __FUNCT__ "SNESSetUp_NCG" 53 PetscErrorCode SNESSetUp_NCG(SNES snes) 54 { 55 PetscErrorCode ierr; 56 SNES_NCG *ncg = (SNES_NCG *)snes->data; 57 const char *optionsprefix; 58 59 PetscFunctionBegin; 60 ierr = SNESDefaultGetWork(snes,2);CHKERRQ(ierr); 61 ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 62 ierr = PetscLineSearchRegisterDynamic(PETSCLINESEARCHNCGLINEAR, PETSC_NULL,"PetscLineSearchCreate_NCGLinear", PetscLineSearchCreate_NCGLinear);CHKERRQ(ierr); 63 ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 64 ierr = PetscLineSearchCreate(((PetscObject)snes)->comm, &ncg->linesearch);CHKERRQ(ierr); 65 ierr = PetscLineSearchSetSNES(ncg->linesearch, snes);CHKERRQ(ierr); 66 if (!snes->pc) { 67 ierr = PetscLineSearchSetType(ncg->linesearch, PETSCLINESEARCHCP);CHKERRQ(ierr); 68 } else { 69 ierr = PetscLineSearchSetType(ncg->linesearch, PETSCLINESEARCHL2);CHKERRQ(ierr); 70 } 71 ierr = PetscLineSearchAppendOptionsPrefix(ncg->linesearch, optionsprefix);CHKERRQ(ierr); 72 ierr = PetscLineSearchSetFromOptions(ncg->linesearch);CHKERRQ(ierr); 73 74 PetscFunctionReturn(0); 75 } 76 /* 77 SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method. 78 79 Input Parameter: 80 . snes - the SNES context 81 82 Application Interface Routine: SNESSetFromOptions() 83 */ 84 #undef __FUNCT__ 85 #define __FUNCT__ "SNESSetFromOptions_NCG" 86 static PetscErrorCode SNESSetFromOptions_NCG(SNES snes) 87 { 88 SNES_NCG *ncg = (SNES_NCG *)snes->data; 89 const char *betas[] = {"FR","PRP","HS", "DY", "CD"}; 90 PetscErrorCode ierr; 91 PetscBool debug, flg; 92 PetscInt indx; 93 94 PetscFunctionBegin; 95 ierr = PetscOptionsHead("SNES NCG options");CHKERRQ(ierr); 96 ierr = PetscOptionsBool("-snes_ncg_monitor", "Monitor NCG iterations", "SNES", ncg->monitor ? PETSC_TRUE: PETSC_FALSE, &debug, PETSC_NULL);CHKERRQ(ierr); 97 ierr = PetscOptionsEList("-snes_ncg_type","Nonlinear CG update used","", betas, 5, "FR",&indx, &flg);CHKERRQ(ierr); 98 if (flg) { 99 ncg->betatype = indx; 100 } 101 if (debug) { 102 ncg->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 103 } 104 ierr = PetscOptionsTail();CHKERRQ(ierr); 105 PetscFunctionReturn(0); 106 } 107 108 /* 109 SNESView_NCG - Prints info from the SNESNCG data structure. 110 111 Input Parameters: 112 + SNES - the SNES context 113 - viewer - visualization context 114 115 Application Interface Routine: SNESView() 116 */ 117 #undef __FUNCT__ 118 #define __FUNCT__ "SNESView_NCG" 119 static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer) 120 { 121 PetscBool iascii; 122 PetscErrorCode ierr; 123 SNES_NCG *ncg = (SNES_NCG *)snes->data; 124 125 PetscFunctionBegin; 126 ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 127 if (iascii) { 128 ierr = PetscViewerASCIIPrintf(viewer," line search type variant: %s\n", ((PetscObject)ncg->linesearch)->type_name);CHKERRQ(ierr); 129 } 130 PetscFunctionReturn(0); 131 } 132 133 #undef __FUNCT__ 134 #define __FUNCT__ "PetscLineSearchApply_NCGLinear" 135 PetscErrorCode PetscLineSearchApply_NCGLinear(PetscLineSearch linesearch) 136 { 137 PetscScalar alpha, ptAp; 138 Vec X, Y, F, W; 139 SNES snes; 140 PetscErrorCode ierr; 141 PetscReal *fnorm, *xnorm, *ynorm; 142 MatStructure flg = DIFFERENT_NONZERO_PATTERN; 143 144 PetscFunctionBegin; 145 146 ierr = PetscLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 147 X = linesearch->vec_sol; 148 W = linesearch->vec_sol_new; 149 F = linesearch->vec_func; 150 Y = linesearch->vec_update; 151 fnorm = &linesearch->fnorm; 152 xnorm = &linesearch->xnorm; 153 ynorm = &linesearch->ynorm; 154 155 /* 156 157 The exact step size for unpreconditioned linear CG is just: 158 alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy) 159 */ 160 ierr = SNESComputeJacobian(snes, X, &snes->jacobian, &snes->jacobian_pre, &flg);CHKERRQ(ierr); 161 ierr = VecDot(F, F, &alpha);CHKERRQ(ierr); 162 ierr = MatMult(snes->jacobian, Y, W);CHKERRQ(ierr); 163 ierr = VecDot(Y, W, &ptAp);CHKERRQ(ierr); 164 alpha = alpha / ptAp; 165 ierr = PetscPrintf(((PetscObject)snes)->comm, "alpha: %G\n", PetscRealPart(alpha));CHKERRQ(ierr); 166 ierr = VecAXPY(X, alpha, Y);CHKERRQ(ierr); 167 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 168 169 ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr); 170 ierr = VecNorm(X, NORM_2, xnorm);CHKERRQ(ierr); 171 ierr = VecNorm(Y, NORM_2, ynorm);CHKERRQ(ierr); 172 173 PetscFunctionReturn(0); 174 } 175 176 EXTERN_C_BEGIN 177 #undef __FUNCT__ 178 #define __FUNCT__ "PetscLineSearchCreate_NCGLinear" 179 180 PetscErrorCode PetscLineSearchCreate_NCGLinear(PetscLineSearch linesearch) 181 { 182 PetscFunctionBegin; 183 linesearch->ops->apply = PetscLineSearchApply_NCGLinear; 184 linesearch->ops->destroy = PETSC_NULL; 185 linesearch->ops->setfromoptions = PETSC_NULL; 186 linesearch->ops->reset = PETSC_NULL; 187 linesearch->ops->view = PETSC_NULL; 188 linesearch->ops->setup = PETSC_NULL; 189 PetscFunctionReturn(0); 190 } 191 EXTERN_C_END 192 193 #undef __FUNCT__ 194 #define __FUNCT__ "ComputeYtJtF_Private" 195 /* 196 197 Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian. 198 199 */ 200 PetscErrorCode ComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf) { 201 PetscErrorCode ierr; 202 PetscScalar ftf, ftg, fty, h; 203 PetscFunctionBegin; 204 ierr = VecDot(F, F, &ftf);CHKERRQ(ierr); 205 ierr = VecDot(F, Y, &fty);CHKERRQ(ierr); 206 h = 1e-5*fty / fty; 207 ierr = VecCopy(X, W);CHKERRQ(ierr); 208 ierr = VecAXPY(W, -h, Y);CHKERRQ(ierr); /* this is arbitrary */ 209 ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr); 210 ierr = VecDot(G, F, &ftg);CHKERRQ(ierr); 211 *ytJtf = (ftg - ftf) / h; 212 PetscFunctionReturn(0); 213 } 214 215 /* 216 SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method. 217 218 Input Parameters: 219 . snes - the SNES context 220 221 Output Parameter: 222 . outits - number of iterations until termination 223 224 Application Interface Routine: SNESSolve() 225 */ 226 #undef __FUNCT__ 227 #define __FUNCT__ "SNESSolve_NCG" 228 PetscErrorCode SNESSolve_NCG(SNES snes) 229 { 230 SNES_NCG *ncg = (SNES_NCG *)snes->data; 231 Vec X, dX, lX, F, B, Fold; 232 PetscReal fnorm, ynorm, xnorm, beta = 0.0; 233 PetscScalar dXdotF, dXolddotFold, dXdotFold, lXdotF, lXdotFold; 234 PetscInt maxits, i; 235 PetscErrorCode ierr; 236 SNESConvergedReason reason; 237 PetscBool lsSuccess = PETSC_TRUE; 238 239 PetscFunctionBegin; 240 snes->reason = SNES_CONVERGED_ITERATING; 241 242 maxits = snes->max_its; /* maximum number of iterations */ 243 X = snes->vec_sol; /* X^n */ 244 Fold = snes->work[0]; /* The previous iterate of X */ 245 dX = snes->work[1]; /* the preconditioned direction */ 246 lX = snes->vec_sol_update; /* the conjugate direction */ 247 F = snes->vec_func; /* residual vector */ 248 B = snes->vec_rhs; /* the right hand side */ 249 250 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 251 snes->iter = 0; 252 snes->norm = 0.; 253 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 254 255 /* compute the initial function and preconditioned update dX */ 256 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 257 if (snes->domainerror) { 258 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 259 PetscFunctionReturn(0); 260 } 261 /* convergence test */ 262 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 263 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 264 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 265 snes->norm = fnorm; 266 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 267 SNESLogConvHistory(snes,fnorm,0); 268 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 269 270 /* set parameter for default relative tolerance convergence test */ 271 snes->ttol = fnorm*snes->rtol; 272 /* test convergence */ 273 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 274 if (snes->reason) PetscFunctionReturn(0); 275 276 /* Call general purpose update function */ 277 if (snes->ops->update) { 278 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 279 } 280 281 /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */ 282 283 if (snes->pc) { 284 ierr = VecCopy(X, dX);CHKERRQ(ierr); 285 ierr = SNESSolve(snes->pc, B, dX);CHKERRQ(ierr); 286 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 287 if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 288 snes->reason = SNES_DIVERGED_INNER; 289 PetscFunctionReturn(0); 290 } 291 ierr = VecAYPX(dX,-1.0,X);CHKERRQ(ierr); 292 } else { 293 ierr = VecCopy(F, dX);CHKERRQ(ierr); 294 } 295 ierr = VecCopy(dX, lX);CHKERRQ(ierr); 296 ierr = VecDot(F, dX, &dXdotF);CHKERRQ(ierr); 297 /* 298 } else { 299 ierr = ComputeYtJtF_Private(snes, X, F, dX, W, G, &dXdotF);CHKERRQ(ierr); 300 } 301 */ 302 for(i = 1; i < maxits + 1; i++) { 303 lsSuccess = PETSC_TRUE; 304 /* some update types require the old update direction or conjugate direction */ 305 if (ncg->betatype == 1 || ncg->betatype == 2 || ncg->betatype == 3 || ncg->betatype == 4) { /* prp, hs, dy, cd need dXold*/ 306 ierr = VecCopy(F, Fold);CHKERRQ(ierr); 307 } 308 ierr = PetscLineSearchApply(ncg->linesearch, X, F, &fnorm, lX);CHKERRQ(ierr); 309 ierr = PetscLineSearchGetSuccess(ncg->linesearch, &lsSuccess);CHKERRQ(ierr); 310 if (!lsSuccess) { 311 if (++snes->numFailures >= snes->maxFailures) { 312 snes->reason = SNES_DIVERGED_LINE_SEARCH; 313 PetscFunctionReturn(0); 314 } 315 } 316 if (snes->nfuncs >= snes->max_funcs) { 317 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 318 PetscFunctionReturn(0); 319 } 320 if (snes->domainerror) { 321 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 322 PetscFunctionReturn(0); 323 } 324 ierr = PetscLineSearchGetNorms(ncg->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 325 /* Monitor convergence */ 326 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 327 snes->iter = i; 328 snes->norm = fnorm; 329 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 330 SNESLogConvHistory(snes,snes->norm,0); 331 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 332 333 /* Test for convergence */ 334 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 335 if (snes->reason) PetscFunctionReturn(0); 336 337 /* Call general purpose update function */ 338 if (snes->ops->update) { 339 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 340 } 341 if (snes->pc) { 342 ierr = VecCopy(X,dX);CHKERRQ(ierr); 343 ierr = SNESSolve(snes->pc, B, dX);CHKERRQ(ierr); 344 ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 345 if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) { 346 snes->reason = SNES_DIVERGED_INNER; 347 PetscFunctionReturn(0); 348 } 349 ierr = VecAYPX(dX,-1.0,X);CHKERRQ(ierr); 350 } else { 351 ierr = VecCopy(F, dX);CHKERRQ(ierr); 352 } 353 354 /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/ 355 switch(ncg->betatype) { 356 case 0: /* Fletcher-Reeves */ 357 dXolddotFold = dXdotF; 358 ierr = VecDot(dX, dX, &dXdotF);CHKERRQ(ierr); 359 beta = PetscRealPart(dXdotF / dXolddotFold); 360 break; 361 case 1: /* Polak-Ribiere-Poylak */ 362 dXolddotFold = dXdotF; 363 ierr = VecDot(F, dX, &dXdotF);CHKERRQ(ierr); 364 ierr = VecDot(Fold, dX, &dXdotFold);CHKERRQ(ierr); 365 beta = PetscRealPart(((dXdotF - dXdotFold) / dXolddotFold)); 366 if (beta < 0.0) beta = 0.0; /* restart */ 367 break; 368 case 2: /* Hestenes-Stiefel */ 369 ierr = VecDot(dX, F, &dXdotF);CHKERRQ(ierr); 370 ierr = VecDot(dX, Fold, &dXdotFold);CHKERRQ(ierr); 371 ierr = VecDot(lX, F, &lXdotF);CHKERRQ(ierr); 372 ierr = VecDot(lX, Fold, &lXdotFold);CHKERRQ(ierr); 373 beta = PetscRealPart((dXdotF - dXdotFold) / (lXdotF - lXdotFold)); 374 break; 375 case 3: /* Dai-Yuan */ 376 ierr = VecDot(dX, F, &dXdotF);CHKERRQ(ierr); 377 ierr = VecDot(lX, F, &lXdotF);CHKERRQ(ierr); 378 ierr = VecDot(lX, Fold, &lXdotFold);CHKERRQ(ierr); 379 beta = PetscRealPart(dXdotF / (lXdotFold - lXdotF));CHKERRQ(ierr); 380 break; 381 case 4: /* Conjugate Descent */ 382 ierr = VecDot(dX, F, &dXdotF);CHKERRQ(ierr); 383 ierr = VecDot(lX, Fold, &lXdotFold);CHKERRQ(ierr); 384 beta = PetscRealPart(dXdotF / lXdotFold);CHKERRQ(ierr); 385 break; 386 } 387 if (ncg->monitor) { 388 ierr = PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", beta);CHKERRQ(ierr); 389 } 390 ierr = VecAYPX(lX, beta, dX);CHKERRQ(ierr); 391 } 392 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 393 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 394 PetscFunctionReturn(0); 395 } 396 397 /*MC 398 SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems. 399 400 Level: beginner 401 402 Options Database: 403 + -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter. 404 . -snes_ls <basic,basicnormnorms,quadratic,critical,test> - Line search type. 405 - -snes_ncg_monitor - Print relevant information about the ncg iteration. 406 407 Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate 408 gradient method. This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise 409 chooses the initial search direction as F(x) for the initial guess x. 410 411 412 .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNGMRES, SNESNQN 413 M*/ 414 EXTERN_C_BEGIN 415 #undef __FUNCT__ 416 #define __FUNCT__ "SNESCreate_NCG" 417 PetscErrorCode SNESCreate_NCG(SNES snes) 418 { 419 PetscErrorCode ierr; 420 SNES_NCG * neP; 421 422 PetscFunctionBegin; 423 snes->ops->destroy = SNESDestroy_NCG; 424 snes->ops->setup = SNESSetUp_NCG; 425 snes->ops->setfromoptions = SNESSetFromOptions_NCG; 426 snes->ops->view = SNESView_NCG; 427 snes->ops->solve = SNESSolve_NCG; 428 snes->ops->reset = SNESReset_NCG; 429 430 snes->usesksp = PETSC_FALSE; 431 snes->usespc = PETSC_TRUE; 432 433 snes->max_funcs = 30000; 434 snes->max_its = 10000; 435 436 ierr = PetscNewLog(snes, SNES_NCG, &neP);CHKERRQ(ierr); 437 snes->data = (void*) neP; 438 neP->monitor = PETSC_NULL; 439 neP->betatype = 1; 440 PetscFunctionReturn(0); 441 } 442 EXTERN_C_END 443