Actual source code: convest.c
petsc-3.13.0 2020-03-29
1: #include <petscconvest.h>
2: #include <petscdmplex.h>
3: #include <petscds.h>
5: #include <petsc/private/petscconvestimpl.h>
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: }
14: /*@
15: PetscConvEstDestroy - Destroys a PetscConvEst object
17: Collective on PetscConvEst
19: Input Parameter:
20: . ce - The PetscConvEst object
22: Level: beginner
24: .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
25: @*/
26: PetscErrorCode PetscConvEstDestroy(PetscConvEst *ce)
27: {
31: if (!*ce) return(0);
33: if (--((PetscObject)(*ce))->refct > 0) {
34: *ce = NULL;
35: return(0);
36: }
37: PetscFree3((*ce)->initGuess, (*ce)->exactSol, (*ce)->ctxs);
38: PetscFree((*ce)->errors);
39: PetscHeaderDestroy(ce);
40: return(0);
41: }
43: /*@
44: PetscConvEstSetFromOptions - Sets a PetscConvEst object from options
46: Collective on PetscConvEst
48: Input Parameters:
49: . ce - The PetscConvEst object
51: Level: beginner
53: .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
54: @*/
55: PetscErrorCode PetscConvEstSetFromOptions(PetscConvEst ce)
56: {
60: PetscOptionsBegin(PetscObjectComm((PetscObject) ce), "", "Convergence Estimator Options", "PetscConvEst");
61: PetscOptionsInt("-convest_num_refine", "The number of refinements for the convergence check", "PetscConvEst", ce->Nr, &ce->Nr, NULL);
62: PetscOptionsReal("-convest_refine_factor", "The increase in resolution in each dimension", "PetscConvEst", ce->r, &ce->r, NULL);
63: PetscOptionsBool("-convest_monitor", "Monitor the error for each convergence check", "PetscConvEst", ce->monitor, &ce->monitor, NULL);
64: PetscOptionsEnd();
65: return(0);
66: }
68: /*@
69: PetscConvEstView - Views a PetscConvEst object
71: Collective on PetscConvEst
73: Input Parameters:
74: + ce - The PetscConvEst object
75: - viewer - The PetscViewer object
77: Level: beginner
79: .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
80: @*/
81: PetscErrorCode PetscConvEstView(PetscConvEst ce, PetscViewer viewer)
82: {
86: PetscObjectPrintClassNamePrefixType((PetscObject) ce, viewer);
87: PetscViewerASCIIPrintf(viewer, "ConvEst with %D levels\n", ce->Nr+1);
88: return(0);
89: }
91: /*@
92: PetscConvEstGetSolver - Gets the solver used to produce discrete solutions
94: Not collective
96: Input Parameter:
97: . ce - The PetscConvEst object
99: Output Parameter:
100: . solver - The solver
102: Level: intermediate
104: .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate()
105: @*/
106: PetscErrorCode PetscConvEstGetSolver(PetscConvEst ce, PetscObject *solver)
107: {
111: *solver = ce->solver;
112: return(0);
113: }
115: /*@
116: PetscConvEstSetSolver - Sets the solver used to produce discrete solutions
118: Not collective
120: Input Parameters:
121: + ce - The PetscConvEst object
122: - solver - The solver
124: Level: intermediate
126: Note: The solver MUST have an attached DM/DS, so that we know the exact solution
128: .seealso: PetscConvEstGetSNES(), PetscConvEstCreate(), PetscConvEstGetConvRate()
129: @*/
130: PetscErrorCode PetscConvEstSetSolver(PetscConvEst ce, PetscObject solver)
131: {
137: ce->solver = solver;
138: (*ce->ops->setsolver)(ce, solver);
139: return(0);
140: }
142: /*@
143: PetscConvEstSetUp - After the solver is specified, we create structures for estimating convergence
145: Collective on PetscConvEst
147: Input Parameters:
148: . ce - The PetscConvEst object
150: Level: beginner
152: .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
153: @*/
154: PetscErrorCode PetscConvEstSetUp(PetscConvEst ce)
155: {
156: PetscDS ds;
157: PetscInt Nf, f;
161: DMGetDS(ce->idm, &ds);
162: PetscDSGetNumFields(ds, &Nf);
163: ce->Nf = PetscMax(Nf, 1);
164: PetscMalloc1((ce->Nr+1)*ce->Nf, &ce->errors);
165: PetscMalloc3(ce->Nf, &ce->initGuess, ce->Nf, &ce->exactSol, ce->Nf, &ce->ctxs);
166: for (f = 0; f < Nf; ++f) ce->initGuess[f] = zero_private;
167: for (f = 0; f < Nf; ++f) {
168: PetscDSGetExactSolution(ds, f, &ce->exactSol[f], &ce->ctxs[f]);
169: 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);
170: }
171: return(0);
172: }
174: PetscErrorCode PetscConvEstComputeInitialGuess(PetscConvEst ce, PetscInt r, DM dm, Vec u)
175: {
182: (*ce->ops->initguess)(ce, r, dm, u);
183: return(0);
184: }
186: PetscErrorCode PetscConvEstComputeError(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[])
187: {
195: (*ce->ops->computeerror)(ce, r, dm, u, errors);
196: return(0);
197: }
199: static PetscErrorCode PetscConvEstMonitor_Private(PetscConvEst ce, PetscInt r)
200: {
201: MPI_Comm comm;
202: PetscInt f;
206: if (ce->monitor) {
207: PetscReal *errors = &ce->errors[r*ce->Nf];
209: PetscObjectGetComm((PetscObject) ce, &comm);
210: PetscPrintf(comm, "L_2 Error: ");
211: if (ce->Nf > 1) {PetscPrintf(comm, "[");}
212: for (f = 0; f < ce->Nf; ++f) {
213: if (f > 0) {PetscPrintf(comm, ", ");}
214: if (errors[f] < 1.0e-11) {PetscPrintf(comm, "< 1e-11");}
215: else {PetscPrintf(comm, "%g", (double) errors[f]);}
216: }
217: if (ce->Nf > 1) {PetscPrintf(comm, "]");}
218: PetscPrintf(comm, "\n");
219: }
220: return(0);
221: }
223: static PetscErrorCode PetscConvEstSetSNES_Private(PetscConvEst ce, PetscObject solver)
224: {
225: PetscClassId id;
229: PetscObjectGetClassId(ce->solver, &id);
230: if (id != SNES_CLASSID) SETERRQ(PetscObjectComm((PetscObject) ce), PETSC_ERR_ARG_WRONG, "Solver was not a SNES");
231: SNESGetDM((SNES) ce->solver, &ce->idm);
232: return(0);
233: }
235: static PetscErrorCode PetscConvEstInitGuessSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u)
236: {
240: DMProjectFunction(dm, 0.0, ce->initGuess, ce->ctxs, INSERT_VALUES, u);
241: return(0);
242: }
244: static PetscErrorCode PetscConvEstComputeErrorSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[])
245: {
249: DMComputeL2FieldDiff(dm, 0.0, ce->exactSol, ce->ctxs, u, errors);
250: return(0);
251: }
253: static PetscErrorCode PetscConvEstGetConvRateSNES_Private(PetscConvEst ce, PetscReal alpha[])
254: {
255: SNES snes = (SNES) ce->solver;
256: DM *dm;
257: PetscObject disc;
258: PetscReal *x, *y, slope, intercept;
259: PetscInt *dof, Nr = ce->Nr, r, f, dim, oldlevel, oldnlev;
260: void *ctx;
264: if (ce->r != 2.0) SETERRQ1(PetscObjectComm((PetscObject) ce), PETSC_ERR_SUP, "Only refinement factor 2 is currently supported (not %g)", (double) ce->r);
265: DMGetDimension(ce->idm, &dim);
266: DMGetApplicationContext(ce->idm, &ctx);
267: DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE);
268: DMGetRefineLevel(ce->idm, &oldlevel);
269: PetscMalloc2((Nr+1), &dm, (Nr+1)*ce->Nf, &dof);
270: /* Loop over meshes */
271: dm[0] = ce->idm;
272: for (r = 0; r <= Nr; ++r) {
273: Vec u;
274: PetscLogStage stage;
275: char stageName[PETSC_MAX_PATH_LEN];
276: const char *dmname, *uname;
278: PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN-1, "ConvEst Refinement Level %D", r);
279: PetscLogStageRegister(stageName, &stage);
280: PetscLogStagePush(stage);
281: if (r > 0) {
282: DMRefine(dm[r-1], MPI_COMM_NULL, &dm[r]);
283: DMSetCoarseDM(dm[r], dm[r-1]);
284: DMCopyDisc(ce->idm, dm[r]);
285: DMCopyTransform(ce->idm, dm[r]);
286: PetscObjectGetName((PetscObject) dm[r-1], &dmname);
287: PetscObjectSetName((PetscObject) dm[r], dmname);
288: for (f = 0; f <= ce->Nf; ++f) {
289: PetscErrorCode (*nspconstr)(DM, PetscInt, MatNullSpace *);
291: DMGetNullSpaceConstructor(dm[r-1], f, &nspconstr);
292: DMSetNullSpaceConstructor(dm[r], f, nspconstr);
293: }
294: }
295: DMViewFromOptions(dm[r], NULL, "-conv_dm_view");
296: /* Create solution */
297: DMCreateGlobalVector(dm[r], &u);
298: DMGetField(dm[r], 0, NULL, &disc);
299: PetscObjectGetName(disc, &uname);
300: PetscObjectSetName((PetscObject) u, uname);
301: /* Setup solver */
302: SNESReset(snes);
303: SNESSetDM(snes, dm[r]);
304: DMPlexSetSNESLocalFEM(dm[r], ctx, ctx, ctx);
305: SNESSetFromOptions(snes);
306: /* Create initial guess */
307: PetscConvEstComputeInitialGuess(ce, r, dm[r], u);
308: SNESSolve(snes, NULL, u);
309: PetscLogEventBegin(ce->event, ce, 0, 0, 0);
310: PetscConvEstComputeError(ce, r, dm[r], u, &ce->errors[r*ce->Nf]);
311: PetscLogEventEnd(ce->event, ce, 0, 0, 0);
312: for (f = 0; f < ce->Nf; ++f) {
313: PetscSection s, fs;
314: PetscInt lsize;
316: /* Could use DMGetOutputDM() to add in Dirichlet dofs */
317: DMGetLocalSection(dm[r], &s);
318: PetscSectionGetField(s, f, &fs);
319: PetscSectionGetConstrainedStorageSize(fs, &lsize);
320: MPI_Allreduce(&lsize, &dof[r*ce->Nf+f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) snes));
321: PetscLogEventSetDof(ce->event, f, dof[r*ce->Nf+f]);
322: PetscLogEventSetError(ce->event, f, ce->errors[r*ce->Nf+f]);
323: }
324: /* Monitor */
325: PetscConvEstMonitor_Private(ce, r);
326: if (!r) {
327: /* PCReset() does not wipe out the level structure */
328: KSP ksp;
329: PC pc;
331: SNESGetKSP(snes, &ksp);
332: KSPGetPC(ksp, &pc);
333: PCMGGetLevels(pc, &oldnlev);
334: }
335: /* Cleanup */
336: VecDestroy(&u);
337: PetscLogStagePop();
338: }
339: for (r = 1; r <= Nr; ++r) {
340: DMDestroy(&dm[r]);
341: }
342: /* Fit convergence rate */
343: PetscMalloc2(Nr+1, &x, Nr+1, &y);
344: for (f = 0; f < ce->Nf; ++f) {
345: for (r = 0; r <= Nr; ++r) {
346: x[r] = PetscLog10Real(dof[r*ce->Nf+f]);
347: y[r] = PetscLog10Real(ce->errors[r*ce->Nf+f]);
348: }
349: PetscLinearRegression(Nr+1, x, y, &slope, &intercept);
350: /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */
351: alpha[f] = -slope * dim;
352: }
353: PetscFree2(x, y);
354: PetscFree2(dm, dof);
355: /* Restore solver */
356: SNESReset(snes);
357: {
358: /* PCReset() does not wipe out the level structure */
359: KSP ksp;
360: PC pc;
362: SNESGetKSP(snes, &ksp);
363: KSPGetPC(ksp, &pc);
364: PCMGSetLevels(pc, oldnlev, NULL);
365: DMSetRefineLevel(ce->idm, oldlevel); /* The damn DMCoarsen() calls in PCMG can reset this */
366: }
367: SNESSetDM(snes, ce->idm);
368: DMPlexSetSNESLocalFEM(ce->idm, ctx, ctx, ctx);
369: SNESSetFromOptions(snes);
370: return(0);
371: }
373: /*@
374: PetscConvEstGetConvRate - Returns an estimate of the convergence rate for the discretization
376: Not collective
378: Input Parameter:
379: . ce - The PetscConvEst object
381: Output Parameter:
382: . alpha - The convergence rate for each field
384: Note: The convergence rate alpha is defined by
385: $ || u_\Delta - u_exact || < C \Delta^alpha
386: where u_\Delta is the discrete solution, and Delta is a measure of the discretization size. We usually use h for the
387: spatial resolution and \Delta t for the temporal resolution.
389: We solve a series of problems using increasing resolution (refined meshes or decreased timesteps), calculate an error
390: based upon the exact solution in the DS, and then fit the result to our model above using linear regression.
392: Options database keys:
393: + -snes_convergence_estimate : Execute convergence estimation inside SNESSolve() and print out the rate
394: - -ts_convergence_estimate : Execute convergence estimation inside TSSolve() and print out the rate
396: Level: intermediate
398: .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate(), SNESSolve(), TSSolve()
399: @*/
400: PetscErrorCode PetscConvEstGetConvRate(PetscConvEst ce, PetscReal alpha[])
401: {
402: PetscInt f;
406: if (ce->event < 0) {PetscLogEventRegister("ConvEst Error", PETSC_OBJECT_CLASSID, &ce->event);}
407: for (f = 0; f < ce->Nf; ++f) alpha[f] = 0.0;
408: (*ce->ops->getconvrate)(ce, alpha);
409: return(0);
410: }
412: /*@
413: PetscConvEstRateView - Displays the convergence rate to a viewer
415: Collective on SNES
417: Parameter:
418: + snes - iterative context obtained from SNESCreate()
419: . alpha - the convergence rate for each field
420: - viewer - the viewer to display the reason
422: Options Database Keys:
423: . -snes_convergence_estimate - print the convergence rate
425: Level: developer
427: .seealso: PetscConvEstGetRate()
428: @*/
429: PetscErrorCode PetscConvEstRateView(PetscConvEst ce, const PetscReal alpha[], PetscViewer viewer)
430: {
431: PetscBool isAscii;
435: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);
436: if (isAscii) {
437: PetscInt Nf = ce->Nf, f;
439: PetscViewerASCIIAddTab(viewer, ((PetscObject) ce)->tablevel);
440: PetscViewerASCIIPrintf(viewer, "L_2 convergence rate: ");
441: if (Nf > 1) {PetscViewerASCIIPrintf(viewer, "[");}
442: for (f = 0; f < Nf; ++f) {
443: if (f > 0) {PetscViewerASCIIPrintf(viewer, ", ");}
444: PetscViewerASCIIPrintf(viewer, "%#.2g", (double) alpha[f]);
445: }
446: if (Nf > 1) {PetscViewerASCIIPrintf(viewer, "]");}
447: PetscViewerASCIIPrintf(viewer, "\n");
448: PetscViewerASCIISubtractTab(viewer, ((PetscObject) ce)->tablevel);
449: }
450: return(0);
451: }
453: /*@
454: PetscConvEstCreate - Create a PetscConvEst object
456: Collective
458: Input Parameter:
459: . comm - The communicator for the PetscConvEst object
461: Output Parameter:
462: . ce - The PetscConvEst object
464: Level: beginner
466: .seealso: PetscConvEstDestroy(), PetscConvEstGetConvRate()
467: @*/
468: PetscErrorCode PetscConvEstCreate(MPI_Comm comm, PetscConvEst *ce)
469: {
474: PetscSysInitializePackage();
475: PetscHeaderCreate(*ce, PETSC_OBJECT_CLASSID, "PetscConvEst", "ConvergenceEstimator", "SNES", comm, PetscConvEstDestroy, PetscConvEstView);
476: (*ce)->monitor = PETSC_FALSE;
477: (*ce)->r = 2.0;
478: (*ce)->Nr = 4;
479: (*ce)->event = -1;
480: (*ce)->ops->setsolver = PetscConvEstSetSNES_Private;
481: (*ce)->ops->initguess = PetscConvEstInitGuessSNES_Private;
482: (*ce)->ops->computeerror = PetscConvEstComputeErrorSNES_Private;
483: (*ce)->ops->getconvrate = PetscConvEstGetConvRateSNES_Private;
484: return(0);
485: }