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