Actual source code: convest.c
petsc-3.14.6 2021-03-30
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: PetscInt Nf, f, Nds, s;
160: DMGetNumFields(ce->idm, &Nf);
161: ce->Nf = PetscMax(Nf, 1);
162: PetscMalloc1((ce->Nr+1)*ce->Nf, &ce->errors);
163: PetscCalloc3(ce->Nf, &ce->initGuess, ce->Nf, &ce->exactSol, ce->Nf, &ce->ctxs);
164: for (f = 0; f < Nf; ++f) ce->initGuess[f] = zero_private;
165: DMGetNumDS(ce->idm, &Nds);
166: for (s = 0; s < Nds; ++s) {
167: PetscDS ds;
168: DMLabel label;
169: IS fieldIS;
170: const PetscInt *fields;
171: PetscInt dsNf;
173: DMGetRegionNumDS(ce->idm, s, &label, &fieldIS, &ds);
174: PetscDSGetNumFields(ds, &dsNf);
175: if (fieldIS) {ISGetIndices(fieldIS, &fields);}
176: for (f = 0; f < dsNf; ++f) {
177: const PetscInt field = fields[f];
178: PetscDSGetExactSolution(ds, field, &ce->exactSol[field], &ce->ctxs[field]);
179: }
180: if (fieldIS) {ISRestoreIndices(fieldIS, &fields);}
181: }
182: for (f = 0; f < Nf; ++f) {
183: 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);
184: }
185: return(0);
186: }
188: PetscErrorCode PetscConvEstComputeInitialGuess(PetscConvEst ce, PetscInt r, DM dm, Vec u)
189: {
196: (*ce->ops->initguess)(ce, r, dm, u);
197: return(0);
198: }
200: PetscErrorCode PetscConvEstComputeError(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[])
201: {
209: (*ce->ops->computeerror)(ce, r, dm, u, errors);
210: return(0);
211: }
213: /*@
214: PetscConvEstMonitorDefault - Monitors the convergence estimation loop
216: Collective on PetscConvEst
218: Input Parameter:
219: + ce - The PetscConvEst object
220: - r - The refinement level
222: Options database keys:
223: . -convest_monitor - Activate the monitor
225: Level: intermediate
227: .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate(), SNESSolve(), TSSolve()
228: @*/
229: PetscErrorCode PetscConvEstMonitorDefault(PetscConvEst ce, PetscInt r)
230: {
231: MPI_Comm comm;
232: PetscInt f;
236: if (ce->monitor) {
237: PetscReal *errors = &ce->errors[r*ce->Nf];
239: PetscObjectGetComm((PetscObject) ce, &comm);
240: PetscPrintf(comm, "L_2 Error: ");
241: if (ce->Nf > 1) {PetscPrintf(comm, "[");}
242: for (f = 0; f < ce->Nf; ++f) {
243: if (f > 0) {PetscPrintf(comm, ", ");}
244: if (errors[f] < 1.0e-11) {PetscPrintf(comm, "< 1e-11");}
245: else {PetscPrintf(comm, "%g", (double) errors[f]);}
246: }
247: if (ce->Nf > 1) {PetscPrintf(comm, "]");}
248: PetscPrintf(comm, "\n");
249: }
250: return(0);
251: }
253: static PetscErrorCode PetscConvEstSetSNES_Private(PetscConvEst ce, PetscObject solver)
254: {
255: PetscClassId id;
259: PetscObjectGetClassId(ce->solver, &id);
260: if (id != SNES_CLASSID) SETERRQ(PetscObjectComm((PetscObject) ce), PETSC_ERR_ARG_WRONG, "Solver was not a SNES");
261: SNESGetDM((SNES) ce->solver, &ce->idm);
262: return(0);
263: }
265: static PetscErrorCode PetscConvEstInitGuessSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u)
266: {
270: DMProjectFunction(dm, 0.0, ce->initGuess, ce->ctxs, INSERT_VALUES, u);
271: return(0);
272: }
274: static PetscErrorCode PetscConvEstComputeErrorSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[])
275: {
279: DMComputeL2FieldDiff(dm, 0.0, ce->exactSol, ce->ctxs, u, errors);
280: return(0);
281: }
283: static PetscErrorCode PetscConvEstGetConvRateSNES_Private(PetscConvEst ce, PetscReal alpha[])
284: {
285: SNES snes = (SNES) ce->solver;
286: DM *dm;
287: PetscObject disc;
288: PetscReal *x, *y, slope, intercept;
289: PetscInt *dof, Nr = ce->Nr, r, f, dim, oldlevel, oldnlev;
290: void *ctx;
294: if (ce->r != 2.0) SETERRQ1(PetscObjectComm((PetscObject) ce), PETSC_ERR_SUP, "Only refinement factor 2 is currently supported (not %g)", (double) ce->r);
295: DMGetDimension(ce->idm, &dim);
296: DMGetApplicationContext(ce->idm, &ctx);
297: DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE);
298: DMGetRefineLevel(ce->idm, &oldlevel);
299: PetscMalloc2((Nr+1), &dm, (Nr+1)*ce->Nf, &dof);
300: /* Loop over meshes */
301: dm[0] = ce->idm;
302: for (r = 0; r <= Nr; ++r) {
303: Vec u;
304: #if defined(PETSC_USE_LOG)
305: PetscLogStage stage;
306: #endif
307: char stageName[PETSC_MAX_PATH_LEN];
308: const char *dmname, *uname;
310: PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN-1, "ConvEst Refinement Level %D", r);
311: PetscLogStageRegister(stageName, &stage);
312: PetscLogStagePush(stage);
313: if (r > 0) {
314: DMRefine(dm[r-1], MPI_COMM_NULL, &dm[r]);
315: DMSetCoarseDM(dm[r], dm[r-1]);
316: DMCopyTransform(ce->idm, dm[r]);
317: PetscObjectGetName((PetscObject) dm[r-1], &dmname);
318: PetscObjectSetName((PetscObject) dm[r], dmname);
319: for (f = 0; f <= ce->Nf; ++f) {
320: PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *);
322: DMGetNullSpaceConstructor(dm[r-1], f, &nspconstr);
323: DMSetNullSpaceConstructor(dm[r], f, nspconstr);
324: }
325: }
326: DMViewFromOptions(dm[r], NULL, "-conv_dm_view");
327: /* Create solution */
328: DMCreateGlobalVector(dm[r], &u);
329: DMGetField(dm[r], 0, NULL, &disc);
330: PetscObjectGetName(disc, &uname);
331: PetscObjectSetName((PetscObject) u, uname);
332: /* Setup solver */
333: SNESReset(snes);
334: SNESSetDM(snes, dm[r]);
335: DMPlexSetSNESLocalFEM(dm[r], ctx, ctx, ctx);
336: SNESSetFromOptions(snes);
337: /* Create initial guess */
338: PetscConvEstComputeInitialGuess(ce, r, dm[r], u);
339: SNESSolve(snes, NULL, u);
340: PetscLogEventBegin(ce->event, ce, 0, 0, 0);
341: PetscConvEstComputeError(ce, r, dm[r], u, &ce->errors[r*ce->Nf]);
342: PetscLogEventEnd(ce->event, ce, 0, 0, 0);
343: for (f = 0; f < ce->Nf; ++f) {
344: PetscSection s, fs;
345: PetscInt lsize;
347: /* Could use DMGetOutputDM() to add in Dirichlet dofs */
348: DMGetLocalSection(dm[r], &s);
349: PetscSectionGetField(s, f, &fs);
350: PetscSectionGetConstrainedStorageSize(fs, &lsize);
351: MPI_Allreduce(&lsize, &dof[r*ce->Nf+f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) snes));
352: PetscLogEventSetDof(ce->event, f, dof[r*ce->Nf+f]);
353: PetscLogEventSetError(ce->event, f, ce->errors[r*ce->Nf+f]);
354: }
355: /* Monitor */
356: PetscConvEstMonitorDefault(ce, r);
357: if (!r) {
358: /* PCReset() does not wipe out the level structure */
359: KSP ksp;
360: PC pc;
362: SNESGetKSP(snes, &ksp);
363: KSPGetPC(ksp, &pc);
364: PCMGGetLevels(pc, &oldnlev);
365: }
366: /* Cleanup */
367: VecDestroy(&u);
368: PetscLogStagePop();
369: }
370: for (r = 1; r <= Nr; ++r) {
371: DMDestroy(&dm[r]);
372: }
373: /* Fit convergence rate */
374: PetscMalloc2(Nr+1, &x, Nr+1, &y);
375: for (f = 0; f < ce->Nf; ++f) {
376: for (r = 0; r <= Nr; ++r) {
377: x[r] = PetscLog10Real(dof[r*ce->Nf+f]);
378: y[r] = PetscLog10Real(ce->errors[r*ce->Nf+f]);
379: }
380: PetscLinearRegression(Nr+1, x, y, &slope, &intercept);
381: /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */
382: alpha[f] = -slope * dim;
383: }
384: PetscFree2(x, y);
385: PetscFree2(dm, dof);
386: /* Restore solver */
387: SNESReset(snes);
388: {
389: /* PCReset() does not wipe out the level structure */
390: KSP ksp;
391: PC pc;
393: SNESGetKSP(snes, &ksp);
394: KSPGetPC(ksp, &pc);
395: PCMGSetLevels(pc, oldnlev, NULL);
396: DMSetRefineLevel(ce->idm, oldlevel); /* The damn DMCoarsen() calls in PCMG can reset this */
397: }
398: SNESSetDM(snes, ce->idm);
399: DMPlexSetSNESLocalFEM(ce->idm, ctx, ctx, ctx);
400: SNESSetFromOptions(snes);
401: return(0);
402: }
404: /*@
405: PetscConvEstGetConvRate - Returns an estimate of the convergence rate for the discretization
407: Not collective
409: Input Parameter:
410: . ce - The PetscConvEst object
412: Output Parameter:
413: . alpha - The convergence rate for each field
415: Note: The convergence rate alpha is defined by
416: $ || u_\Delta - u_exact || < C \Delta^alpha
417: where u_\Delta is the discrete solution, and Delta is a measure of the discretization size. We usually use h for the
418: spatial resolution and \Delta t for the temporal resolution.
420: We solve a series of problems using increasing resolution (refined meshes or decreased timesteps), calculate an error
421: based upon the exact solution in the DS, and then fit the result to our model above using linear regression.
423: Options database keys:
424: + -snes_convergence_estimate : Execute convergence estimation inside SNESSolve() and print out the rate
425: - -ts_convergence_estimate : Execute convergence estimation inside TSSolve() and print out the rate
427: Level: intermediate
429: .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate(), SNESSolve(), TSSolve()
430: @*/
431: PetscErrorCode PetscConvEstGetConvRate(PetscConvEst ce, PetscReal alpha[])
432: {
433: PetscInt f;
437: if (ce->event < 0) {PetscLogEventRegister("ConvEst Error", PETSC_OBJECT_CLASSID, &ce->event);}
438: for (f = 0; f < ce->Nf; ++f) alpha[f] = 0.0;
439: (*ce->ops->getconvrate)(ce, alpha);
440: return(0);
441: }
443: /*@
444: PetscConvEstRateView - Displays the convergence rate to a viewer
446: Collective on SNES
448: Parameter:
449: + snes - iterative context obtained from SNESCreate()
450: . alpha - the convergence rate for each field
451: - viewer - the viewer to display the reason
453: Options Database Keys:
454: . -snes_convergence_estimate - print the convergence rate
456: Level: developer
458: .seealso: PetscConvEstGetRate()
459: @*/
460: PetscErrorCode PetscConvEstRateView(PetscConvEst ce, const PetscReal alpha[], PetscViewer viewer)
461: {
462: PetscBool isAscii;
466: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);
467: if (isAscii) {
468: PetscInt Nf = ce->Nf, f;
470: PetscViewerASCIIAddTab(viewer, ((PetscObject) ce)->tablevel);
471: PetscViewerASCIIPrintf(viewer, "L_2 convergence rate: ");
472: if (Nf > 1) {PetscViewerASCIIPrintf(viewer, "[");}
473: for (f = 0; f < Nf; ++f) {
474: if (f > 0) {PetscViewerASCIIPrintf(viewer, ", ");}
475: PetscViewerASCIIPrintf(viewer, "%#.2g", (double) alpha[f]);
476: }
477: if (Nf > 1) {PetscViewerASCIIPrintf(viewer, "]");}
478: PetscViewerASCIIPrintf(viewer, "\n");
479: PetscViewerASCIISubtractTab(viewer, ((PetscObject) ce)->tablevel);
480: }
481: return(0);
482: }
484: /*@
485: PetscConvEstCreate - Create a PetscConvEst object
487: Collective
489: Input Parameter:
490: . comm - The communicator for the PetscConvEst object
492: Output Parameter:
493: . ce - The PetscConvEst object
495: Level: beginner
497: .seealso: PetscConvEstDestroy(), PetscConvEstGetConvRate()
498: @*/
499: PetscErrorCode PetscConvEstCreate(MPI_Comm comm, PetscConvEst *ce)
500: {
505: PetscSysInitializePackage();
506: PetscHeaderCreate(*ce, PETSC_OBJECT_CLASSID, "PetscConvEst", "ConvergenceEstimator", "SNES", comm, PetscConvEstDestroy, PetscConvEstView);
507: (*ce)->monitor = PETSC_FALSE;
508: (*ce)->r = 2.0;
509: (*ce)->Nr = 4;
510: (*ce)->event = -1;
511: (*ce)->ops->setsolver = PetscConvEstSetSNES_Private;
512: (*ce)->ops->initguess = PetscConvEstInitGuessSNES_Private;
513: (*ce)->ops->computeerror = PetscConvEstComputeErrorSNES_Private;
514: (*ce)->ops->getconvrate = PetscConvEstGetConvRateSNES_Private;
515: return(0);
516: }