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