16af0ca60SMatthew G. Knepley #include <petscconvest.h> /*I "petscconvest.h" I*/ 26af0ca60SMatthew G. Knepley #include <petscdmplex.h> 36af0ca60SMatthew G. Knepley #include <petscds.h> 46af0ca60SMatthew G. Knepley 56af0ca60SMatthew G. Knepley #include <petsc/private/petscconvestimpl.h> 66af0ca60SMatthew G. Knepley 76af0ca60SMatthew G. Knepley static PetscErrorCode zero_private(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 86af0ca60SMatthew G. Knepley { 96af0ca60SMatthew G. Knepley PetscInt c; 106af0ca60SMatthew G. Knepley for (c = 0; c < Nc; ++c) u[c] = 0.0; 116af0ca60SMatthew G. Knepley return 0; 126af0ca60SMatthew G. Knepley } 136af0ca60SMatthew G. Knepley 146af0ca60SMatthew G. Knepley 156af0ca60SMatthew G. Knepley /*@ 166af0ca60SMatthew G. Knepley PetscConvEstCreate - Create a PetscConvEst object 176af0ca60SMatthew G. Knepley 18d083f849SBarry Smith Collective 196af0ca60SMatthew G. Knepley 206af0ca60SMatthew G. Knepley Input Parameter: 216af0ca60SMatthew G. Knepley . comm - The communicator for the PetscConvEst object 226af0ca60SMatthew G. Knepley 236af0ca60SMatthew G. Knepley Output Parameter: 246af0ca60SMatthew G. Knepley . ce - The PetscConvEst object 256af0ca60SMatthew G. Knepley 266af0ca60SMatthew G. Knepley Level: beginner 276af0ca60SMatthew G. Knepley 286af0ca60SMatthew G. Knepley .seealso: PetscConvEstDestroy(), PetscConvEstGetConvRate() 296af0ca60SMatthew G. Knepley @*/ 306af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstCreate(MPI_Comm comm, PetscConvEst *ce) 316af0ca60SMatthew G. Knepley { 326af0ca60SMatthew G. Knepley PetscErrorCode ierr; 336af0ca60SMatthew G. Knepley 346af0ca60SMatthew G. Knepley PetscFunctionBegin; 356af0ca60SMatthew G. Knepley PetscValidPointer(ce, 2); 366af0ca60SMatthew G. Knepley ierr = PetscSysInitializePackage();CHKERRQ(ierr); 376af0ca60SMatthew G. Knepley ierr = PetscHeaderCreate(*ce, PETSC_OBJECT_CLASSID, "PetscConvEst", "ConvergenceEstimator", "SNES", comm, PetscConvEstDestroy, PetscConvEstView);CHKERRQ(ierr); 386af0ca60SMatthew G. Knepley (*ce)->monitor = PETSC_FALSE; 396af0ca60SMatthew G. Knepley (*ce)->Nr = 4; 406af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 416af0ca60SMatthew G. Knepley } 426af0ca60SMatthew G. Knepley 436af0ca60SMatthew G. Knepley /*@ 446af0ca60SMatthew G. Knepley PetscConvEstDestroy - Destroys a PetscConvEst object 456af0ca60SMatthew G. Knepley 466af0ca60SMatthew G. Knepley Collective on PetscConvEst 476af0ca60SMatthew G. Knepley 486af0ca60SMatthew G. Knepley Input Parameter: 496af0ca60SMatthew G. Knepley . ce - The PetscConvEst object 506af0ca60SMatthew G. Knepley 516af0ca60SMatthew G. Knepley Level: beginner 526af0ca60SMatthew G. Knepley 536af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate() 546af0ca60SMatthew G. Knepley @*/ 556af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstDestroy(PetscConvEst *ce) 566af0ca60SMatthew G. Knepley { 576af0ca60SMatthew G. Knepley PetscErrorCode ierr; 586af0ca60SMatthew G. Knepley 596af0ca60SMatthew G. Knepley PetscFunctionBegin; 606af0ca60SMatthew G. Knepley if (!*ce) PetscFunctionReturn(0); 616af0ca60SMatthew G. Knepley PetscValidHeaderSpecific((*ce),PETSC_OBJECT_CLASSID,1); 626af0ca60SMatthew G. Knepley if (--((PetscObject)(*ce))->refct > 0) { 636af0ca60SMatthew G. Knepley *ce = NULL; 646af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 656af0ca60SMatthew G. Knepley } 6695cbbfd3SMatthew G. Knepley ierr = PetscFree3((*ce)->initGuess, (*ce)->exactSol, (*ce)->ctxs);CHKERRQ(ierr); 676af0ca60SMatthew G. Knepley ierr = PetscFree((*ce)->errors);CHKERRQ(ierr); 686af0ca60SMatthew G. Knepley ierr = PetscHeaderDestroy(ce);CHKERRQ(ierr); 696af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 706af0ca60SMatthew G. Knepley } 716af0ca60SMatthew G. Knepley 726af0ca60SMatthew G. Knepley /*@ 736af0ca60SMatthew G. Knepley PetscConvEstSetFromOptions - Sets a PetscConvEst object from options 746af0ca60SMatthew G. Knepley 756af0ca60SMatthew G. Knepley Collective on PetscConvEst 766af0ca60SMatthew G. Knepley 776af0ca60SMatthew G. Knepley Input Parameters: 786af0ca60SMatthew G. Knepley . ce - The PetscConvEst object 796af0ca60SMatthew G. Knepley 806af0ca60SMatthew G. Knepley Level: beginner 816af0ca60SMatthew G. Knepley 826af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate() 836af0ca60SMatthew G. Knepley @*/ 846af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstSetFromOptions(PetscConvEst ce) 856af0ca60SMatthew G. Knepley { 866af0ca60SMatthew G. Knepley PetscErrorCode ierr; 876af0ca60SMatthew G. Knepley 886af0ca60SMatthew G. Knepley PetscFunctionBegin; 896af0ca60SMatthew G. Knepley ierr = PetscOptionsBegin(PetscObjectComm((PetscObject) ce), "", "Convergence Estimator Options", "PetscConvEst");CHKERRQ(ierr); 9046079b62SMatthew G. Knepley ierr = PetscOptionsInt("-convest_num_refine", "The number of refinements for the convergence check", "PetscConvEst", ce->Nr, &ce->Nr, NULL);CHKERRQ(ierr); 9146079b62SMatthew G. Knepley ierr = PetscOptionsBool("-convest_monitor", "Monitor the error for each convergence check", "PetscConvEst", ce->monitor, &ce->monitor, NULL);CHKERRQ(ierr); 926af0ca60SMatthew G. Knepley ierr = PetscOptionsEnd(); 936af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 946af0ca60SMatthew G. Knepley } 956af0ca60SMatthew G. Knepley 966af0ca60SMatthew G. Knepley /*@ 976af0ca60SMatthew G. Knepley PetscConvEstView - Views a PetscConvEst object 986af0ca60SMatthew G. Knepley 996af0ca60SMatthew G. Knepley Collective on PetscConvEst 1006af0ca60SMatthew G. Knepley 1016af0ca60SMatthew G. Knepley Input Parameters: 1026af0ca60SMatthew G. Knepley + ce - The PetscConvEst object 1036af0ca60SMatthew G. Knepley - viewer - The PetscViewer object 1046af0ca60SMatthew G. Knepley 1056af0ca60SMatthew G. Knepley Level: beginner 1066af0ca60SMatthew G. Knepley 1076af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate() 1086af0ca60SMatthew G. Knepley @*/ 1096af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstView(PetscConvEst ce, PetscViewer viewer) 1106af0ca60SMatthew G. Knepley { 1116af0ca60SMatthew G. Knepley PetscErrorCode ierr; 1126af0ca60SMatthew G. Knepley 1136af0ca60SMatthew G. Knepley PetscFunctionBegin; 1146af0ca60SMatthew G. Knepley ierr = PetscObjectPrintClassNamePrefixType((PetscObject) ce, viewer);CHKERRQ(ierr); 1156af0ca60SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "ConvEst with %D levels\n", ce->Nr+1);CHKERRQ(ierr); 1166af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 1176af0ca60SMatthew G. Knepley } 1186af0ca60SMatthew G. Knepley 1196af0ca60SMatthew G. Knepley /*@ 1206af0ca60SMatthew G. Knepley PetscConvEstGetSolver - Gets the solver used to produce discrete solutions 1216af0ca60SMatthew G. Knepley 1226af0ca60SMatthew G. Knepley Not collective 1236af0ca60SMatthew G. Knepley 1246af0ca60SMatthew G. Knepley Input Parameter: 1256af0ca60SMatthew G. Knepley . ce - The PetscConvEst object 1266af0ca60SMatthew G. Knepley 1276af0ca60SMatthew G. Knepley Output Parameter: 1286af0ca60SMatthew G. Knepley . snes - The solver 1296af0ca60SMatthew G. Knepley 1306af0ca60SMatthew G. Knepley Level: intermediate 1316af0ca60SMatthew G. Knepley 1326af0ca60SMatthew G. Knepley .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate() 1336af0ca60SMatthew G. Knepley @*/ 1346af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstGetSolver(PetscConvEst ce, SNES *snes) 1356af0ca60SMatthew G. Knepley { 1366af0ca60SMatthew G. Knepley PetscFunctionBegin; 1376af0ca60SMatthew G. Knepley PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 1386af0ca60SMatthew G. Knepley PetscValidPointer(snes, 2); 1396af0ca60SMatthew G. Knepley *snes = ce->snes; 1406af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 1416af0ca60SMatthew G. Knepley } 1426af0ca60SMatthew G. Knepley 1436af0ca60SMatthew G. Knepley /*@ 1446af0ca60SMatthew G. Knepley PetscConvEstSetSolver - Sets the solver used to produce discrete solutions 1456af0ca60SMatthew G. Knepley 1466af0ca60SMatthew G. Knepley Not collective 1476af0ca60SMatthew G. Knepley 1486af0ca60SMatthew G. Knepley Input Parameters: 1496af0ca60SMatthew G. Knepley + ce - The PetscConvEst object 1506af0ca60SMatthew G. Knepley - snes - The solver 1516af0ca60SMatthew G. Knepley 1526af0ca60SMatthew G. Knepley Level: intermediate 1536af0ca60SMatthew G. Knepley 1546af0ca60SMatthew G. Knepley Note: The solver MUST have an attached DM/DS, so that we know the exact solution 1556af0ca60SMatthew G. Knepley 1566af0ca60SMatthew G. Knepley .seealso: PetscConvEstGetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate() 1576af0ca60SMatthew G. Knepley @*/ 1586af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstSetSolver(PetscConvEst ce, SNES snes) 1596af0ca60SMatthew G. Knepley { 1606af0ca60SMatthew G. Knepley PetscErrorCode ierr; 1616af0ca60SMatthew G. Knepley 1626af0ca60SMatthew G. Knepley PetscFunctionBegin; 1636af0ca60SMatthew G. Knepley PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 1646af0ca60SMatthew G. Knepley PetscValidHeaderSpecific(snes, SNES_CLASSID, 2); 1656af0ca60SMatthew G. Knepley ce->snes = snes; 1666af0ca60SMatthew G. Knepley ierr = SNESGetDM(ce->snes, &ce->idm);CHKERRQ(ierr); 1676af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 1686af0ca60SMatthew G. Knepley } 1696af0ca60SMatthew G. Knepley 1706af0ca60SMatthew G. Knepley /*@ 1716af0ca60SMatthew G. Knepley PetscConvEstSetUp - After the solver is specified, we create structures for estimating convergence 1726af0ca60SMatthew G. Knepley 1736af0ca60SMatthew G. Knepley Collective on PetscConvEst 1746af0ca60SMatthew G. Knepley 1756af0ca60SMatthew G. Knepley Input Parameters: 1766af0ca60SMatthew G. Knepley . ce - The PetscConvEst object 1776af0ca60SMatthew G. Knepley 1786af0ca60SMatthew G. Knepley Level: beginner 1796af0ca60SMatthew G. Knepley 1806af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate() 1816af0ca60SMatthew G. Knepley @*/ 1820955ed61SMatthew G. Knepley PetscErrorCode PetscConvEstSetUp(PetscConvEst ce) 1836af0ca60SMatthew G. Knepley { 1846af0ca60SMatthew G. Knepley PetscDS prob; 1856af0ca60SMatthew G. Knepley PetscInt f; 1866af0ca60SMatthew G. Knepley PetscErrorCode ierr; 1876af0ca60SMatthew G. Knepley 1886af0ca60SMatthew G. Knepley PetscFunctionBegin; 1896af0ca60SMatthew G. Knepley ierr = DMGetDS(ce->idm, &prob);CHKERRQ(ierr); 1906af0ca60SMatthew G. Knepley ierr = PetscDSGetNumFields(prob, &ce->Nf);CHKERRQ(ierr); 1916af0ca60SMatthew G. Knepley ierr = PetscMalloc1((ce->Nr+1)*ce->Nf, &ce->errors);CHKERRQ(ierr); 19295cbbfd3SMatthew G. Knepley ierr = PetscMalloc3(ce->Nf, &ce->initGuess, ce->Nf, &ce->exactSol, ce->Nf, &ce->ctxs);CHKERRQ(ierr); 1936af0ca60SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) ce->initGuess[f] = zero_private; 1946af0ca60SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 19595cbbfd3SMatthew G. Knepley ierr = PetscDSGetExactSolution(prob, f, &ce->exactSol[f], &ce->ctxs[f]);CHKERRQ(ierr); 1966af0ca60SMatthew G. Knepley if (!ce->exactSol[f]) SETERRQ1(PetscObjectComm((PetscObject) ce), PETSC_ERR_ARG_WRONG, "DS must contain exact solution functions in order to estimate convergence, missing for field %D", f); 1976af0ca60SMatthew G. Knepley } 1986af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 1996af0ca60SMatthew G. Knepley } 2006af0ca60SMatthew G. Knepley 2016af0ca60SMatthew G. Knepley /*@ 2026af0ca60SMatthew G. Knepley PetscConvEstGetConvRate - Returns an estimate of the convergence rate for the discretization 2036af0ca60SMatthew G. Knepley 2046af0ca60SMatthew G. Knepley Not collective 2056af0ca60SMatthew G. Knepley 2066af0ca60SMatthew G. Knepley Input Parameter: 2076af0ca60SMatthew G. Knepley . ce - The PetscConvEst object 2086af0ca60SMatthew G. Knepley 2096af0ca60SMatthew G. Knepley Output Parameter: 21046079b62SMatthew G. Knepley . alpha - The convergence rate for each field 2116af0ca60SMatthew G. Knepley 2126af0ca60SMatthew G. Knepley Note: The convergence rate alpha is defined by 2136af0ca60SMatthew G. Knepley $ || u_h - u_exact || < C h^alpha 2146af0ca60SMatthew G. Knepley where u_h is the discrete solution, and h is a measure of the discretization size. 2156af0ca60SMatthew G. Knepley 2166af0ca60SMatthew G. Knepley We solve a series of problems on refined meshes, calculate an error based upon the exact solution in the DS, 2176af0ca60SMatthew G. Knepley and then fit the result to our model above using linear regression. 2186af0ca60SMatthew G. Knepley 2196af0ca60SMatthew G. Knepley Options database keys: 2206af0ca60SMatthew G. Knepley . -snes_convergence_estimate : Execute convergence estimation and print out the rate 2216af0ca60SMatthew G. Knepley 2226af0ca60SMatthew G. Knepley Level: intermediate 2236af0ca60SMatthew G. Knepley 2246af0ca60SMatthew G. Knepley .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate() 2256af0ca60SMatthew G. Knepley @*/ 22646079b62SMatthew G. Knepley PetscErrorCode PetscConvEstGetConvRate(PetscConvEst ce, PetscReal alpha[]) 2276af0ca60SMatthew G. Knepley { 2286af0ca60SMatthew G. Knepley DM *dm; 2296af0ca60SMatthew G. Knepley PetscObject disc; 2306af0ca60SMatthew G. Knepley MPI_Comm comm; 2316af0ca60SMatthew G. Knepley const char *uname, *dmname; 2326af0ca60SMatthew G. Knepley void *ctx; 2336af0ca60SMatthew G. Knepley Vec u; 2346af0ca60SMatthew G. Knepley PetscReal t = 0.0, *x, *y, slope, intercept; 23573269098SMatthew G. Knepley PetscInt *dof, dim, Nr = ce->Nr, r, f, oldlevel, oldnlev; 236adfa7136SMatthew G. Knepley PetscLogEvent event; 2376af0ca60SMatthew G. Knepley PetscErrorCode ierr; 2386af0ca60SMatthew G. Knepley 2396af0ca60SMatthew G. Knepley PetscFunctionBegin; 2406af0ca60SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) ce, &comm);CHKERRQ(ierr); 2416af0ca60SMatthew G. Knepley ierr = DMGetDimension(ce->idm, &dim);CHKERRQ(ierr); 2426af0ca60SMatthew G. Knepley ierr = DMGetApplicationContext(ce->idm, &ctx);CHKERRQ(ierr); 2436af0ca60SMatthew G. Knepley ierr = DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE);CHKERRQ(ierr); 24473269098SMatthew G. Knepley ierr = DMGetRefineLevel(ce->idm, &oldlevel);CHKERRQ(ierr); 245adfa7136SMatthew G. Knepley ierr = PetscMalloc2((Nr+1), &dm, (Nr+1)*ce->Nf, &dof);CHKERRQ(ierr); 2466af0ca60SMatthew G. Knepley dm[0] = ce->idm; 24746079b62SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) alpha[f] = 0.0; 2486af0ca60SMatthew G. Knepley /* Loop over meshes */ 249adfa7136SMatthew G. Knepley ierr = PetscLogEventRegister("ConvEst Error", PETSC_OBJECT_CLASSID, &event);CHKERRQ(ierr); 2506af0ca60SMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 251adfa7136SMatthew G. Knepley PetscLogStage stage; 252adfa7136SMatthew G. Knepley char stageName[PETSC_MAX_PATH_LEN]; 253adfa7136SMatthew G. Knepley 254adfa7136SMatthew G. Knepley ierr = PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN-1, "ConvEst Refinement Level %D", r);CHKERRQ(ierr); 255adfa7136SMatthew G. Knepley ierr = PetscLogStageRegister(stageName, &stage);CHKERRQ(ierr); 256adfa7136SMatthew G. Knepley ierr = PetscLogStagePush(stage);CHKERRQ(ierr); 2576af0ca60SMatthew G. Knepley if (r > 0) { 2586af0ca60SMatthew G. Knepley ierr = DMRefine(dm[r-1], MPI_COMM_NULL, &dm[r]);CHKERRQ(ierr); 2596af0ca60SMatthew G. Knepley ierr = DMSetCoarseDM(dm[r], dm[r-1]);CHKERRQ(ierr); 260e5e52638SMatthew G. Knepley ierr = DMCopyDisc(ce->idm, dm[r]);CHKERRQ(ierr); 261ca3d3a14SMatthew G. Knepley ierr = DMCopyTransform(ce->idm, dm[r]);CHKERRQ(ierr); 2626af0ca60SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) dm[r-1], &dmname);CHKERRQ(ierr); 2636af0ca60SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) dm[r], dmname);CHKERRQ(ierr); 2648f892730SMatthew G. Knepley for (f = 0; f <= ce->Nf; ++f) { 2658f892730SMatthew G. Knepley PetscErrorCode (*nspconstr)(DM, PetscInt, MatNullSpace *); 2668f892730SMatthew G. Knepley ierr = DMGetNullSpaceConstructor(dm[r-1], f, &nspconstr);CHKERRQ(ierr); 2678f892730SMatthew G. Knepley ierr = DMSetNullSpaceConstructor(dm[r], f, nspconstr);CHKERRQ(ierr); 2688f892730SMatthew G. Knepley } 2696af0ca60SMatthew G. Knepley } 2706af0ca60SMatthew G. Knepley ierr = DMViewFromOptions(dm[r], NULL, "-conv_dm_view");CHKERRQ(ierr); 2716af0ca60SMatthew G. Knepley /* Create solution */ 2726af0ca60SMatthew G. Knepley ierr = DMCreateGlobalVector(dm[r], &u);CHKERRQ(ierr); 273e5e52638SMatthew G. Knepley ierr = DMGetField(dm[r], 0, NULL, &disc);CHKERRQ(ierr); 274239a0881SMatthew G. Knepley ierr = PetscObjectGetName(disc, &uname);CHKERRQ(ierr); 2756af0ca60SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) u, uname);CHKERRQ(ierr); 2766af0ca60SMatthew G. Knepley /* Setup solver */ 2776af0ca60SMatthew G. Knepley ierr = SNESReset(ce->snes);CHKERRQ(ierr); 2786af0ca60SMatthew G. Knepley ierr = SNESSetDM(ce->snes, dm[r]);CHKERRQ(ierr); 2796af0ca60SMatthew G. Knepley ierr = DMPlexSetSNESLocalFEM(dm[r], ctx, ctx, ctx);CHKERRQ(ierr); 2806af0ca60SMatthew G. Knepley ierr = SNESSetFromOptions(ce->snes);CHKERRQ(ierr); 2816af0ca60SMatthew G. Knepley /* Create initial guess */ 28295cbbfd3SMatthew G. Knepley ierr = DMProjectFunction(dm[r], t, ce->initGuess, ce->ctxs, INSERT_VALUES, u);CHKERRQ(ierr); 2836af0ca60SMatthew G. Knepley ierr = SNESSolve(ce->snes, NULL, u);CHKERRQ(ierr); 2845d68e14cSMatthew G. Knepley ierr = PetscLogEventBegin(event, ce, 0, 0, 0);CHKERRQ(ierr); 28595cbbfd3SMatthew G. Knepley ierr = DMComputeL2FieldDiff(dm[r], t, ce->exactSol, ce->ctxs, u, &ce->errors[r*ce->Nf]);CHKERRQ(ierr); 286adfa7136SMatthew G. Knepley ierr = PetscLogEventEnd(event, ce, 0, 0, 0);CHKERRQ(ierr); 287adfa7136SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 288a56d3bf6SMatthew G. Knepley PetscSection s, fs; 289a56d3bf6SMatthew G. Knepley PetscInt lsize; 290a56d3bf6SMatthew G. Knepley 291a56d3bf6SMatthew G. Knepley /* Could use DMGetOutputDM() to add in Dirichlet dofs */ 292*92fd8e1eSJed Brown ierr = DMGetLocalSection(dm[r], &s);CHKERRQ(ierr); 293a56d3bf6SMatthew G. Knepley ierr = PetscSectionGetField(s, f, &fs);CHKERRQ(ierr); 294a56d3bf6SMatthew G. Knepley ierr = PetscSectionGetConstrainedStorageSize(fs, &lsize);CHKERRQ(ierr); 295a56d3bf6SMatthew G. Knepley ierr = MPI_Allreduce(&lsize, &dof[r*ce->Nf+f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) ce->snes));CHKERRQ(ierr); 296a56d3bf6SMatthew G. Knepley ierr = PetscLogEventSetDof(event, f, dof[r*ce->Nf+f]);CHKERRQ(ierr); 297a56d3bf6SMatthew G. Knepley ierr = PetscLogEventSetError(event, f, ce->errors[r*ce->Nf+f]);CHKERRQ(ierr); 298adfa7136SMatthew G. Knepley } 2996af0ca60SMatthew G. Knepley /* Monitor */ 3006af0ca60SMatthew G. Knepley if (ce->monitor) { 3016af0ca60SMatthew G. Knepley PetscReal *errors = &ce->errors[r*ce->Nf]; 3026af0ca60SMatthew G. Knepley 303adfa7136SMatthew G. Knepley ierr = PetscPrintf(comm, "L_2 Error: ");CHKERRQ(ierr); 304adfa7136SMatthew G. Knepley if (ce->Nf > 1) {ierr = PetscPrintf(comm, "[");CHKERRQ(ierr);} 3056af0ca60SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 3066af0ca60SMatthew G. Knepley if (f > 0) {ierr = PetscPrintf(comm, ", ");CHKERRQ(ierr);} 3076af0ca60SMatthew G. Knepley if (errors[f] < 1.0e-11) {ierr = PetscPrintf(comm, "< 1e-11");CHKERRQ(ierr);} 3088dfade3bSBarry Smith else {ierr = PetscPrintf(comm, "%g", (double)errors[f]);CHKERRQ(ierr);} 3096af0ca60SMatthew G. Knepley } 310adfa7136SMatthew G. Knepley if (ce->Nf > 1) {ierr = PetscPrintf(comm, "]");CHKERRQ(ierr);} 311adfa7136SMatthew G. Knepley ierr = PetscPrintf(comm, "\n");CHKERRQ(ierr); 3126af0ca60SMatthew G. Knepley } 31373269098SMatthew G. Knepley if (!r) { 31473269098SMatthew G. Knepley /* PCReset() does not wipe out the level structure */ 31573269098SMatthew G. Knepley KSP ksp; 31673269098SMatthew G. Knepley PC pc; 31773269098SMatthew G. Knepley 31873269098SMatthew G. Knepley ierr = SNESGetKSP(ce->snes, &ksp);CHKERRQ(ierr); 31973269098SMatthew G. Knepley ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 32073269098SMatthew G. Knepley ierr = PCMGGetLevels(pc, &oldnlev);CHKERRQ(ierr); 32173269098SMatthew G. Knepley } 3226af0ca60SMatthew G. Knepley /* Cleanup */ 3236af0ca60SMatthew G. Knepley ierr = VecDestroy(&u);CHKERRQ(ierr); 324adfa7136SMatthew G. Knepley ierr = PetscLogStagePop();CHKERRQ(ierr); 3256af0ca60SMatthew G. Knepley } 3266af0ca60SMatthew G. Knepley for (r = 1; r <= Nr; ++r) { 3276af0ca60SMatthew G. Knepley ierr = DMDestroy(&dm[r]);CHKERRQ(ierr); 3286af0ca60SMatthew G. Knepley } 3296af0ca60SMatthew G. Knepley /* Fit convergence rate */ 3306af0ca60SMatthew G. Knepley ierr = PetscMalloc2(Nr+1, &x, Nr+1, &y);CHKERRQ(ierr); 33146079b62SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 3326af0ca60SMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 333adfa7136SMatthew G. Knepley x[r] = PetscLog10Real(dof[r*ce->Nf+f]); 33446079b62SMatthew G. Knepley y[r] = PetscLog10Real(ce->errors[r*ce->Nf+f]); 3356af0ca60SMatthew G. Knepley } 336bebf13c0SMatthew G. Knepley ierr = PetscLinearRegression(Nr+1, x, y, &slope, &intercept);CHKERRQ(ierr); 3376af0ca60SMatthew G. Knepley /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */ 33846079b62SMatthew G. Knepley alpha[f] = -slope * dim; 33946079b62SMatthew G. Knepley } 34046079b62SMatthew G. Knepley ierr = PetscFree2(x, y);CHKERRQ(ierr); 3416af0ca60SMatthew G. Knepley ierr = PetscFree2(dm, dof);CHKERRQ(ierr); 3422cae373cSMatthew G. Knepley /* Restore solver */ 3432cae373cSMatthew G. Knepley ierr = SNESReset(ce->snes);CHKERRQ(ierr); 34473269098SMatthew G. Knepley { 34573269098SMatthew G. Knepley /* PCReset() does not wipe out the level structure */ 34673269098SMatthew G. Knepley KSP ksp; 34773269098SMatthew G. Knepley PC pc; 34873269098SMatthew G. Knepley 34973269098SMatthew G. Knepley ierr = SNESGetKSP(ce->snes, &ksp);CHKERRQ(ierr); 35073269098SMatthew G. Knepley ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 35173269098SMatthew G. Knepley ierr = PCMGSetLevels(pc, oldnlev, NULL);CHKERRQ(ierr); 35273269098SMatthew G. Knepley ierr = DMSetRefineLevel(ce->idm, oldlevel);CHKERRQ(ierr); /* The damn DMCoarsen() calls in PCMG can reset this */ 35373269098SMatthew G. Knepley } 3542cae373cSMatthew G. Knepley ierr = SNESSetDM(ce->snes, ce->idm);CHKERRQ(ierr); 3552cae373cSMatthew G. Knepley ierr = DMPlexSetSNESLocalFEM(ce->idm, ctx, ctx, ctx);CHKERRQ(ierr); 3562cae373cSMatthew G. Knepley ierr = SNESSetFromOptions(ce->snes);CHKERRQ(ierr); 3576af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 3586af0ca60SMatthew G. Knepley } 3596af0ca60SMatthew G. Knepley 3606af0ca60SMatthew G. Knepley /*@ 3616af0ca60SMatthew G. Knepley PetscConvEstRateView - Displays the convergence rate to a viewer 3626af0ca60SMatthew G. Knepley 3636af0ca60SMatthew G. Knepley Collective on SNES 3646af0ca60SMatthew G. Knepley 3656af0ca60SMatthew G. Knepley Parameter: 3666af0ca60SMatthew G. Knepley + snes - iterative context obtained from SNESCreate() 36746079b62SMatthew G. Knepley . alpha - the convergence rate for each field 3686af0ca60SMatthew G. Knepley - viewer - the viewer to display the reason 3696af0ca60SMatthew G. Knepley 3706af0ca60SMatthew G. Knepley Options Database Keys: 3716af0ca60SMatthew G. Knepley . -snes_convergence_estimate - print the convergence rate 3726af0ca60SMatthew G. Knepley 3736af0ca60SMatthew G. Knepley Level: developer 3746af0ca60SMatthew G. Knepley 3756af0ca60SMatthew G. Knepley .seealso: PetscConvEstGetRate() 3766af0ca60SMatthew G. Knepley @*/ 377a56d3bf6SMatthew G. Knepley PetscErrorCode PetscConvEstRateView(PetscConvEst ce, const PetscReal alpha[], PetscViewer viewer) 3786af0ca60SMatthew G. Knepley { 3796af0ca60SMatthew G. Knepley PetscBool isAscii; 3806af0ca60SMatthew G. Knepley PetscErrorCode ierr; 3816af0ca60SMatthew G. Knepley 3826af0ca60SMatthew G. Knepley PetscFunctionBegin; 3836af0ca60SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr); 3846af0ca60SMatthew G. Knepley if (isAscii) { 385a56d3bf6SMatthew G. Knepley DM dm; 386a56d3bf6SMatthew G. Knepley PetscInt Nf, f; 387a56d3bf6SMatthew G. Knepley 388a56d3bf6SMatthew G. Knepley ierr = SNESGetDM(ce->snes, &dm);CHKERRQ(ierr); 389a56d3bf6SMatthew G. Knepley ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 3906af0ca60SMatthew G. Knepley ierr = PetscViewerASCIIAddTab(viewer, ((PetscObject) ce)->tablevel);CHKERRQ(ierr); 391adfa7136SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "L_2 convergence rate: ");CHKERRQ(ierr); 392a56d3bf6SMatthew G. Knepley if (Nf > 1) {ierr = PetscViewerASCIIPrintf(viewer, "[");CHKERRQ(ierr);} 393a56d3bf6SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 39446079b62SMatthew G. Knepley if (f > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);} 3955f9b3039SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "%#.2g", (double) alpha[f]);CHKERRQ(ierr); 39646079b62SMatthew G. Knepley } 397a56d3bf6SMatthew G. Knepley if (Nf > 1) {ierr = PetscViewerASCIIPrintf(viewer, "]");CHKERRQ(ierr);} 398adfa7136SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 3996af0ca60SMatthew G. Knepley ierr = PetscViewerASCIISubtractTab(viewer, ((PetscObject) ce)->tablevel);CHKERRQ(ierr); 4006af0ca60SMatthew G. Knepley } 4016af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 4026af0ca60SMatthew G. Knepley } 403