1 2 #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 3 #include <petscblaslapack.h> 4 5 #undef __FUNCT__ 6 #define __FUNCT__ "SNESMonitorSolution" 7 /*@C 8 SNESMonitorSolution - Monitors progress of the SNES solvers by calling 9 VecView() for the approximate solution at each iteration. 10 11 Collective on SNES 12 13 Input Parameters: 14 + snes - the SNES context 15 . its - iteration number 16 . fgnorm - 2-norm of residual 17 - dummy - either a viewer or PETSC_NULL 18 19 Level: intermediate 20 21 .keywords: SNES, nonlinear, vector, monitor, view 22 23 .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView() 24 @*/ 25 PetscErrorCode SNESMonitorSolution(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 26 { 27 PetscErrorCode ierr; 28 Vec x; 29 PetscViewer viewer = (PetscViewer) dummy; 30 31 PetscFunctionBegin; 32 ierr = SNESGetSolution(snes,&x);CHKERRQ(ierr); 33 if (!viewer) { 34 MPI_Comm comm; 35 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 36 viewer = PETSC_VIEWER_DRAW_(comm); 37 } 38 ierr = VecView(x,viewer);CHKERRQ(ierr); 39 40 PetscFunctionReturn(0); 41 } 42 43 #undef __FUNCT__ 44 #define __FUNCT__ "SNESMonitorResidual" 45 /*@C 46 SNESMonitorResidual - Monitors progress of the SNES solvers by calling 47 VecView() for the residual at each iteration. 48 49 Collective on SNES 50 51 Input Parameters: 52 + snes - the SNES context 53 . its - iteration number 54 . fgnorm - 2-norm of residual 55 - dummy - either a viewer or PETSC_NULL 56 57 Level: intermediate 58 59 .keywords: SNES, nonlinear, vector, monitor, view 60 61 .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView() 62 @*/ 63 PetscErrorCode SNESMonitorResidual(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 64 { 65 PetscErrorCode ierr; 66 Vec x; 67 PetscViewer viewer = (PetscViewer) dummy; 68 69 PetscFunctionBegin; 70 ierr = SNESGetFunction(snes,&x,0,0);CHKERRQ(ierr); 71 if (!viewer) { 72 MPI_Comm comm; 73 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 74 viewer = PETSC_VIEWER_DRAW_(comm); 75 } 76 ierr = VecView(x,viewer);CHKERRQ(ierr); 77 78 PetscFunctionReturn(0); 79 } 80 81 #undef __FUNCT__ 82 #define __FUNCT__ "SNESMonitorSolutionUpdate" 83 /*@C 84 SNESMonitorSolutionUpdate - Monitors progress of the SNES solvers by calling 85 VecView() for the UPDATE to the solution at each iteration. 86 87 Collective on SNES 88 89 Input Parameters: 90 + snes - the SNES context 91 . its - iteration number 92 . fgnorm - 2-norm of residual 93 - dummy - either a viewer or PETSC_NULL 94 95 Level: intermediate 96 97 .keywords: SNES, nonlinear, vector, monitor, view 98 99 .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView() 100 @*/ 101 PetscErrorCode SNESMonitorSolutionUpdate(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 102 { 103 PetscErrorCode ierr; 104 Vec x; 105 PetscViewer viewer = (PetscViewer) dummy; 106 107 PetscFunctionBegin; 108 ierr = SNESGetSolutionUpdate(snes,&x);CHKERRQ(ierr); 109 if (!viewer) { 110 MPI_Comm comm; 111 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 112 viewer = PETSC_VIEWER_DRAW_(comm); 113 } 114 ierr = VecView(x,viewer);CHKERRQ(ierr); 115 116 PetscFunctionReturn(0); 117 } 118 119 #undef __FUNCT__ 120 #define __FUNCT__ "SNESMonitorDefault" 121 /*@C 122 SNESMonitorDefault - Monitors progress of the SNES solvers (default). 123 124 Collective on SNES 125 126 Input Parameters: 127 + snes - the SNES context 128 . its - iteration number 129 . fgnorm - 2-norm of residual 130 - dummy - unused context 131 132 Notes: 133 This routine prints the residual norm at each iteration. 134 135 Level: intermediate 136 137 .keywords: SNES, nonlinear, default, monitor, norm 138 139 .seealso: SNESMonitorSet(), SNESMonitorSolution() 140 @*/ 141 PetscErrorCode SNESMonitorDefault(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 142 { 143 PetscErrorCode ierr; 144 PetscViewer viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm); 145 146 PetscFunctionBegin; 147 ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 148 ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);CHKERRQ(ierr); 149 ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 150 PetscFunctionReturn(0); 151 } 152 153 #undef __FUNCT__ 154 #define __FUNCT__ "SNESMonitorJacUpdateSpectrum" 155 PetscErrorCode SNESMonitorJacUpdateSpectrum(SNES snes,PetscInt it,PetscReal fnorm,void *ctx) { 156 157 Vec X; 158 Mat J,dJ,dJdense; 159 PetscErrorCode ierr; 160 PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*); 161 PetscInt n,i; 162 PetscBLASInt nb,lwork; 163 PetscReal *eigr,*eigi; 164 MatStructure flg = DIFFERENT_NONZERO_PATTERN; 165 PetscScalar *work; 166 PetscScalar *a; 167 168 PetscFunctionBegin; 169 if (it == 0) PetscFunctionReturn(0); 170 /* create the difference between the current update and the current jacobian */ 171 ierr = SNESGetSolution(snes,&X);CHKERRQ(ierr); 172 ierr = SNESGetJacobian(snes,&J,PETSC_NULL,&func,PETSC_NULL);CHKERRQ(ierr); 173 ierr = MatDuplicate(J,MAT_COPY_VALUES,&dJ);CHKERRQ(ierr); 174 ierr = SNESComputeJacobian(snes,X,&dJ,&dJ,&flg);CHKERRQ(ierr); 175 ierr = MatAXPY(dJ,-1.0,J,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 176 /* compute the spectrum directly */ 177 ierr = MatConvert(dJ,MATSEQDENSE,MAT_INITIAL_MATRIX,&dJdense);CHKERRQ(ierr); 178 ierr = MatGetSize(dJ,&n,PETSC_NULL);CHKERRQ(ierr); 179 nb = PetscBLASIntCast(n); 180 lwork = 3*nb; 181 ierr = PetscMalloc(n*sizeof(PetscReal),&eigr);CHKERRQ(ierr); 182 ierr = PetscMalloc(n*sizeof(PetscReal),&eigi);CHKERRQ(ierr); 183 ierr = PetscMalloc(lwork*sizeof(PetscScalar),&work);CHKERRQ(ierr); 184 ierr = MatGetArray(dJdense,&a);CHKERRQ(ierr); 185 #if !defined(PETSC_USE_COMPLEX) 186 { 187 PetscBLASInt lierr; 188 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 189 LAPACKgeev_("N","N",&nb,a,&nb,eigr,eigi,PETSC_NULL,&nb,PETSC_NULL,&nb,work,&lwork,&lierr); 190 if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"geev() error %d",lierr); 191 ierr = PetscFPTrapPop();CHKERRQ(ierr); 192 } 193 #else 194 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not coded for complex"); 195 #endif 196 PetscPrintf(((PetscObject)snes)->comm,"Eigenvalues of J_%d - J_%d:\n",it,it-1);CHKERRQ(ierr); 197 for (i=0;i<n;i++) { 198 PetscPrintf(((PetscObject)snes)->comm,"%5d: %20.5g + %20.5gi\n",i,eigr[i],eigi[i]);CHKERRQ(ierr); 199 } 200 ierr = MatRestoreArray(dJdense,&a);CHKERRQ(ierr); 201 ierr = MatDestroy(&dJ);CHKERRQ(ierr); 202 ierr = MatDestroy(&dJdense);CHKERRQ(ierr); 203 ierr = PetscFree(eigr);CHKERRQ(ierr); 204 ierr = PetscFree(eigi);CHKERRQ(ierr); 205 ierr = PetscFree(work);CHKERRQ(ierr); 206 PetscFunctionReturn(0); 207 } 208 209 #undef __FUNCT__ 210 #define __FUNCT__ "SNESMonitorRange_Private" 211 PetscErrorCode SNESMonitorRange_Private(SNES snes,PetscInt it,PetscReal *per) 212 { 213 PetscErrorCode ierr; 214 Vec resid; 215 PetscReal rmax,pwork; 216 PetscInt i,n,N; 217 PetscScalar *r; 218 219 PetscFunctionBegin; 220 ierr = SNESGetFunction(snes,&resid,0,0);CHKERRQ(ierr); 221 ierr = VecNorm(resid,NORM_INFINITY,&rmax);CHKERRQ(ierr); 222 ierr = VecGetLocalSize(resid,&n);CHKERRQ(ierr); 223 ierr = VecGetSize(resid,&N);CHKERRQ(ierr); 224 ierr = VecGetArray(resid,&r);CHKERRQ(ierr); 225 pwork = 0.0; 226 for (i=0; i<n; i++) { 227 pwork += (PetscAbsScalar(r[i]) > .20*rmax); 228 } 229 ierr = MPI_Allreduce(&pwork,per,1,MPIU_REAL,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr); 230 ierr = VecRestoreArray(resid,&r);CHKERRQ(ierr); 231 *per = *per/N; 232 PetscFunctionReturn(0); 233 } 234 235 #undef __FUNCT__ 236 #define __FUNCT__ "SNESMonitorRange" 237 /*@C 238 SNESMonitorRange - Prints the percentage of residual elements that are more then 10 percent of the maximum value. 239 240 Collective on SNES 241 242 Input Parameters: 243 + snes - iterative context 244 . it - iteration number 245 . rnorm - 2-norm (preconditioned) residual value (may be estimated). 246 - dummy - unused monitor context 247 248 Options Database Key: 249 . -snes_monitor_range - Activates SNESMonitorRange() 250 251 Level: intermediate 252 253 .keywords: SNES, default, monitor, residual 254 255 .seealso: SNESMonitorSet(), SNESMonitorDefault(), SNESMonitorLGCreate() 256 @*/ 257 PetscErrorCode SNESMonitorRange(SNES snes,PetscInt it,PetscReal rnorm,void *dummy) 258 { 259 PetscErrorCode ierr; 260 PetscReal perc,rel; 261 PetscViewer viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm); 262 /* should be in a MonitorRangeContext */ 263 static PetscReal prev; 264 265 PetscFunctionBegin; 266 if (!it) prev = rnorm; 267 ierr = SNESMonitorRange_Private(snes,it,&perc);CHKERRQ(ierr); 268 269 rel = (prev - rnorm)/prev; 270 prev = rnorm; 271 ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 272 ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES preconditioned resid norm %14.12e Percent values above 20 percent of maximum %5.2f relative decrease %5.2e ratio %5.2e \n",it,(double)rnorm,(double)100.0*perc,(double)rel,(double)rel/perc);CHKERRQ(ierr); 273 ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 274 PetscFunctionReturn(0); 275 } 276 277 typedef struct { 278 PetscViewer viewer; 279 PetscReal *history; 280 } SNESMonitorRatioContext; 281 282 #undef __FUNCT__ 283 #define __FUNCT__ "SNESMonitorRatio" 284 /*@C 285 SNESMonitorRatio - Monitors progress of the SNES solvers by printing the ratio 286 of residual norm at each iteration to the previous. 287 288 Collective on SNES 289 290 Input Parameters: 291 + snes - the SNES context 292 . its - iteration number 293 . fgnorm - 2-norm of residual (or gradient) 294 - dummy - context of monitor 295 296 Level: intermediate 297 298 .keywords: SNES, nonlinear, monitor, norm 299 300 .seealso: SNESMonitorSet(), SNESMonitorSolution() 301 @*/ 302 PetscErrorCode SNESMonitorRatio(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 303 { 304 PetscErrorCode ierr; 305 PetscInt len; 306 PetscReal *history; 307 SNESMonitorRatioContext *ctx = (SNESMonitorRatioContext*)dummy; 308 309 PetscFunctionBegin; 310 ierr = SNESGetConvergenceHistory(snes,&history,PETSC_NULL,&len);CHKERRQ(ierr); 311 ierr = PetscViewerASCIIAddTab(ctx->viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 312 if (!its || !history || its > len) { 313 ierr = PetscViewerASCIIPrintf(ctx->viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);CHKERRQ(ierr); 314 } else { 315 PetscReal ratio = fgnorm/history[its-1]; 316 ierr = PetscViewerASCIIPrintf(ctx->viewer,"%3D SNES Function norm %14.12e %14.12e \n",its,(double)fgnorm,(double)ratio);CHKERRQ(ierr); 317 } 318 ierr = PetscViewerASCIISubtractTab(ctx->viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 319 PetscFunctionReturn(0); 320 } 321 322 /* 323 If the we set the history monitor space then we must destroy it 324 */ 325 #undef __FUNCT__ 326 #define __FUNCT__ "SNESMonitorRatioDestroy" 327 PetscErrorCode SNESMonitorRatioDestroy(void **ct) 328 { 329 PetscErrorCode ierr; 330 SNESMonitorRatioContext *ctx = *(SNESMonitorRatioContext**)ct; 331 332 PetscFunctionBegin; 333 ierr = PetscFree(ctx->history);CHKERRQ(ierr); 334 ierr = PetscViewerDestroy(&ctx->viewer);CHKERRQ(ierr); 335 ierr = PetscFree(ctx);CHKERRQ(ierr); 336 PetscFunctionReturn(0); 337 } 338 339 #undef __FUNCT__ 340 #define __FUNCT__ "SNESMonitorSetRatio" 341 /*@C 342 SNESMonitorSetRatio - Sets SNES to use a monitor that prints the 343 ratio of the function norm at each iteration. 344 345 Collective on SNES 346 347 Input Parameters: 348 + snes - the SNES context 349 - viewer - ASCII viewer to print output 350 351 Level: intermediate 352 353 .keywords: SNES, nonlinear, monitor, norm 354 355 .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault() 356 @*/ 357 PetscErrorCode SNESMonitorSetRatio(SNES snes,PetscViewer viewer) 358 { 359 PetscErrorCode ierr; 360 SNESMonitorRatioContext *ctx; 361 PetscReal *history; 362 363 PetscFunctionBegin; 364 if (!viewer) { 365 ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,"stdout",&viewer);CHKERRQ(ierr); 366 ierr = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr); 367 } 368 ierr = PetscNewLog(snes,SNESMonitorRatioContext,&ctx); 369 ierr = SNESGetConvergenceHistory(snes,&history,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 370 if (!history) { 371 ierr = PetscMalloc(100*sizeof(PetscReal),&ctx->history);CHKERRQ(ierr); 372 ierr = SNESSetConvergenceHistory(snes,ctx->history,0,100,PETSC_TRUE);CHKERRQ(ierr); 373 } 374 ctx->viewer = viewer; 375 ierr = SNESMonitorSet(snes,SNESMonitorRatio,ctx,SNESMonitorRatioDestroy);CHKERRQ(ierr); 376 PetscFunctionReturn(0); 377 } 378 379 /* ---------------------------------------------------------------- */ 380 #undef __FUNCT__ 381 #define __FUNCT__ "SNESMonitorDefaultShort" 382 /* 383 Default (short) SNES Monitor, same as SNESMonitorDefault() except 384 it prints fewer digits of the residual as the residual gets smaller. 385 This is because the later digits are meaningless and are often 386 different on different machines; by using this routine different 387 machines will usually generate the same output. 388 */ 389 PetscErrorCode SNESMonitorDefaultShort(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 390 { 391 PetscErrorCode ierr; 392 PetscViewer viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm); 393 394 PetscFunctionBegin; 395 ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 396 if (fgnorm > 1.e-9) { 397 ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %G \n",its,fgnorm);CHKERRQ(ierr); 398 } else if (fgnorm > 1.e-11){ 399 ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %5.3e \n",its,fgnorm);CHKERRQ(ierr); 400 } else { 401 ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm < 1.e-11\n",its);CHKERRQ(ierr); 402 } 403 ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 404 PetscFunctionReturn(0); 405 } 406 /* ---------------------------------------------------------------- */ 407 #undef __FUNCT__ 408 #define __FUNCT__ "SNESDefaultConverged" 409 /*@C 410 SNESDefaultConverged - Convergence test of the solvers for 411 systems of nonlinear equations (default). 412 413 Collective on SNES 414 415 Input Parameters: 416 + snes - the SNES context 417 . it - the iteration (0 indicates before any Newton steps) 418 . xnorm - 2-norm of current iterate 419 . snorm - 2-norm of current step 420 . fnorm - 2-norm of function at current iterate 421 - dummy - unused context 422 423 Output Parameter: 424 . reason - one of 425 $ SNES_CONVERGED_FNORM_ABS - (fnorm < abstol), 426 $ SNES_CONVERGED_SNORM_RELATIVE - (snorm < stol*xnorm), 427 $ SNES_CONVERGED_FNORM_RELATIVE - (fnorm < rtol*fnorm0), 428 $ SNES_DIVERGED_FUNCTION_COUNT - (nfct > maxf), 429 $ SNES_DIVERGED_FNORM_NAN - (fnorm == NaN), 430 $ SNES_CONVERGED_ITERATING - (otherwise), 431 432 where 433 + maxf - maximum number of function evaluations, 434 set with SNESSetTolerances() 435 . nfct - number of function evaluations, 436 . abstol - absolute function norm tolerance, 437 set with SNESSetTolerances() 438 - rtol - relative function norm tolerance, set with SNESSetTolerances() 439 440 Level: intermediate 441 442 .keywords: SNES, nonlinear, default, converged, convergence 443 444 .seealso: SNESSetConvergenceTest() 445 @*/ 446 PetscErrorCode SNESDefaultConverged(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 447 { 448 PetscErrorCode ierr; 449 450 PetscFunctionBegin; 451 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 452 PetscValidPointer(reason,6); 453 454 *reason = SNES_CONVERGED_ITERATING; 455 456 if (!it) { 457 /* set parameter for default relative tolerance convergence test */ 458 snes->ttol = fnorm*snes->rtol; 459 } 460 if (PetscIsInfOrNanReal(fnorm)) { 461 ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr); 462 *reason = SNES_DIVERGED_FNORM_NAN; 463 } else if (fnorm < snes->abstol) { 464 ierr = PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e\n",(double)fnorm,(double)snes->abstol);CHKERRQ(ierr); 465 *reason = SNES_CONVERGED_FNORM_ABS; 466 } else if (snes->nfuncs >= snes->max_funcs) { 467 ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr); 468 *reason = SNES_DIVERGED_FUNCTION_COUNT; 469 } 470 471 if (it && !*reason) { 472 if (fnorm <= snes->ttol) { 473 ierr = PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e (relative tolerance)\n",(double)fnorm,(double)snes->ttol);CHKERRQ(ierr); 474 *reason = SNES_CONVERGED_FNORM_RELATIVE; 475 } else if (snorm < snes->stol*xnorm) { 476 ierr = PetscInfo3(snes,"Converged due to small update length: %14.12e < %14.12e * %14.12e\n",(double)snorm,(double)snes->stol,(double)xnorm);CHKERRQ(ierr); 477 *reason = SNES_CONVERGED_SNORM_RELATIVE; 478 } 479 } 480 PetscFunctionReturn(0); 481 } 482 483 #undef __FUNCT__ 484 #define __FUNCT__ "SNESSkipConverged" 485 /*@C 486 SNESSkipConverged - Convergence test for SNES that NEVER returns as 487 converged, UNLESS the maximum number of iteration have been reached. 488 489 Logically Collective on SNES 490 491 Input Parameters: 492 + snes - the SNES context 493 . it - the iteration (0 indicates before any Newton steps) 494 . xnorm - 2-norm of current iterate 495 . snorm - 2-norm of current step 496 . fnorm - 2-norm of function at current iterate 497 - dummy - unused context 498 499 Output Parameter: 500 . reason - SNES_CONVERGED_ITERATING, SNES_CONVERGED_ITS, or SNES_DIVERGED_FNORM_NAN 501 502 Notes: 503 Convergence is then declared after a fixed number of iterations have been used. 504 505 Level: advanced 506 507 .keywords: SNES, nonlinear, skip, converged, convergence 508 509 .seealso: SNESSetConvergenceTest() 510 @*/ 511 PetscErrorCode SNESSkipConverged(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 512 { 513 PetscErrorCode ierr; 514 515 PetscFunctionBegin; 516 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 517 PetscValidPointer(reason,6); 518 519 *reason = SNES_CONVERGED_ITERATING; 520 521 if (fnorm != fnorm) { 522 ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr); 523 *reason = SNES_DIVERGED_FNORM_NAN; 524 } else if(it == snes->max_its) { 525 *reason = SNES_CONVERGED_ITS; 526 } 527 PetscFunctionReturn(0); 528 } 529 530 #undef __FUNCT__ 531 #define __FUNCT__ "SNESDefaultGetWork" 532 /*@ 533 SNESDefaultGetWork - Gets a number of work vectors. 534 535 Input Parameters: 536 . snes - the SNES context 537 . nw - number of work vectors to allocate 538 539 Level: developer 540 541 Notes: 542 Call this only if no work vectors have been allocated 543 @*/ 544 PetscErrorCode SNESDefaultGetWork(SNES snes,PetscInt nw) 545 { 546 PetscErrorCode ierr; 547 548 PetscFunctionBegin; 549 if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);} 550 snes->nwork = nw; 551 ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 552 ierr = PetscLogObjectParents(snes,nw,snes->work);CHKERRQ(ierr); 553 PetscFunctionReturn(0); 554 } 555