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 #include <petscblaslapack.h> 56af0ca60SMatthew G. Knepley 66af0ca60SMatthew G. Knepley #include <petsc/private/petscconvestimpl.h> 76af0ca60SMatthew G. Knepley 86af0ca60SMatthew G. Knepley static PetscErrorCode zero_private(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 96af0ca60SMatthew G. Knepley { 106af0ca60SMatthew G. Knepley PetscInt c; 116af0ca60SMatthew G. Knepley for (c = 0; c < Nc; ++c) u[c] = 0.0; 126af0ca60SMatthew G. Knepley return 0; 136af0ca60SMatthew G. Knepley } 146af0ca60SMatthew G. Knepley 156af0ca60SMatthew G. Knepley 166af0ca60SMatthew G. Knepley /*@ 176af0ca60SMatthew G. Knepley PetscConvEstCreate - Create a PetscConvEst object 186af0ca60SMatthew G. Knepley 196af0ca60SMatthew G. Knepley Collective on MPI_Comm 206af0ca60SMatthew G. Knepley 216af0ca60SMatthew G. Knepley Input Parameter: 226af0ca60SMatthew G. Knepley . comm - The communicator for the PetscConvEst object 236af0ca60SMatthew G. Knepley 246af0ca60SMatthew G. Knepley Output Parameter: 256af0ca60SMatthew G. Knepley . ce - The PetscConvEst object 266af0ca60SMatthew G. Knepley 276af0ca60SMatthew G. Knepley Level: beginner 286af0ca60SMatthew G. Knepley 296af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence, create 306af0ca60SMatthew G. Knepley .seealso: PetscConvEstDestroy(), PetscConvEstGetConvRate() 316af0ca60SMatthew G. Knepley @*/ 326af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstCreate(MPI_Comm comm, PetscConvEst *ce) 336af0ca60SMatthew G. Knepley { 346af0ca60SMatthew G. Knepley PetscErrorCode ierr; 356af0ca60SMatthew G. Knepley 366af0ca60SMatthew G. Knepley PetscFunctionBegin; 376af0ca60SMatthew G. Knepley PetscValidPointer(ce, 2); 386af0ca60SMatthew G. Knepley ierr = PetscSysInitializePackage();CHKERRQ(ierr); 396af0ca60SMatthew G. Knepley ierr = PetscHeaderCreate(*ce, PETSC_OBJECT_CLASSID, "PetscConvEst", "ConvergenceEstimator", "SNES", comm, PetscConvEstDestroy, PetscConvEstView);CHKERRQ(ierr); 406af0ca60SMatthew G. Knepley (*ce)->monitor = PETSC_FALSE; 416af0ca60SMatthew G. Knepley (*ce)->Nr = 4; 426af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 436af0ca60SMatthew G. Knepley } 446af0ca60SMatthew G. Knepley 456af0ca60SMatthew G. Knepley /*@ 466af0ca60SMatthew G. Knepley PetscConvEstDestroy - Destroys a PetscConvEst object 476af0ca60SMatthew G. Knepley 486af0ca60SMatthew G. Knepley Collective on PetscConvEst 496af0ca60SMatthew G. Knepley 506af0ca60SMatthew G. Knepley Input Parameter: 516af0ca60SMatthew G. Knepley . ce - The PetscConvEst object 526af0ca60SMatthew G. Knepley 536af0ca60SMatthew G. Knepley Level: beginner 546af0ca60SMatthew G. Knepley 556af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence, destroy 566af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate() 576af0ca60SMatthew G. Knepley @*/ 586af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstDestroy(PetscConvEst *ce) 596af0ca60SMatthew G. Knepley { 606af0ca60SMatthew G. Knepley PetscErrorCode ierr; 616af0ca60SMatthew G. Knepley 626af0ca60SMatthew G. Knepley PetscFunctionBegin; 636af0ca60SMatthew G. Knepley if (!*ce) PetscFunctionReturn(0); 646af0ca60SMatthew G. Knepley PetscValidHeaderSpecific((*ce),PETSC_OBJECT_CLASSID,1); 656af0ca60SMatthew G. Knepley if (--((PetscObject)(*ce))->refct > 0) { 666af0ca60SMatthew G. Knepley *ce = NULL; 676af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 686af0ca60SMatthew G. Knepley } 696af0ca60SMatthew G. Knepley ierr = PetscFree2((*ce)->initGuess, (*ce)->exactSol);CHKERRQ(ierr); 706af0ca60SMatthew G. Knepley ierr = PetscFree((*ce)->errors);CHKERRQ(ierr); 716af0ca60SMatthew G. Knepley ierr = PetscHeaderDestroy(ce);CHKERRQ(ierr); 726af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 736af0ca60SMatthew G. Knepley } 746af0ca60SMatthew G. Knepley 756af0ca60SMatthew G. Knepley /*@ 766af0ca60SMatthew G. Knepley PetscConvEstSetFromOptions - Sets a PetscConvEst object from options 776af0ca60SMatthew G. Knepley 786af0ca60SMatthew G. Knepley Collective on PetscConvEst 796af0ca60SMatthew G. Knepley 806af0ca60SMatthew G. Knepley Input Parameters: 816af0ca60SMatthew G. Knepley . ce - The PetscConvEst object 826af0ca60SMatthew G. Knepley 836af0ca60SMatthew G. Knepley Level: beginner 846af0ca60SMatthew G. Knepley 856af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence, options 866af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate() 876af0ca60SMatthew G. Knepley @*/ 886af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstSetFromOptions(PetscConvEst ce) 896af0ca60SMatthew G. Knepley { 906af0ca60SMatthew G. Knepley PetscErrorCode ierr; 916af0ca60SMatthew G. Knepley 926af0ca60SMatthew G. Knepley PetscFunctionBegin; 936af0ca60SMatthew G. Knepley ierr = PetscOptionsBegin(PetscObjectComm((PetscObject) ce), "", "Convergence Estimator Options", "PetscConvEst");CHKERRQ(ierr); 9446079b62SMatthew G. Knepley ierr = PetscOptionsInt("-convest_num_refine", "The number of refinements for the convergence check", "PetscConvEst", ce->Nr, &ce->Nr, NULL);CHKERRQ(ierr); 9546079b62SMatthew G. Knepley ierr = PetscOptionsBool("-convest_monitor", "Monitor the error for each convergence check", "PetscConvEst", ce->monitor, &ce->monitor, NULL);CHKERRQ(ierr); 966af0ca60SMatthew G. Knepley ierr = PetscOptionsEnd(); 976af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 986af0ca60SMatthew G. Knepley } 996af0ca60SMatthew G. Knepley 1006af0ca60SMatthew G. Knepley /*@ 1016af0ca60SMatthew G. Knepley PetscConvEstView - Views a PetscConvEst object 1026af0ca60SMatthew G. Knepley 1036af0ca60SMatthew G. Knepley Collective on PetscConvEst 1046af0ca60SMatthew G. Knepley 1056af0ca60SMatthew G. Knepley Input Parameters: 1066af0ca60SMatthew G. Knepley + ce - The PetscConvEst object 1076af0ca60SMatthew G. Knepley - viewer - The PetscViewer object 1086af0ca60SMatthew G. Knepley 1096af0ca60SMatthew G. Knepley Level: beginner 1106af0ca60SMatthew G. Knepley 1116af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence, view 1126af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate() 1136af0ca60SMatthew G. Knepley @*/ 1146af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstView(PetscConvEst ce, PetscViewer viewer) 1156af0ca60SMatthew G. Knepley { 1166af0ca60SMatthew G. Knepley PetscErrorCode ierr; 1176af0ca60SMatthew G. Knepley 1186af0ca60SMatthew G. Knepley PetscFunctionBegin; 1196af0ca60SMatthew G. Knepley ierr = PetscObjectPrintClassNamePrefixType((PetscObject) ce, viewer);CHKERRQ(ierr); 1206af0ca60SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "ConvEst with %D levels\n", ce->Nr+1);CHKERRQ(ierr); 1216af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 1226af0ca60SMatthew G. Knepley } 1236af0ca60SMatthew G. Knepley 1246af0ca60SMatthew G. Knepley /*@ 1256af0ca60SMatthew G. Knepley PetscConvEstGetSolver - Gets the solver used to produce discrete solutions 1266af0ca60SMatthew G. Knepley 1276af0ca60SMatthew G. Knepley Not collective 1286af0ca60SMatthew G. Knepley 1296af0ca60SMatthew G. Knepley Input Parameter: 1306af0ca60SMatthew G. Knepley . ce - The PetscConvEst object 1316af0ca60SMatthew G. Knepley 1326af0ca60SMatthew G. Knepley Output Parameter: 1336af0ca60SMatthew G. Knepley . snes - The solver 1346af0ca60SMatthew G. Knepley 1356af0ca60SMatthew G. Knepley Level: intermediate 1366af0ca60SMatthew G. Knepley 1376af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence 1386af0ca60SMatthew G. Knepley .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate() 1396af0ca60SMatthew G. Knepley @*/ 1406af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstGetSolver(PetscConvEst ce, SNES *snes) 1416af0ca60SMatthew G. Knepley { 1426af0ca60SMatthew G. Knepley PetscFunctionBegin; 1436af0ca60SMatthew G. Knepley PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 1446af0ca60SMatthew G. Knepley PetscValidPointer(snes, 2); 1456af0ca60SMatthew G. Knepley *snes = ce->snes; 1466af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 1476af0ca60SMatthew G. Knepley } 1486af0ca60SMatthew G. Knepley 1496af0ca60SMatthew G. Knepley /*@ 1506af0ca60SMatthew G. Knepley PetscConvEstSetSolver - Sets the solver used to produce discrete solutions 1516af0ca60SMatthew G. Knepley 1526af0ca60SMatthew G. Knepley Not collective 1536af0ca60SMatthew G. Knepley 1546af0ca60SMatthew G. Knepley Input Parameters: 1556af0ca60SMatthew G. Knepley + ce - The PetscConvEst object 1566af0ca60SMatthew G. Knepley - snes - The solver 1576af0ca60SMatthew G. Knepley 1586af0ca60SMatthew G. Knepley Level: intermediate 1596af0ca60SMatthew G. Knepley 1606af0ca60SMatthew G. Knepley Note: The solver MUST have an attached DM/DS, so that we know the exact solution 1616af0ca60SMatthew G. Knepley 1626af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence 1636af0ca60SMatthew G. Knepley .seealso: PetscConvEstGetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate() 1646af0ca60SMatthew G. Knepley @*/ 1656af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstSetSolver(PetscConvEst ce, SNES snes) 1666af0ca60SMatthew G. Knepley { 1676af0ca60SMatthew G. Knepley PetscErrorCode ierr; 1686af0ca60SMatthew G. Knepley 1696af0ca60SMatthew G. Knepley PetscFunctionBegin; 1706af0ca60SMatthew G. Knepley PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 1716af0ca60SMatthew G. Knepley PetscValidHeaderSpecific(snes, SNES_CLASSID, 2); 1726af0ca60SMatthew G. Knepley ce->snes = snes; 1736af0ca60SMatthew G. Knepley ierr = SNESGetDM(ce->snes, &ce->idm);CHKERRQ(ierr); 1746af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 1756af0ca60SMatthew G. Knepley } 1766af0ca60SMatthew G. Knepley 1776af0ca60SMatthew G. Knepley /*@ 1786af0ca60SMatthew G. Knepley PetscConvEstSetUp - After the solver is specified, we create structures for estimating convergence 1796af0ca60SMatthew G. Knepley 1806af0ca60SMatthew G. Knepley Collective on PetscConvEst 1816af0ca60SMatthew G. Knepley 1826af0ca60SMatthew G. Knepley Input Parameters: 1836af0ca60SMatthew G. Knepley . ce - The PetscConvEst object 1846af0ca60SMatthew G. Knepley 1856af0ca60SMatthew G. Knepley Level: beginner 1866af0ca60SMatthew G. Knepley 1876af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence, setup 1886af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate() 1896af0ca60SMatthew G. Knepley @*/ 1900955ed61SMatthew G. Knepley PetscErrorCode PetscConvEstSetUp(PetscConvEst ce) 1916af0ca60SMatthew G. Knepley { 1926af0ca60SMatthew G. Knepley PetscDS prob; 1936af0ca60SMatthew G. Knepley PetscInt f; 1946af0ca60SMatthew G. Knepley PetscErrorCode ierr; 1956af0ca60SMatthew G. Knepley 1966af0ca60SMatthew G. Knepley PetscFunctionBegin; 1976af0ca60SMatthew G. Knepley ierr = DMGetDS(ce->idm, &prob);CHKERRQ(ierr); 1986af0ca60SMatthew G. Knepley ierr = PetscDSGetNumFields(prob, &ce->Nf);CHKERRQ(ierr); 1996af0ca60SMatthew G. Knepley ierr = PetscMalloc1((ce->Nr+1)*ce->Nf, &ce->errors);CHKERRQ(ierr); 2006af0ca60SMatthew G. Knepley ierr = PetscMalloc2(ce->Nf, &ce->initGuess, ce->Nf, &ce->exactSol);CHKERRQ(ierr); 2016af0ca60SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) ce->initGuess[f] = zero_private; 2026af0ca60SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 2036af0ca60SMatthew G. Knepley ierr = PetscDSGetExactSolution(prob, f, &ce->exactSol[f]);CHKERRQ(ierr); 2046af0ca60SMatthew 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); 2056af0ca60SMatthew G. Knepley } 2066af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 2076af0ca60SMatthew G. Knepley } 2086af0ca60SMatthew G. Knepley 2096af0ca60SMatthew G. Knepley static PetscErrorCode PetscConvEstLinearRegression_Private(PetscConvEst ce, PetscInt n, const PetscReal x[], const PetscReal y[], PetscReal *slope, PetscReal *intercept) 2106af0ca60SMatthew G. Knepley { 21128cc7b48SMatthew G. Knepley PetscScalar H[4]; 21228cc7b48SMatthew G. Knepley PetscReal *X, *Y, beta[2]; 2136af0ca60SMatthew G. Knepley PetscInt i, j, k; 2146af0ca60SMatthew G. Knepley PetscErrorCode ierr; 2156af0ca60SMatthew G. Knepley 2166af0ca60SMatthew G. Knepley PetscFunctionBegin; 2178fa78361SMatthew G. Knepley *slope = *intercept = 0.0; 2186af0ca60SMatthew G. Knepley ierr = PetscMalloc2(n*2, &X, n*2, &Y);CHKERRQ(ierr); 2196af0ca60SMatthew G. Knepley for (k = 0; k < n; ++k) { 2206af0ca60SMatthew G. Knepley /* X[n,2] = [1, x] */ 2216af0ca60SMatthew G. Knepley X[k*2+0] = 1.0; 2226af0ca60SMatthew G. Knepley X[k*2+1] = x[k]; 2236af0ca60SMatthew G. Knepley } 2246af0ca60SMatthew G. Knepley /* H = X^T X */ 2256af0ca60SMatthew G. Knepley for (i = 0; i < 2; ++i) { 2266af0ca60SMatthew G. Knepley for (j = 0; j < 2; ++j) { 2276af0ca60SMatthew G. Knepley H[i*2+j] = 0.0; 2286af0ca60SMatthew G. Knepley for (k = 0; k < n; ++k) { 2296af0ca60SMatthew G. Knepley H[i*2+j] += X[k*2+i] * X[k*2+j]; 2306af0ca60SMatthew G. Knepley } 2316af0ca60SMatthew G. Knepley } 2326af0ca60SMatthew G. Knepley } 2336af0ca60SMatthew G. Knepley /* H = (X^T X)^{-1} */ 2346af0ca60SMatthew G. Knepley { 2356af0ca60SMatthew G. Knepley PetscBLASInt two = 2, ipiv[2], info; 2360955ed61SMatthew G. Knepley PetscScalar work[2]; 2376af0ca60SMatthew G. Knepley 2386af0ca60SMatthew G. Knepley ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 2396af0ca60SMatthew G. Knepley PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&two, &two, H, &two, ipiv, &info)); 2406af0ca60SMatthew G. Knepley PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&two, H, &two, ipiv, work, &two, &info)); 2416af0ca60SMatthew G. Knepley ierr = PetscFPTrapPop();CHKERRQ(ierr); 2426af0ca60SMatthew G. Knepley } 2436af0ca60SMatthew G. Knepley /* Y = H X^T */ 2446af0ca60SMatthew G. Knepley for (i = 0; i < 2; ++i) { 2456af0ca60SMatthew G. Knepley for (k = 0; k < n; ++k) { 2466af0ca60SMatthew G. Knepley Y[i*n+k] = 0.0; 2476af0ca60SMatthew G. Knepley for (j = 0; j < 2; ++j) { 24828cc7b48SMatthew G. Knepley Y[i*n+k] += PetscRealPart(H[i*2+j]) * X[k*2+j]; 2496af0ca60SMatthew G. Knepley } 2506af0ca60SMatthew G. Knepley } 2516af0ca60SMatthew G. Knepley } 2526af0ca60SMatthew G. Knepley /* beta = Y error = [y-intercept, slope] */ 2536af0ca60SMatthew G. Knepley for (i = 0; i < 2; ++i) { 2546af0ca60SMatthew G. Knepley beta[i] = 0.0; 2556af0ca60SMatthew G. Knepley for (k = 0; k < n; ++k) { 2566af0ca60SMatthew G. Knepley beta[i] += Y[i*n+k] * y[k]; 2576af0ca60SMatthew G. Knepley } 2586af0ca60SMatthew G. Knepley } 2596af0ca60SMatthew G. Knepley ierr = PetscFree2(X, Y);CHKERRQ(ierr); 2606af0ca60SMatthew G. Knepley *intercept = beta[0]; 2616af0ca60SMatthew G. Knepley *slope = beta[1]; 2626af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 2636af0ca60SMatthew G. Knepley } 2646af0ca60SMatthew G. Knepley 2656af0ca60SMatthew G. Knepley /*@ 2666af0ca60SMatthew G. Knepley PetscConvEstGetConvRate - Returns an estimate of the convergence rate for the discretization 2676af0ca60SMatthew G. Knepley 2686af0ca60SMatthew G. Knepley Not collective 2696af0ca60SMatthew G. Knepley 2706af0ca60SMatthew G. Knepley Input Parameter: 2716af0ca60SMatthew G. Knepley . ce - The PetscConvEst object 2726af0ca60SMatthew G. Knepley 2736af0ca60SMatthew G. Knepley Output Parameter: 274a56d3bf6SMatthew G. Knepley . alpha - The convergence rate for each field 2756af0ca60SMatthew G. Knepley 2766af0ca60SMatthew G. Knepley Note: The convergence rate alpha is defined by 2776af0ca60SMatthew G. Knepley $ || u_h - u_exact || < C h^alpha 2786af0ca60SMatthew G. Knepley where u_h is the discrete solution, and h is a measure of the discretization size. 2796af0ca60SMatthew G. Knepley 2806af0ca60SMatthew G. Knepley We solve a series of problems on refined meshes, calculate an error based upon the exact solution in the DS, 2816af0ca60SMatthew G. Knepley and then fit the result to our model above using linear regression. 2826af0ca60SMatthew G. Knepley 2836af0ca60SMatthew G. Knepley Options database keys: 2846af0ca60SMatthew G. Knepley . -snes_convergence_estimate : Execute convergence estimation and print out the rate 2856af0ca60SMatthew G. Knepley 2866af0ca60SMatthew G. Knepley Level: intermediate 2876af0ca60SMatthew G. Knepley 2886af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence 2896af0ca60SMatthew G. Knepley .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate() 2906af0ca60SMatthew G. Knepley @*/ 29146079b62SMatthew G. Knepley PetscErrorCode PetscConvEstGetConvRate(PetscConvEst ce, PetscReal alpha[]) 2926af0ca60SMatthew G. Knepley { 2936af0ca60SMatthew G. Knepley DM *dm; 2946af0ca60SMatthew G. Knepley PetscDS prob; 2956af0ca60SMatthew G. Knepley PetscObject disc; 2966af0ca60SMatthew G. Knepley MPI_Comm comm; 2976af0ca60SMatthew G. Knepley const char *uname, *dmname; 2986af0ca60SMatthew G. Knepley void *ctx; 2996af0ca60SMatthew G. Knepley Vec u; 3006af0ca60SMatthew G. Knepley PetscReal t = 0.0, *x, *y, slope, intercept; 30173269098SMatthew G. Knepley PetscInt *dof, dim, Nr = ce->Nr, r, f, oldlevel, oldnlev; 3025d68e14cSMatthew G. Knepley PetscLogEvent event; 3036af0ca60SMatthew G. Knepley PetscErrorCode ierr; 3046af0ca60SMatthew G. Knepley 3056af0ca60SMatthew G. Knepley PetscFunctionBegin; 3066af0ca60SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) ce, &comm);CHKERRQ(ierr); 3076af0ca60SMatthew G. Knepley ierr = DMGetDimension(ce->idm, &dim);CHKERRQ(ierr); 3086af0ca60SMatthew G. Knepley ierr = DMGetApplicationContext(ce->idm, &ctx);CHKERRQ(ierr); 3096af0ca60SMatthew G. Knepley ierr = DMGetDS(ce->idm, &prob);CHKERRQ(ierr); 3106af0ca60SMatthew G. Knepley ierr = DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE);CHKERRQ(ierr); 31173269098SMatthew G. Knepley ierr = DMGetRefineLevel(ce->idm, &oldlevel);CHKERRQ(ierr); 3120b755f55SMatthew G. Knepley ierr = PetscMalloc2((Nr+1), &dm, (Nr+1)*ce->Nf, &dof);CHKERRQ(ierr); 3136af0ca60SMatthew G. Knepley dm[0] = ce->idm; 3140b755f55SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) alpha[f] = 0.0; 3156af0ca60SMatthew G. Knepley /* Loop over meshes */ 3165d68e14cSMatthew G. Knepley ierr = PetscLogEventRegister("ConvEst Error", PETSC_OBJECT_CLASSID, &event);CHKERRQ(ierr); 3176af0ca60SMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 3185d68e14cSMatthew G. Knepley PetscLogStage stage; 3195d68e14cSMatthew G. Knepley char stageName[PETSC_MAX_PATH_LEN]; 3205d68e14cSMatthew G. Knepley 3215d68e14cSMatthew G. Knepley ierr = PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN-1, "ConvEst Refinement Level %D", r);CHKERRQ(ierr); 3225d68e14cSMatthew G. Knepley ierr = PetscLogStageRegister(stageName, &stage);CHKERRQ(ierr); 3235d68e14cSMatthew G. Knepley ierr = PetscLogStagePush(stage);CHKERRQ(ierr); 3246af0ca60SMatthew G. Knepley if (r > 0) { 3256af0ca60SMatthew G. Knepley ierr = DMRefine(dm[r-1], MPI_COMM_NULL, &dm[r]);CHKERRQ(ierr); 3266af0ca60SMatthew G. Knepley ierr = DMSetCoarseDM(dm[r], dm[r-1]);CHKERRQ(ierr); 3276af0ca60SMatthew G. Knepley ierr = DMSetDS(dm[r], prob);CHKERRQ(ierr); 3286af0ca60SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) dm[r-1], &dmname);CHKERRQ(ierr); 3296af0ca60SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) dm[r], dmname);CHKERRQ(ierr); 3308f892730SMatthew G. Knepley for (f = 0; f <= ce->Nf; ++f) { 3318f892730SMatthew G. Knepley PetscErrorCode (*nspconstr)(DM, PetscInt, MatNullSpace *); 3328f892730SMatthew G. Knepley ierr = DMGetNullSpaceConstructor(dm[r-1], f, &nspconstr);CHKERRQ(ierr); 3338f892730SMatthew G. Knepley ierr = DMSetNullSpaceConstructor(dm[r], f, nspconstr);CHKERRQ(ierr); 3348f892730SMatthew G. Knepley } 3356af0ca60SMatthew G. Knepley } 3366af0ca60SMatthew G. Knepley ierr = DMViewFromOptions(dm[r], NULL, "-conv_dm_view");CHKERRQ(ierr); 3376af0ca60SMatthew G. Knepley /* Create solution */ 3386af0ca60SMatthew G. Knepley ierr = DMCreateGlobalVector(dm[r], &u);CHKERRQ(ierr); 339239a0881SMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, 0, &disc);CHKERRQ(ierr); 340239a0881SMatthew G. Knepley ierr = PetscObjectGetName(disc, &uname);CHKERRQ(ierr); 3416af0ca60SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) u, uname);CHKERRQ(ierr); 3426af0ca60SMatthew G. Knepley /* Setup solver */ 3436af0ca60SMatthew G. Knepley ierr = SNESReset(ce->snes);CHKERRQ(ierr); 3446af0ca60SMatthew G. Knepley ierr = SNESSetDM(ce->snes, dm[r]);CHKERRQ(ierr); 3456af0ca60SMatthew G. Knepley ierr = DMPlexSetSNESLocalFEM(dm[r], ctx, ctx, ctx);CHKERRQ(ierr); 3466af0ca60SMatthew G. Knepley ierr = SNESSetFromOptions(ce->snes);CHKERRQ(ierr); 3476af0ca60SMatthew G. Knepley /* Create initial guess */ 3486af0ca60SMatthew G. Knepley ierr = DMProjectFunction(dm[r], t, ce->initGuess, NULL, INSERT_VALUES, u);CHKERRQ(ierr); 3496af0ca60SMatthew G. Knepley ierr = SNESSolve(ce->snes, NULL, u);CHKERRQ(ierr); 3505d68e14cSMatthew G. Knepley ierr = PetscLogEventBegin(event, ce, 0, 0, 0);CHKERRQ(ierr); 3516af0ca60SMatthew G. Knepley ierr = DMComputeL2FieldDiff(dm[r], t, ce->exactSol, NULL, u, &ce->errors[r*ce->Nf]);CHKERRQ(ierr); 3525d68e14cSMatthew G. Knepley ierr = PetscLogEventEnd(event, ce, 0, 0, 0);CHKERRQ(ierr); 353a56d3bf6SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 354a56d3bf6SMatthew G. Knepley PetscSection s, fs; 355a56d3bf6SMatthew G. Knepley PetscInt lsize; 356a56d3bf6SMatthew G. Knepley 357a56d3bf6SMatthew G. Knepley /* Could use DMGetOutputDM() to add in Dirichlet dofs */ 358a56d3bf6SMatthew G. Knepley ierr = DMGetSection(dm[r], &s);CHKERRQ(ierr); 359a56d3bf6SMatthew G. Knepley ierr = PetscSectionGetField(s, f, &fs);CHKERRQ(ierr); 360a56d3bf6SMatthew G. Knepley ierr = PetscSectionGetConstrainedStorageSize(fs, &lsize);CHKERRQ(ierr); 361a56d3bf6SMatthew G. Knepley ierr = MPI_Allreduce(&lsize, &dof[r*ce->Nf+f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) ce->snes));CHKERRQ(ierr); 362a56d3bf6SMatthew G. Knepley ierr = PetscLogEventSetDof(event, f, dof[r*ce->Nf+f]);CHKERRQ(ierr); 363a56d3bf6SMatthew G. Knepley ierr = PetscLogEventSetError(event, f, ce->errors[r*ce->Nf+f]);CHKERRQ(ierr); 364a56d3bf6SMatthew G. Knepley } 3656af0ca60SMatthew G. Knepley /* Monitor */ 3666af0ca60SMatthew G. Knepley if (ce->monitor) { 3676af0ca60SMatthew G. Knepley PetscReal *errors = &ce->errors[r*ce->Nf]; 3686af0ca60SMatthew G. Knepley 3690b755f55SMatthew G. Knepley ierr = PetscPrintf(comm, "L_2 Error: ");CHKERRQ(ierr); 3700b755f55SMatthew G. Knepley if (ce->Nf > 1) {ierr = PetscPrintf(comm, "[");CHKERRQ(ierr);} 3716af0ca60SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 3726af0ca60SMatthew G. Knepley if (f > 0) {ierr = PetscPrintf(comm, ", ");CHKERRQ(ierr);} 3736af0ca60SMatthew G. Knepley if (errors[f] < 1.0e-11) {ierr = PetscPrintf(comm, "< 1e-11");CHKERRQ(ierr);} 3748dfade3bSBarry Smith else {ierr = PetscPrintf(comm, "%g", (double)errors[f]);CHKERRQ(ierr);} 3756af0ca60SMatthew G. Knepley } 3760b755f55SMatthew G. Knepley if (ce->Nf > 1) {ierr = PetscPrintf(comm, "]");CHKERRQ(ierr);} 3770b755f55SMatthew G. Knepley ierr = PetscPrintf(comm, "\n");CHKERRQ(ierr); 3786af0ca60SMatthew G. Knepley } 37973269098SMatthew G. Knepley if (!r) { 38073269098SMatthew G. Knepley /* PCReset() does not wipe out the level structure */ 38173269098SMatthew G. Knepley KSP ksp; 38273269098SMatthew G. Knepley PC pc; 38373269098SMatthew G. Knepley 38473269098SMatthew G. Knepley ierr = SNESGetKSP(ce->snes, &ksp);CHKERRQ(ierr); 38573269098SMatthew G. Knepley ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 38673269098SMatthew G. Knepley ierr = PCMGGetLevels(pc, &oldnlev);CHKERRQ(ierr); 38773269098SMatthew G. Knepley } 3886af0ca60SMatthew G. Knepley /* Cleanup */ 3896af0ca60SMatthew G. Knepley ierr = VecDestroy(&u);CHKERRQ(ierr); 3905d68e14cSMatthew G. Knepley ierr = PetscLogStagePop();CHKERRQ(ierr); 3916af0ca60SMatthew G. Knepley } 3926af0ca60SMatthew G. Knepley for (r = 1; r <= Nr; ++r) { 3936af0ca60SMatthew G. Knepley ierr = DMDestroy(&dm[r]);CHKERRQ(ierr); 3946af0ca60SMatthew G. Knepley } 3956af0ca60SMatthew G. Knepley /* Fit convergence rate */ 3966af0ca60SMatthew G. Knepley ierr = PetscMalloc2(Nr+1, &x, Nr+1, &y);CHKERRQ(ierr); 397a56d3bf6SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 3986af0ca60SMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 399a56d3bf6SMatthew G. Knepley x[r] = PetscLog10Real(dof[r*ce->Nf+f]); 400a56d3bf6SMatthew G. Knepley y[r] = PetscLog10Real(ce->errors[r*ce->Nf+f]); 4016af0ca60SMatthew G. Knepley } 4026af0ca60SMatthew G. Knepley ierr = PetscConvEstLinearRegression_Private(ce, Nr+1, x, y, &slope, &intercept);CHKERRQ(ierr); 4036af0ca60SMatthew G. Knepley /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */ 404a56d3bf6SMatthew G. Knepley alpha[f] = -slope * dim; 405a56d3bf6SMatthew G. Knepley } 406a56d3bf6SMatthew G. Knepley ierr = PetscFree2(x, y);CHKERRQ(ierr); 4076af0ca60SMatthew G. Knepley ierr = PetscFree2(dm, dof);CHKERRQ(ierr); 4082cae373cSMatthew G. Knepley /* Restore solver */ 4092cae373cSMatthew G. Knepley ierr = SNESReset(ce->snes);CHKERRQ(ierr); 41073269098SMatthew G. Knepley { 41173269098SMatthew G. Knepley /* PCReset() does not wipe out the level structure */ 41273269098SMatthew G. Knepley KSP ksp; 41373269098SMatthew G. Knepley PC pc; 41473269098SMatthew G. Knepley 41573269098SMatthew G. Knepley ierr = SNESGetKSP(ce->snes, &ksp);CHKERRQ(ierr); 41673269098SMatthew G. Knepley ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 41773269098SMatthew G. Knepley ierr = PCMGSetLevels(pc, oldnlev, NULL);CHKERRQ(ierr); 41873269098SMatthew G. Knepley ierr = DMSetRefineLevel(ce->idm, oldlevel);CHKERRQ(ierr); /* The damn DMCoarsen() calls in PCMG can reset this */ 41973269098SMatthew G. Knepley } 4202cae373cSMatthew G. Knepley ierr = SNESSetDM(ce->snes, ce->idm);CHKERRQ(ierr); 4212cae373cSMatthew G. Knepley ierr = DMPlexSetSNESLocalFEM(ce->idm, ctx, ctx, ctx);CHKERRQ(ierr); 4222cae373cSMatthew G. Knepley ierr = SNESSetFromOptions(ce->snes);CHKERRQ(ierr); 4236af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 4246af0ca60SMatthew G. Knepley } 4256af0ca60SMatthew G. Knepley 4266af0ca60SMatthew G. Knepley /*@ 4276af0ca60SMatthew G. Knepley PetscConvEstRateView - Displays the convergence rate to a viewer 4286af0ca60SMatthew G. Knepley 4296af0ca60SMatthew G. Knepley Collective on SNES 4306af0ca60SMatthew G. Knepley 4316af0ca60SMatthew G. Knepley Parameter: 4326af0ca60SMatthew G. Knepley + snes - iterative context obtained from SNESCreate() 433a56d3bf6SMatthew G. Knepley . alpha - the convergence rate for each field 4346af0ca60SMatthew G. Knepley - viewer - the viewer to display the reason 4356af0ca60SMatthew G. Knepley 4366af0ca60SMatthew G. Knepley Options Database Keys: 4376af0ca60SMatthew G. Knepley . -snes_convergence_estimate - print the convergence rate 4386af0ca60SMatthew G. Knepley 4396af0ca60SMatthew G. Knepley Level: developer 4406af0ca60SMatthew G. Knepley 4416af0ca60SMatthew G. Knepley .seealso: PetscConvEstGetRate() 4426af0ca60SMatthew G. Knepley @*/ 443a56d3bf6SMatthew G. Knepley PetscErrorCode PetscConvEstRateView(PetscConvEst ce, const PetscReal alpha[], PetscViewer viewer) 4446af0ca60SMatthew G. Knepley { 4456af0ca60SMatthew G. Knepley PetscBool isAscii; 4466af0ca60SMatthew G. Knepley PetscErrorCode ierr; 4476af0ca60SMatthew G. Knepley 4486af0ca60SMatthew G. Knepley PetscFunctionBegin; 4496af0ca60SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr); 4506af0ca60SMatthew G. Knepley if (isAscii) { 451a56d3bf6SMatthew G. Knepley DM dm; 452a56d3bf6SMatthew G. Knepley PetscInt Nf, f; 453a56d3bf6SMatthew G. Knepley 454a56d3bf6SMatthew G. Knepley ierr = SNESGetDM(ce->snes, &dm);CHKERRQ(ierr); 455a56d3bf6SMatthew G. Knepley ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 4566af0ca60SMatthew G. Knepley ierr = PetscViewerASCIIAddTab(viewer, ((PetscObject) ce)->tablevel);CHKERRQ(ierr); 457a56d3bf6SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "L_2 convergence rate: ");CHKERRQ(ierr); 458a56d3bf6SMatthew G. Knepley if (Nf > 1) {ierr = PetscViewerASCIIPrintf(viewer, "[");CHKERRQ(ierr);} 459a56d3bf6SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 460a56d3bf6SMatthew G. Knepley if (f > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);} 461*5f9b3039SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "%#.2g", (double) alpha[f]);CHKERRQ(ierr); 462a56d3bf6SMatthew G. Knepley } 463a56d3bf6SMatthew G. Knepley if (Nf > 1) {ierr = PetscViewerASCIIPrintf(viewer, "]");CHKERRQ(ierr);} 464a56d3bf6SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 4656af0ca60SMatthew G. Knepley ierr = PetscViewerASCIISubtractTab(viewer, ((PetscObject) ce)->tablevel);CHKERRQ(ierr); 4666af0ca60SMatthew G. Knepley } 4676af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 4686af0ca60SMatthew G. Knepley } 469