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