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