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