1 #include <petscconvest.h> /*I "petscconvest.h" I*/ 2 #include <petscdmplex.h> 3 #include <petscds.h> 4 5 #include <petsc/private/petscconvestimpl.h> 6 7 static PetscErrorCode zero_private(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 8 { 9 PetscInt c; 10 for (c = 0; c < Nc; ++c) u[c] = 0.0; 11 return PETSC_SUCCESS; 12 } 13 14 /*@ 15 PetscConvEstDestroy - Destroys a PetscConvEst object 16 17 Collective 18 19 Input Parameter: 20 . ce - The PetscConvEst object 21 22 Level: beginner 23 24 .seealso: `PetscConvEstCreate()`, `PetscConvEstGetConvRate()` 25 @*/ 26 PetscErrorCode PetscConvEstDestroy(PetscConvEst *ce) 27 { 28 PetscFunctionBegin; 29 if (!*ce) PetscFunctionReturn(PETSC_SUCCESS); 30 PetscValidHeaderSpecific((*ce), PETSC_OBJECT_CLASSID, 1); 31 if (--((PetscObject)(*ce))->refct > 0) { 32 *ce = NULL; 33 PetscFunctionReturn(PETSC_SUCCESS); 34 } 35 PetscCall(PetscFree3((*ce)->initGuess, (*ce)->exactSol, (*ce)->ctxs)); 36 PetscCall(PetscFree2((*ce)->dofs, (*ce)->errors)); 37 PetscCall(PetscHeaderDestroy(ce)); 38 PetscFunctionReturn(PETSC_SUCCESS); 39 } 40 41 /*@ 42 PetscConvEstSetFromOptions - Sets a PetscConvEst object from options 43 44 Collective 45 46 Input Parameters: 47 . ce - The PetscConvEst object 48 49 Level: beginner 50 51 .seealso: `PetscConvEstCreate()`, `PetscConvEstGetConvRate()` 52 @*/ 53 PetscErrorCode PetscConvEstSetFromOptions(PetscConvEst ce) 54 { 55 PetscFunctionBegin; 56 PetscOptionsBegin(PetscObjectComm((PetscObject)ce), "", "Convergence Estimator Options", "PetscConvEst"); 57 PetscCall(PetscOptionsInt("-convest_num_refine", "The number of refinements for the convergence check", "PetscConvEst", ce->Nr, &ce->Nr, NULL)); 58 PetscCall(PetscOptionsReal("-convest_refine_factor", "The increase in resolution in each dimension", "PetscConvEst", ce->r, &ce->r, NULL)); 59 PetscCall(PetscOptionsBool("-convest_monitor", "Monitor the error for each convergence check", "PetscConvEst", ce->monitor, &ce->monitor, NULL)); 60 PetscCall(PetscOptionsBool("-convest_no_refine", "Debugging flag to run on the same mesh each time", "PetscConvEst", ce->noRefine, &ce->noRefine, NULL)); 61 PetscOptionsEnd(); 62 PetscFunctionReturn(PETSC_SUCCESS); 63 } 64 65 /*@ 66 PetscConvEstView - Views a PetscConvEst object 67 68 Collective 69 70 Input Parameters: 71 + ce - The PetscConvEst object 72 - viewer - The PetscViewer object 73 74 Level: beginner 75 76 .seealso: `PetscConvEstCreate()`, `PetscConvEstGetConvRate()` 77 @*/ 78 PetscErrorCode PetscConvEstView(PetscConvEst ce, PetscViewer viewer) 79 { 80 PetscFunctionBegin; 81 PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)ce, viewer)); 82 PetscCall(PetscViewerASCIIPrintf(viewer, "ConvEst with %" PetscInt_FMT " levels\n", ce->Nr + 1)); 83 PetscFunctionReturn(PETSC_SUCCESS); 84 } 85 86 /*@ 87 PetscConvEstGetSolver - Gets the solver used to produce discrete solutions 88 89 Not collective 90 91 Input Parameter: 92 . ce - The PetscConvEst object 93 94 Output Parameter: 95 . solver - The solver 96 97 Level: intermediate 98 99 .seealso: `PetscConvEstSetSolver()`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()` 100 @*/ 101 PetscErrorCode PetscConvEstGetSolver(PetscConvEst ce, PetscObject *solver) 102 { 103 PetscFunctionBegin; 104 PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 105 PetscValidPointer(solver, 2); 106 *solver = ce->solver; 107 PetscFunctionReturn(PETSC_SUCCESS); 108 } 109 110 /*@ 111 PetscConvEstSetSolver - Sets the solver used to produce discrete solutions 112 113 Not collective 114 115 Input Parameters: 116 + ce - The PetscConvEst object 117 - solver - The solver 118 119 Level: intermediate 120 121 Note: The solver MUST have an attached DM/DS, so that we know the exact solution 122 123 .seealso: `PetscConvEstGetSNES()`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()` 124 @*/ 125 PetscErrorCode PetscConvEstSetSolver(PetscConvEst ce, PetscObject solver) 126 { 127 PetscFunctionBegin; 128 PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 129 PetscValidHeader(solver, 2); 130 ce->solver = solver; 131 PetscUseTypeMethod(ce, setsolver, solver); 132 PetscFunctionReturn(PETSC_SUCCESS); 133 } 134 135 /*@ 136 PetscConvEstSetUp - After the solver is specified, we create structures for estimating convergence 137 138 Collective 139 140 Input Parameters: 141 . ce - The PetscConvEst object 142 143 Level: beginner 144 145 .seealso: `PetscConvEstCreate()`, `PetscConvEstGetConvRate()` 146 @*/ 147 PetscErrorCode PetscConvEstSetUp(PetscConvEst ce) 148 { 149 PetscInt Nf, f, Nds, s; 150 151 PetscFunctionBegin; 152 PetscCall(DMGetNumFields(ce->idm, &Nf)); 153 ce->Nf = PetscMax(Nf, 1); 154 PetscCall(PetscMalloc2((ce->Nr + 1) * ce->Nf, &ce->dofs, (ce->Nr + 1) * ce->Nf, &ce->errors)); 155 PetscCall(PetscCalloc3(ce->Nf, &ce->initGuess, ce->Nf, &ce->exactSol, ce->Nf, &ce->ctxs)); 156 for (f = 0; f < Nf; ++f) ce->initGuess[f] = zero_private; 157 PetscCall(DMGetNumDS(ce->idm, &Nds)); 158 for (s = 0; s < Nds; ++s) { 159 PetscDS ds; 160 DMLabel label; 161 IS fieldIS; 162 const PetscInt *fields; 163 PetscInt dsNf; 164 165 PetscCall(DMGetRegionNumDS(ce->idm, s, &label, &fieldIS, &ds)); 166 PetscCall(PetscDSGetNumFields(ds, &dsNf)); 167 if (fieldIS) PetscCall(ISGetIndices(fieldIS, &fields)); 168 for (f = 0; f < dsNf; ++f) { 169 const PetscInt field = fields[f]; 170 PetscCall(PetscDSGetExactSolution(ds, field, &ce->exactSol[field], &ce->ctxs[field])); 171 } 172 if (fieldIS) PetscCall(ISRestoreIndices(fieldIS, &fields)); 173 } 174 for (f = 0; f < Nf; ++f) PetscCheck(ce->exactSol[f], PetscObjectComm((PetscObject)ce), PETSC_ERR_ARG_WRONG, "DS must contain exact solution functions in order to estimate convergence, missing for field %" PetscInt_FMT, f); 175 PetscFunctionReturn(PETSC_SUCCESS); 176 } 177 178 PetscErrorCode PetscConvEstComputeInitialGuess(PetscConvEst ce, PetscInt r, DM dm, Vec u) 179 { 180 PetscFunctionBegin; 181 PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 182 if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3); 183 PetscValidHeaderSpecific(u, VEC_CLASSID, 4); 184 PetscUseTypeMethod(ce, initguess, r, dm, u); 185 PetscFunctionReturn(PETSC_SUCCESS); 186 } 187 188 PetscErrorCode PetscConvEstComputeError(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[]) 189 { 190 PetscFunctionBegin; 191 PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 192 if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3); 193 PetscValidHeaderSpecific(u, VEC_CLASSID, 4); 194 PetscValidRealPointer(errors, 5); 195 PetscUseTypeMethod(ce, computeerror, r, dm, u, errors); 196 PetscFunctionReturn(PETSC_SUCCESS); 197 } 198 199 /*@ 200 PetscConvEstMonitorDefault - Monitors the convergence estimation loop 201 202 Collective 203 204 Input Parameters: 205 + ce - The `PetscConvEst` object 206 - r - The refinement level 207 208 Options database keys: 209 . -convest_monitor - Activate the monitor 210 211 Level: intermediate 212 213 .seealso: `PetscConvEst`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`, `SNESSolve()`, `TSSolve()` 214 @*/ 215 PetscErrorCode PetscConvEstMonitorDefault(PetscConvEst ce, PetscInt r) 216 { 217 MPI_Comm comm; 218 PetscInt f; 219 220 PetscFunctionBegin; 221 if (ce->monitor) { 222 PetscInt *dofs = &ce->dofs[r * ce->Nf]; 223 PetscReal *errors = &ce->errors[r * ce->Nf]; 224 225 PetscCall(PetscObjectGetComm((PetscObject)ce, &comm)); 226 PetscCall(PetscPrintf(comm, "N: ")); 227 if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "[")); 228 for (f = 0; f < ce->Nf; ++f) { 229 if (f > 0) PetscCall(PetscPrintf(comm, ", ")); 230 PetscCall(PetscPrintf(comm, "%7" PetscInt_FMT, dofs[f])); 231 } 232 if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "]")); 233 PetscCall(PetscPrintf(comm, " ")); 234 PetscCall(PetscPrintf(comm, "L_2 Error: ")); 235 if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "[")); 236 for (f = 0; f < ce->Nf; ++f) { 237 if (f > 0) PetscCall(PetscPrintf(comm, ", ")); 238 if (errors[f] < 1.0e-11) PetscCall(PetscPrintf(comm, "< 1e-11")); 239 else PetscCall(PetscPrintf(comm, "%g", (double)errors[f])); 240 } 241 if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "]")); 242 PetscCall(PetscPrintf(comm, "\n")); 243 } 244 PetscFunctionReturn(PETSC_SUCCESS); 245 } 246 247 static PetscErrorCode PetscConvEstSetSNES_Private(PetscConvEst ce, PetscObject solver) 248 { 249 PetscClassId id; 250 251 PetscFunctionBegin; 252 PetscCall(PetscObjectGetClassId(ce->solver, &id)); 253 PetscCheck(id == SNES_CLASSID, PetscObjectComm((PetscObject)ce), PETSC_ERR_ARG_WRONG, "Solver was not a SNES"); 254 PetscCall(SNESGetDM((SNES)ce->solver, &ce->idm)); 255 PetscFunctionReturn(PETSC_SUCCESS); 256 } 257 258 static PetscErrorCode PetscConvEstInitGuessSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u) 259 { 260 PetscFunctionBegin; 261 PetscCall(DMProjectFunction(dm, 0.0, ce->initGuess, ce->ctxs, INSERT_VALUES, u)); 262 PetscFunctionReturn(PETSC_SUCCESS); 263 } 264 265 static PetscErrorCode PetscConvEstComputeErrorSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[]) 266 { 267 PetscFunctionBegin; 268 PetscCall(DMComputeL2FieldDiff(dm, 0.0, ce->exactSol, ce->ctxs, u, errors)); 269 PetscFunctionReturn(PETSC_SUCCESS); 270 } 271 272 static PetscErrorCode PetscConvEstSetJacobianNullSpace_Private(PetscConvEst ce, SNES snes) 273 { 274 DM dm; 275 PetscInt f; 276 277 PetscFunctionBegin; 278 PetscCall(SNESGetDM(snes, &dm)); 279 for (f = 0; f < ce->Nf; ++f) { 280 PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *); 281 282 PetscCall(DMGetNullSpaceConstructor(dm, f, &nspconstr)); 283 if (nspconstr) { 284 MatNullSpace nullsp; 285 Mat J; 286 287 PetscCall((*nspconstr)(dm, f, f, &nullsp)); 288 PetscCall(SNESSetUp(snes)); 289 PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL)); 290 PetscCall(MatSetNullSpace(J, nullsp)); 291 PetscCall(MatNullSpaceDestroy(&nullsp)); 292 break; 293 } 294 } 295 PetscFunctionReturn(PETSC_SUCCESS); 296 } 297 298 static PetscErrorCode PetscConvEstGetConvRateSNES_Private(PetscConvEst ce, PetscReal alpha[]) 299 { 300 SNES snes = (SNES)ce->solver; 301 DM *dm; 302 PetscObject disc; 303 PetscReal *x, *y, slope, intercept; 304 PetscInt Nr = ce->Nr, r, f, dim, oldlevel, oldnlev; 305 void *ctx; 306 307 PetscFunctionBegin; 308 PetscCheck(ce->r == 2.0, PetscObjectComm((PetscObject)ce), PETSC_ERR_SUP, "Only refinement factor 2 is currently supported (not %g)", (double)ce->r); 309 PetscCall(DMGetDimension(ce->idm, &dim)); 310 PetscCall(DMGetApplicationContext(ce->idm, &ctx)); 311 PetscCall(DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE)); 312 PetscCall(DMGetRefineLevel(ce->idm, &oldlevel)); 313 PetscCall(PetscMalloc1((Nr + 1), &dm)); 314 /* Loop over meshes */ 315 dm[0] = ce->idm; 316 for (r = 0; r <= Nr; ++r) { 317 Vec u; 318 #if defined(PETSC_USE_LOG) 319 PetscLogStage stage; 320 #endif 321 char stageName[PETSC_MAX_PATH_LEN]; 322 const char *dmname, *uname; 323 324 PetscCall(PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN - 1, "ConvEst Refinement Level %" PetscInt_FMT, r)); 325 #if defined(PETSC_USE_LOG) 326 PetscCall(PetscLogStageGetId(stageName, &stage)); 327 if (stage < 0) PetscCall(PetscLogStageRegister(stageName, &stage)); 328 #endif 329 PetscCall(PetscLogStagePush(stage)); 330 if (r > 0) { 331 if (!ce->noRefine) { 332 PetscCall(DMRefine(dm[r - 1], MPI_COMM_NULL, &dm[r])); 333 PetscCall(DMSetCoarseDM(dm[r], dm[r - 1])); 334 } else { 335 DM cdm, rcdm; 336 337 PetscCall(DMClone(dm[r - 1], &dm[r])); 338 PetscCall(DMCopyDisc(dm[r - 1], dm[r])); 339 PetscCall(DMGetCoordinateDM(dm[r - 1], &cdm)); 340 PetscCall(DMGetCoordinateDM(dm[r], &rcdm)); 341 PetscCall(DMCopyDisc(cdm, rcdm)); 342 } 343 PetscCall(DMCopyTransform(ce->idm, dm[r])); 344 PetscCall(PetscObjectGetName((PetscObject)dm[r - 1], &dmname)); 345 PetscCall(PetscObjectSetName((PetscObject)dm[r], dmname)); 346 for (f = 0; f < ce->Nf; ++f) { 347 PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *); 348 349 PetscCall(DMGetNullSpaceConstructor(dm[r - 1], f, &nspconstr)); 350 PetscCall(DMSetNullSpaceConstructor(dm[r], f, nspconstr)); 351 } 352 } 353 PetscCall(DMViewFromOptions(dm[r], NULL, "-conv_dm_view")); 354 /* Create solution */ 355 PetscCall(DMCreateGlobalVector(dm[r], &u)); 356 PetscCall(DMGetField(dm[r], 0, NULL, &disc)); 357 PetscCall(PetscObjectGetName(disc, &uname)); 358 PetscCall(PetscObjectSetName((PetscObject)u, uname)); 359 /* Setup solver */ 360 PetscCall(SNESReset(snes)); 361 PetscCall(SNESSetDM(snes, dm[r])); 362 PetscCall(DMPlexSetSNESLocalFEM(dm[r], ctx, ctx, ctx)); 363 PetscCall(SNESSetFromOptions(snes)); 364 /* Set nullspace for Jacobian */ 365 PetscCall(PetscConvEstSetJacobianNullSpace_Private(ce, snes)); 366 /* Create initial guess */ 367 PetscCall(PetscConvEstComputeInitialGuess(ce, r, dm[r], u)); 368 PetscCall(SNESSolve(snes, NULL, u)); 369 PetscCall(PetscLogEventBegin(ce->event, ce, 0, 0, 0)); 370 PetscCall(PetscConvEstComputeError(ce, r, dm[r], u, &ce->errors[r * ce->Nf])); 371 PetscCall(PetscLogEventEnd(ce->event, ce, 0, 0, 0)); 372 for (f = 0; f < ce->Nf; ++f) { 373 PetscSection s, fs; 374 PetscInt lsize; 375 376 /* Could use DMGetOutputDM() to add in Dirichlet dofs */ 377 PetscCall(DMGetLocalSection(dm[r], &s)); 378 PetscCall(PetscSectionGetField(s, f, &fs)); 379 PetscCall(PetscSectionGetConstrainedStorageSize(fs, &lsize)); 380 PetscCallMPI(MPI_Allreduce(&lsize, &ce->dofs[r * ce->Nf + f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes))); 381 PetscCall(PetscLogEventSetDof(ce->event, f, ce->dofs[r * ce->Nf + f])); 382 PetscCall(PetscLogEventSetError(ce->event, f, ce->errors[r * ce->Nf + f])); 383 } 384 /* Monitor */ 385 PetscCall(PetscConvEstMonitorDefault(ce, r)); 386 if (!r) { 387 /* PCReset() does not wipe out the level structure */ 388 KSP ksp; 389 PC pc; 390 391 PetscCall(SNESGetKSP(snes, &ksp)); 392 PetscCall(KSPGetPC(ksp, &pc)); 393 PetscCall(PCMGGetLevels(pc, &oldnlev)); 394 } 395 /* Cleanup */ 396 PetscCall(VecDestroy(&u)); 397 PetscCall(PetscLogStagePop()); 398 } 399 for (r = 1; r <= Nr; ++r) PetscCall(DMDestroy(&dm[r])); 400 /* Fit convergence rate */ 401 PetscCall(PetscMalloc2(Nr + 1, &x, Nr + 1, &y)); 402 for (f = 0; f < ce->Nf; ++f) { 403 for (r = 0; r <= Nr; ++r) { 404 x[r] = PetscLog10Real(ce->dofs[r * ce->Nf + f]); 405 y[r] = PetscLog10Real(ce->errors[r * ce->Nf + f]); 406 } 407 PetscCall(PetscLinearRegression(Nr + 1, x, y, &slope, &intercept)); 408 /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */ 409 alpha[f] = -slope * dim; 410 } 411 PetscCall(PetscFree2(x, y)); 412 PetscCall(PetscFree(dm)); 413 /* Restore solver */ 414 PetscCall(SNESReset(snes)); 415 { 416 /* PCReset() does not wipe out the level structure */ 417 KSP ksp; 418 PC pc; 419 420 PetscCall(SNESGetKSP(snes, &ksp)); 421 PetscCall(KSPGetPC(ksp, &pc)); 422 PetscCall(PCMGSetLevels(pc, oldnlev, NULL)); 423 PetscCall(DMSetRefineLevel(ce->idm, oldlevel)); /* The damn DMCoarsen() calls in PCMG can reset this */ 424 } 425 PetscCall(SNESSetDM(snes, ce->idm)); 426 PetscCall(DMPlexSetSNESLocalFEM(ce->idm, ctx, ctx, ctx)); 427 PetscCall(SNESSetFromOptions(snes)); 428 PetscCall(PetscConvEstSetJacobianNullSpace_Private(ce, snes)); 429 PetscFunctionReturn(PETSC_SUCCESS); 430 } 431 432 /*@ 433 PetscConvEstGetConvRate - Returns an estimate of the convergence rate for the discretization 434 435 Not collective 436 437 Input Parameter: 438 . ce - The `PetscConvEst` object 439 440 Output Parameter: 441 . alpha - The convergence rate for each field 442 443 Options database keys: 444 + -snes_convergence_estimate - Execute convergence estimation inside `SNESSolve()` and print out the rate 445 - -ts_convergence_estimate - Execute convergence estimation inside `TSSolve()` and print out the rate 446 447 Notes: 448 The convergence rate alpha is defined by 449 $ || u_\Delta - u_exact || < C \Delta^alpha 450 where u_\Delta is the discrete solution, and Delta is a measure of the discretization size. We usually use h for the 451 spatial resolution and \Delta t for the temporal resolution. 452 453 We solve a series of problems using increasing resolution (refined meshes or decreased timesteps), calculate an error 454 based upon the exact solution in the DS, and then fit the result to our model above using linear regression. 455 456 Level: intermediate 457 458 .seealso: `PetscConvEstSetSolver()`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`, `SNESSolve()`, `TSSolve()` 459 @*/ 460 PetscErrorCode PetscConvEstGetConvRate(PetscConvEst ce, PetscReal alpha[]) 461 { 462 PetscInt f; 463 464 PetscFunctionBegin; 465 if (ce->event < 0) PetscCall(PetscLogEventRegister("ConvEst Error", PETSC_OBJECT_CLASSID, &ce->event)); 466 for (f = 0; f < ce->Nf; ++f) alpha[f] = 0.0; 467 PetscUseTypeMethod(ce, getconvrate, alpha); 468 PetscFunctionReturn(PETSC_SUCCESS); 469 } 470 471 /*@ 472 PetscConvEstRateView - Displays the convergence rate to a viewer 473 474 Collective 475 476 Parameter: 477 + snes - iterative context obtained from `SNESCreate()` 478 . alpha - the convergence rate for each field 479 - viewer - the viewer to display the reason 480 481 Options Database Keys: 482 . -snes_convergence_estimate - print the convergence rate 483 484 Level: developer 485 486 .seealso: `PetscConvEst`, `PetscConvEstGetRate()` 487 @*/ 488 PetscErrorCode PetscConvEstRateView(PetscConvEst ce, const PetscReal alpha[], PetscViewer viewer) 489 { 490 PetscBool isAscii; 491 492 PetscFunctionBegin; 493 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isAscii)); 494 if (isAscii) { 495 PetscInt Nf = ce->Nf, f; 496 497 PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ce)->tablevel)); 498 PetscCall(PetscViewerASCIIPrintf(viewer, "L_2 convergence rate: ")); 499 if (Nf > 1) PetscCall(PetscViewerASCIIPrintf(viewer, "[")); 500 for (f = 0; f < Nf; ++f) { 501 if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", ")); 502 PetscCall(PetscViewerASCIIPrintf(viewer, "%#.2g", (double)alpha[f])); 503 } 504 if (Nf > 1) PetscCall(PetscViewerASCIIPrintf(viewer, "]")); 505 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 506 PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ce)->tablevel)); 507 } 508 PetscFunctionReturn(PETSC_SUCCESS); 509 } 510 511 /*@ 512 PetscConvEstCreate - Create a `PetscConvEst` object 513 514 Collective 515 516 Input Parameter: 517 . comm - The communicator for the `PetscConvEst` object 518 519 Output Parameter: 520 . ce - The `PetscConvEst` object 521 522 Level: beginner 523 524 .seealso: `PetscConvEst`, `PetscConvEstDestroy()`, `PetscConvEstGetConvRate()` 525 @*/ 526 PetscErrorCode PetscConvEstCreate(MPI_Comm comm, PetscConvEst *ce) 527 { 528 PetscFunctionBegin; 529 PetscValidPointer(ce, 2); 530 PetscCall(PetscSysInitializePackage()); 531 PetscCall(PetscHeaderCreate(*ce, PETSC_OBJECT_CLASSID, "PetscConvEst", "ConvergenceEstimator", "SNES", comm, PetscConvEstDestroy, PetscConvEstView)); 532 (*ce)->monitor = PETSC_FALSE; 533 (*ce)->r = 2.0; 534 (*ce)->Nr = 4; 535 (*ce)->event = -1; 536 (*ce)->ops->setsolver = PetscConvEstSetSNES_Private; 537 (*ce)->ops->initguess = PetscConvEstInitGuessSNES_Private; 538 (*ce)->ops->computeerror = PetscConvEstComputeErrorSNES_Private; 539 (*ce)->ops->getconvrate = PetscConvEstGetConvRateSNES_Private; 540 PetscFunctionReturn(PETSC_SUCCESS); 541 } 542