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