PetscErrorCode KSPMonitorSNESLGResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,PetscObject *objs) { SNES snes = (SNES) objs[0]; PetscDrawLG lg = (PetscDrawLG) objs[1]; PetscErrorCode ierr; PetscReal y[2]; Vec snes_solution,work1,work2; PetscFunctionBegin; if (rnorm > 0.0) y[0] = PetscLog10Real(rnorm); else y[0] = -15.0; ierr = SNESGetSolution(snes,&snes_solution);CHKERRQ(ierr); ierr = VecDuplicate(snes_solution,&work1);CHKERRQ(ierr); ierr = VecDuplicate(snes_solution,&work2);CHKERRQ(ierr); ierr = KSPBuildSolution(ksp,work1,NULL);CHKERRQ(ierr); ierr = VecAYPX(work1,-1.0,snes_solution);CHKERRQ(ierr); ierr = SNESComputeFunction(snes,work1,work2);CHKERRQ(ierr); ierr = VecNorm(work2,NORM_2,y+1);CHKERRQ(ierr); if (y[1] > 0.0) y[1] = PetscLog10Real(y[1]); else y[1] = -15.0; ierr = VecDestroy(&work1);CHKERRQ(ierr); ierr = VecDestroy(&work2);CHKERRQ(ierr); ierr = PetscDrawLGAddPoint(lg,NULL,y);CHKERRQ(ierr); if (n < 20 || !(n % 5)) { ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); } PetscFunctionReturn(0); }
PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,PetscObject *objs) { PetscDrawLG lg = (PetscDrawLG) objs[0]; PetscReal x[2],y[2],scnorm; PetscErrorCode ierr; PetscMPIInt rank; Vec resid,work; PetscFunctionBegin; ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)ksp),&rank);CHKERRQ(ierr); if (!rank) { if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} x[0] = x[1] = (PetscReal) n; if (rnorm > 0.0) y[0] = PetscLog10Real(rnorm); else y[0] = -15.0; } ierr = VecDuplicate(ksp->vec_rhs,&work);CHKERRQ(ierr); ierr = KSPBuildResidual(ksp,0,work,&resid);CHKERRQ(ierr); ierr = VecNorm(resid,NORM_2,&scnorm);CHKERRQ(ierr); ierr = VecDestroy(&work);CHKERRQ(ierr); if (!rank) { if (scnorm > 0.0) y[1] = PetscLog10Real(scnorm); else y[1] = -15.0; ierr = PetscDrawLGAddPoint(lg,x,y);CHKERRQ(ierr); if (n <= 20 || (n % 3)) { ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); } } PetscFunctionReturn(0); }
PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,void *ctx) { PetscDrawLG lg = (PetscDrawLG) ctx; PetscReal x[2],y[2],scnorm; Vec resid,work; PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(ksp,KSP_CLASSID,1); PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,4); if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} x[0] = x[1] = (PetscReal) n; if (rnorm > 0.0) y[0] = PetscLog10Real(rnorm); else y[0] = -15.0; ierr = VecDuplicate(ksp->vec_rhs,&work);CHKERRQ(ierr); ierr = KSPBuildResidual(ksp,NULL,work,&resid);CHKERRQ(ierr); ierr = VecNorm(resid,NORM_2,&scnorm);CHKERRQ(ierr); ierr = VecDestroy(&work);CHKERRQ(ierr); if (scnorm > 0.0) y[1] = PetscLog10Real(scnorm); else y[1] = -15.0; ierr = PetscDrawLGAddPoint(lg,x,y);CHKERRQ(ierr); if (n <= 20 || !(n % 5)) { ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); } PetscFunctionReturn(0); }
PetscErrorCode PetscAGetBase(PetscReal vmin,PetscReal vmax,int num,PetscReal *Base,int *power) { PetscReal base,ftemp,e10; static PetscReal base_try[5] = {10.0,5.0,2.0,1.0,0.5}; PetscErrorCode ierr; int i; PetscFunctionBegin; /* labels of the form n * BASE */ /* get an approximate value for BASE */ base = (vmax - vmin) / (double)(num + 1); /* make it of form m x 10^power, m in [1.0, 10) */ if (base <= 0.0) { base = PetscAbsReal(vmin); if (base < 1.0) base = 1.0; } ftemp = PetscLog10Real((1.0 + EPS) * base); if (ftemp < 0.0) ftemp -= 1.0; *power = (int)ftemp; ierr = PetscExp10((double)-*power,&e10);CHKERRQ(ierr); base = base * e10; if (base < 1.0) base = 1.0; /* now reduce it to one of 1, 2, or 5 */ for (i=1; i<5; i++) { if (base >= base_try[i]) { ierr = PetscExp10((double)*power,&e10);CHKERRQ(ierr); base = base_try[i-1] * e10; if (i == 1) *power = *power + 1; break; } } *Base = base; PetscFunctionReturn(0); }
PetscErrorCode KSPMonitorLGResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,PetscObject *objs) { PetscErrorCode ierr; PetscReal x,y; PetscDrawLG lg = (PetscDrawLG) objs[0]; PetscFunctionBegin; if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} x = (PetscReal) n; if (rnorm > 0.0) y = PetscLog10Real(rnorm); else y = -15.0; ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); if (n < 20 || !(n % 5)) { ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); } PetscFunctionReturn(0); }
PetscErrorCode KSPMonitorLGResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,void *ctx) { PetscDrawLG lg = (PetscDrawLG) ctx; PetscReal x,y; PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,4); if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} x = (PetscReal) n; if (rnorm > 0.0) y = PetscLog10Real(rnorm); else y = -15.0; ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); if (n <= 20 || !(n % 5)) { ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); } PetscFunctionReturn(0); }
PetscErrorCode KSPMonitorLGRange(KSP ksp,PetscInt n,PetscReal rnorm,void *monctx) { PetscDrawLG lg; PetscErrorCode ierr; PetscReal x,y,per; PetscViewer v = (PetscViewer)monctx; static PetscReal prev; /* should be in the context */ PetscDraw draw; PetscFunctionBegin; ierr = PetscViewerDrawGetDrawLG(v,0,&lg);CHKERRQ(ierr); if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); ierr = PetscDrawSetTitle(draw,"Residual norm");CHKERRQ(ierr); x = (PetscReal) n; if (rnorm > 0.0) y = PetscLog10Real(rnorm); else y = -15.0; ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); if (n < 20 || !(n % 5)) { ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); } ierr = PetscViewerDrawGetDrawLG(v,1,&lg);CHKERRQ(ierr); ierr = KSPMonitorRange_Private(ksp,n,&per);CHKERRQ(ierr); if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); ierr = PetscDrawSetTitle(draw,"% elemts > .2*max elemt");CHKERRQ(ierr); x = (PetscReal) n; y = 100.0*per; ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); if (n < 20 || !(n % 5)) { ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); } ierr = PetscViewerDrawGetDrawLG(v,2,&lg);CHKERRQ(ierr); if (!n) {prev = rnorm;ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");CHKERRQ(ierr); ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");CHKERRQ(ierr); x = (PetscReal) n; y = (prev - rnorm)/prev; ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); if (n < 20 || !(n % 5)) { ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); } ierr = PetscViewerDrawGetDrawLG(v,3,&lg);CHKERRQ(ierr); if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr); ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");CHKERRQ(ierr); x = (PetscReal) n; y = (prev - rnorm)/(prev*per); if (n > 5) { ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); } if (n < 20 || !(n % 5)) { ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); } prev = rnorm; PetscFunctionReturn(0); }
/*@ PetscConvEstGetConvRate - Returns an estimate of the convergence rate for the discretization Not collective Input Parameter: . ce - The PetscConvEst object Output Parameter: . alpha - The convergence rate for each field Note: The convergence rate alpha is defined by $ || u_h - u_exact || < C h^alpha where u_h is the discrete solution, and h is a measure of the discretization size. We solve a series of problems on refined meshes, calculate an error based upon the exact solution in the DS, and then fit the result to our model above using linear regression. Options database keys: . -snes_convergence_estimate : Execute convergence estimation and print out the rate Level: intermediate .keywords: PetscConvEst, convergence .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate() @*/ PetscErrorCode PetscConvEstGetConvRate(PetscConvEst ce, PetscReal alpha[]) { DM *dm; PetscObject disc; MPI_Comm comm; const char *uname, *dmname; void *ctx; Vec u; PetscReal t = 0.0, *x, *y, slope, intercept; PetscInt *dof, dim, Nr = ce->Nr, r, f, oldlevel, oldnlev; PetscLogEvent event; PetscErrorCode ierr; PetscFunctionBegin; ierr = PetscObjectGetComm((PetscObject) ce, &comm);CHKERRQ(ierr); ierr = DMGetDimension(ce->idm, &dim);CHKERRQ(ierr); ierr = DMGetApplicationContext(ce->idm, &ctx);CHKERRQ(ierr); ierr = DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE);CHKERRQ(ierr); ierr = DMGetRefineLevel(ce->idm, &oldlevel);CHKERRQ(ierr); ierr = PetscMalloc2((Nr+1), &dm, (Nr+1)*ce->Nf, &dof);CHKERRQ(ierr); dm[0] = ce->idm; for (f = 0; f < ce->Nf; ++f) alpha[f] = 0.0; /* Loop over meshes */ ierr = PetscLogEventRegister("ConvEst Error", PETSC_OBJECT_CLASSID, &event);CHKERRQ(ierr); for (r = 0; r <= Nr; ++r) { PetscLogStage stage; char stageName[PETSC_MAX_PATH_LEN]; ierr = PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN-1, "ConvEst Refinement Level %D", r);CHKERRQ(ierr); ierr = PetscLogStageRegister(stageName, &stage);CHKERRQ(ierr); ierr = PetscLogStagePush(stage);CHKERRQ(ierr); if (r > 0) { ierr = DMRefine(dm[r-1], MPI_COMM_NULL, &dm[r]);CHKERRQ(ierr); ierr = DMSetCoarseDM(dm[r], dm[r-1]);CHKERRQ(ierr); ierr = DMCopyDisc(ce->idm, dm[r]);CHKERRQ(ierr); ierr = DMCopyTransform(ce->idm, dm[r]);CHKERRQ(ierr); ierr = PetscObjectGetName((PetscObject) dm[r-1], &dmname);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject) dm[r], dmname);CHKERRQ(ierr); for (f = 0; f <= ce->Nf; ++f) { PetscErrorCode (*nspconstr)(DM, PetscInt, MatNullSpace *); ierr = DMGetNullSpaceConstructor(dm[r-1], f, &nspconstr);CHKERRQ(ierr); ierr = DMSetNullSpaceConstructor(dm[r], f, nspconstr);CHKERRQ(ierr); } } ierr = DMViewFromOptions(dm[r], NULL, "-conv_dm_view");CHKERRQ(ierr); /* Create solution */ ierr = DMCreateGlobalVector(dm[r], &u);CHKERRQ(ierr); ierr = DMGetField(dm[r], 0, NULL, &disc);CHKERRQ(ierr); ierr = PetscObjectGetName(disc, &uname);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject) u, uname);CHKERRQ(ierr); /* Setup solver */ ierr = SNESReset(ce->snes);CHKERRQ(ierr); ierr = SNESSetDM(ce->snes, dm[r]);CHKERRQ(ierr); ierr = DMPlexSetSNESLocalFEM(dm[r], ctx, ctx, ctx);CHKERRQ(ierr); ierr = SNESSetFromOptions(ce->snes);CHKERRQ(ierr); /* Create initial guess */ ierr = DMProjectFunction(dm[r], t, ce->initGuess, ce->ctxs, INSERT_VALUES, u);CHKERRQ(ierr); ierr = SNESSolve(ce->snes, NULL, u);CHKERRQ(ierr); ierr = PetscLogEventBegin(event, ce, 0, 0, 0);CHKERRQ(ierr); ierr = DMComputeL2FieldDiff(dm[r], t, ce->exactSol, ce->ctxs, u, &ce->errors[r*ce->Nf]);CHKERRQ(ierr); ierr = PetscLogEventEnd(event, ce, 0, 0, 0);CHKERRQ(ierr); for (f = 0; f < ce->Nf; ++f) { PetscSection s, fs; PetscInt lsize; /* Could use DMGetOutputDM() to add in Dirichlet dofs */ ierr = DMGetSection(dm[r], &s);CHKERRQ(ierr); ierr = PetscSectionGetField(s, f, &fs);CHKERRQ(ierr); ierr = PetscSectionGetConstrainedStorageSize(fs, &lsize);CHKERRQ(ierr); ierr = MPI_Allreduce(&lsize, &dof[r*ce->Nf+f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) ce->snes));CHKERRQ(ierr); ierr = PetscLogEventSetDof(event, f, dof[r*ce->Nf+f]);CHKERRQ(ierr); ierr = PetscLogEventSetError(event, f, ce->errors[r*ce->Nf+f]);CHKERRQ(ierr); } /* Monitor */ if (ce->monitor) { PetscReal *errors = &ce->errors[r*ce->Nf]; ierr = PetscPrintf(comm, "L_2 Error: ");CHKERRQ(ierr); if (ce->Nf > 1) {ierr = PetscPrintf(comm, "[");CHKERRQ(ierr);} for (f = 0; f < ce->Nf; ++f) { if (f > 0) {ierr = PetscPrintf(comm, ", ");CHKERRQ(ierr);} if (errors[f] < 1.0e-11) {ierr = PetscPrintf(comm, "< 1e-11");CHKERRQ(ierr);} else {ierr = PetscPrintf(comm, "%g", (double)errors[f]);CHKERRQ(ierr);} } if (ce->Nf > 1) {ierr = PetscPrintf(comm, "]");CHKERRQ(ierr);} ierr = PetscPrintf(comm, "\n");CHKERRQ(ierr); } if (!r) { /* PCReset() does not wipe out the level structure */ KSP ksp; PC pc; ierr = SNESGetKSP(ce->snes, &ksp);CHKERRQ(ierr); ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); ierr = PCMGGetLevels(pc, &oldnlev);CHKERRQ(ierr); } /* Cleanup */ ierr = VecDestroy(&u);CHKERRQ(ierr); ierr = PetscLogStagePop();CHKERRQ(ierr); } for (r = 1; r <= Nr; ++r) { ierr = DMDestroy(&dm[r]);CHKERRQ(ierr); } /* Fit convergence rate */ ierr = PetscMalloc2(Nr+1, &x, Nr+1, &y);CHKERRQ(ierr); for (f = 0; f < ce->Nf; ++f) { for (r = 0; r <= Nr; ++r) { x[r] = PetscLog10Real(dof[r*ce->Nf+f]); y[r] = PetscLog10Real(ce->errors[r*ce->Nf+f]); } ierr = PetscLinearRegression(Nr+1, x, y, &slope, &intercept);CHKERRQ(ierr); /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */ alpha[f] = -slope * dim; } ierr = PetscFree2(x, y);CHKERRQ(ierr); ierr = PetscFree2(dm, dof);CHKERRQ(ierr); /* Restore solver */ ierr = SNESReset(ce->snes);CHKERRQ(ierr); { /* PCReset() does not wipe out the level structure */ KSP ksp; PC pc; ierr = SNESGetKSP(ce->snes, &ksp);CHKERRQ(ierr); ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); ierr = PCMGSetLevels(pc, oldnlev, NULL);CHKERRQ(ierr); ierr = DMSetRefineLevel(ce->idm, oldlevel);CHKERRQ(ierr); /* The damn DMCoarsen() calls in PCMG can reset this */ } ierr = SNESSetDM(ce->snes, ce->idm);CHKERRQ(ierr); ierr = DMPlexSetSNESLocalFEM(ce->idm, ctx, ctx, ctx);CHKERRQ(ierr); ierr = SNESSetFromOptions(ce->snes);CHKERRQ(ierr); PetscFunctionReturn(0); }
/* Run with -build_twosided allreduce -pc_type bjacobi -sub_pc_type lu -q 16 -ksp_rtol 1.e-34 (or 1.e-14 for double precision) -q <q> number of spectral elements to use -N <N> maximum number of GLL points per element */ int main(int argc,char **args) { PetscErrorCode ierr; PetscGLL gll; PetscInt N = 80,n,q = 8,xs,xn,j,l; PetscReal **A; Mat K; KSP ksp; PC pc; Vec x,b; PetscInt *rows; PetscReal norm,xc,yc,h; PetscScalar *f; PetscDraw draw; PetscDrawLG lg; PetscDrawAxis axis; DM da; PetscMPIInt rank,size; ierr = PetscInitialize(&argc,&args,NULL,NULL);if (ierr) return ierr; ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-q",&q,NULL);CHKERRQ(ierr); ierr = PetscDrawCreate(PETSC_COMM_WORLD,NULL,"Log(Error norm) vs Number of GLL points",0,0,500,500,&draw);CHKERRQ(ierr); ierr = PetscDrawSetFromOptions(draw);CHKERRQ(ierr); ierr = PetscDrawLGCreate(draw,1,&lg);CHKERRQ(ierr); ierr = PetscDrawLGSetUseMarkers(lg,PETSC_TRUE);CHKERRQ(ierr); ierr = PetscDrawLGGetAxis(lg,&axis);CHKERRQ(ierr); ierr = PetscDrawAxisSetLabels(axis,NULL,"Number of GLL points","Log(Error Norm)");CHKERRQ(ierr); for (n=4; n<N; n+=2) { /* da contains the information about the parallel layout of the elements */ ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,q*(n-1)+1,1,1,NULL,&da);CHKERRQ(ierr); ierr = DMSetFromOptions(da);CHKERRQ(ierr); ierr = DMSetUp(da);CHKERRQ(ierr); ierr = DMDAGetInfo(da,NULL,&q,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); q = (q-1)/(n-1); /* number of spectral elements */ /* gll simply contains the GLL node and weight values */ ierr = PetscGLLCreate(n,PETSCGLL_VIA_LINEARALGEBRA,&gll);CHKERRQ(ierr); ierr = DMDASetGLLCoordinates(da,&gll);CHKERRQ(ierr); /* Creates the element stiffness matrix for the given gll */ ierr = PetscGLLElementLaplacianCreate(&gll,&A);CHKERRQ(ierr); /* Scale the element stiffness and weights by the size of the element */ h = 2.0/q; for (j=0; j<n; j++) { gll.weights[j] *= .5*h; for (l=0; l<n; l++) { A[j][l] = 2.*A[j][l]/h; } } /* Create the global stiffness matrix and add the element stiffness for each local element */ ierr = DMCreateMatrix(da,&K);CHKERRQ(ierr); ierr = MatSetOption(K,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xn,NULL,NULL);CHKERRQ(ierr); xs = xs/(n-1); xn = xn/(n-1); ierr = PetscMalloc1(n,&rows);CHKERRQ(ierr); /* loop over local elements */ for (j=xs; j<xs+xn; j++) { for (l=0; l<n; l++) rows[l] = j*(n-1)+l; ierr = MatSetValues(K,n,rows,n,rows,&A[0][0],ADD_VALUES);CHKERRQ(ierr); } ierr = MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatCreateVecs(K,&x,&b);CHKERRQ(ierr); ierr = ComputeRhs(da,&gll,b);CHKERRQ(ierr); /* Replace the first and last rows/columns of the matrix with the identity to obtain the zero Dirichlet boundary conditions */ rows[0] = 0; rows[1] = q*(n-1); ierr = MatZeroRowsColumns(K,2,rows,1.0,x,b);CHKERRQ(ierr); ierr = PetscFree(rows);CHKERRQ(ierr); ierr = KSPCreate(PETSC_COMM_WORLD,&ksp);CHKERRQ(ierr); ierr = KSPSetOperators(ksp,K,K);CHKERRQ(ierr); ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); ierr = PCSetType(pc,PCLU);CHKERRQ(ierr); ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); ierr = KSPSolve(ksp,b,x);CHKERRQ(ierr); /* compute the error to the continium problem */ ierr = ComputeSolution(da,&gll,b);CHKERRQ(ierr); ierr = VecAXPY(x,-1.0,b);CHKERRQ(ierr); /* compute the L^2 norm of the error */ ierr = VecGetArray(x,&f);CHKERRQ(ierr); ierr = PetscGLLIntegrate(&gll,f,&norm);CHKERRQ(ierr); ierr = VecRestoreArray(x,&f);CHKERRQ(ierr); norm = PetscSqrtReal(norm); ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"L^2 norm of the error %D %g\n",n,(double)norm);CHKERRQ(ierr); if (n > 10 && norm > 1.e-8) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Slower convergence than expected"); xc = (PetscReal)n; yc = PetscLog10Real(norm); ierr = PetscDrawLGAddPoint(lg,&xc,&yc);CHKERRQ(ierr); ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); ierr = VecDestroy(&b);CHKERRQ(ierr); ierr = VecDestroy(&x);CHKERRQ(ierr); ierr = KSPDestroy(&ksp);CHKERRQ(ierr); ierr = MatDestroy(&K);CHKERRQ(ierr); ierr = PetscGLLElementLaplacianDestroy(&gll,&A);CHKERRQ(ierr); ierr = PetscGLLDestroy(&gll);CHKERRQ(ierr); ierr = DMDestroy(&da);CHKERRQ(ierr); } ierr = PetscDrawLGDestroy(&lg);CHKERRQ(ierr); ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr); ierr = PetscFinalize(); return ierr; }