xref: /petsc/src/snes/tutorials/ex36.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
1 static char help[] = "First example in homogenization book\n\n";
2 
3 #include <petscsnes.h>
4 #include <petscdmplex.h>
5 #include <petscds.h>
6 #include <petscconvest.h>
7 #include <petscbag.h>
8 
9 /*
10   To control the refinement use -dm_plex_box_faces <n> or -dm_refine <k>, or both
11 
12   To see the exact and computed solutions
13 
14     -compare_view draw -draw_size 500,500 -draw_pause -1
15 
16   To see the delay in convergence of the discretization use
17 
18     -snes_convergence_estimate -convest_num_refine 7 -convest_monitor
19 
20   and to see the proper rate use
21 
22     -dm_refine 5 -snes_convergence_estimate -convest_num_refine 2 -convest_monitor
23 */
24 
25 typedef enum {MOD_CONSTANT, MOD_OSCILLATORY, MOD_TANH, NUM_MOD_TYPES} ModType;
26 const char *modTypes[NUM_MOD_TYPES+1] = {"constant", "oscillatory", "tanh", "unknown"};
27 
28 /* Constants */
29 enum {EPSILON, NUM_CONSTANTS};
30 
31 typedef struct {
32   PetscReal epsilon; /* Wavelength of fine scale oscillation */
33 } Parameter;
34 
35 typedef struct {
36   PetscBag bag;      /* Holds problem parameters */
37   ModType  modType;  /* Model type */
38 } AppCtx;
39 
40 static PetscErrorCode trig_homogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
41 {
42   PetscInt d;
43   *u = 1.0;
44   for (d = 0; d < dim; ++d) *u *= PetscSinReal(2.0*PETSC_PI*x[d]);
45   return 0;
46 }
47 
48 static PetscErrorCode oscillatory_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
49 {
50   Parameter      *param = (Parameter *) ctx;
51   const PetscReal eps   = param->epsilon;
52 
53   u[0] = x[0] - x[0]*x[0] + (eps / (2.*PETSC_PI))*(0.5 - x[0])*PetscSinReal(2.*PETSC_PI*x[0]/eps) + PetscSqr(eps / (2.*PETSC_PI))*(1. - PetscCosReal(2.*PETSC_PI*x[0]/eps));
54   return 0;
55 }
56 
57 static void f0_trig_homogeneous_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
58                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
59                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
60                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
61 {
62   PetscInt d;
63   for (d = 0; d < dim; ++d) {
64     PetscScalar v = 1.;
65     for (PetscInt e = 0; e < dim; e++) {
66       if (e == d) {
67         v *= -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]);
68       } else {
69         v *= PetscSinReal(2.0*PETSC_PI*x[d]);
70       }
71     }
72     f0[0] += v;
73   }
74 }
75 
76 static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
77                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
78                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
79                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
80 {
81   PetscInt d;
82   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
83 }
84 
85 static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
86                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
87                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
88                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
89 {
90   PetscInt d;
91   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
92 }
93 
94 static void f0_oscillatory_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
95                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
96                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
97                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
98 {
99   f0[0] = -1.;
100 }
101 
102 static void f1_oscillatory_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
103                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
104                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
105                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
106 {
107   const PetscReal eps = PetscRealPart(constants[EPSILON]);
108 
109   f1[0] = u_x[0] / (2. + PetscCosReal(2.*PETSC_PI*x[0]/eps));
110 }
111 
112 static void g3_oscillatory_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
113                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
114                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
115                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
116 {
117   const PetscReal eps = PetscRealPart(constants[EPSILON]);
118 
119   g3[0] = 1. / (2. + PetscCosReal(2.*PETSC_PI*x[0]/eps));
120 }
121 
122 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
123 {
124   PetscInt       mod;
125 
126   PetscFunctionBeginUser;
127   options->modType = MOD_CONSTANT;
128   PetscOptionsBegin(comm, "", "Homogenization Problem Options", "DMPLEX");
129   mod = options->modType;
130   PetscCall(PetscOptionsEList("-mod_type", "The model type", "ex36.c", modTypes, NUM_MOD_TYPES, modTypes[options->modType], &mod, NULL));
131   options->modType = (ModType) mod;
132   PetscOptionsEnd();
133   PetscFunctionReturn(0);
134 }
135 
136 static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *user)
137 {
138   PetscBag       bag;
139   Parameter     *p;
140 
141   PetscFunctionBeginUser;
142   PetscCall(PetscBagCreate(comm, sizeof(Parameter), &user->bag));
143   PetscCall(PetscBagGetData(user->bag, (void **) &p));
144   PetscCall(PetscBagSetName(user->bag, "par", "Homogenization parameters"));
145   bag  = user->bag;
146   PetscCall(PetscBagRegisterReal(bag, &p->epsilon, 1.0, "epsilon", "Wavelength of fine scale oscillation"));
147   PetscFunctionReturn(0);
148 }
149 
150 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
151 {
152   PetscFunctionBeginUser;
153   PetscCall(DMCreate(comm, dm));
154   PetscCall(DMSetType(*dm, DMPLEX));
155   PetscCall(DMSetFromOptions(*dm));
156   PetscCall(DMSetApplicationContext(*dm, user));
157   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
158   PetscFunctionReturn(0);
159 }
160 
161 static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
162 {
163   PetscDS              ds;
164   DMLabel              label;
165   PetscSimplePointFunc ex;
166   const PetscInt       id = 1;
167   void                *ctx;
168 
169   PetscFunctionBeginUser;
170   PetscCall(DMGetDS(dm, &ds));
171   PetscCall(PetscBagGetData(user->bag, (void **) &ctx));
172   switch (user->modType) {
173     case MOD_CONSTANT:
174       PetscCall(PetscDSSetResidual(ds, 0, f0_trig_homogeneous_u, f1_u));
175       PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu));
176       PetscCall(DMGetLabel(dm, "marker", &label));
177       ex   = trig_homogeneous_u;
178       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) ex, NULL, ctx, NULL));
179       break;
180     case MOD_OSCILLATORY:
181       PetscCall(PetscDSSetResidual(ds, 0, f0_oscillatory_u, f1_oscillatory_u));
182       PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_oscillatory_uu));
183       PetscCall(DMGetLabel(dm, "marker", &label));
184       ex   = oscillatory_u;
185       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) ex, NULL, ctx, NULL));
186       break;
187     default: SETERRQ(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Unsupported model type: %s (%d)", modTypes[PetscMin(user->modType, NUM_MOD_TYPES)], user->modType);
188   }
189   PetscCall(PetscDSSetExactSolution(ds, 0, ex, ctx));
190   /* Setup constants */
191   {
192     Parameter  *param;
193     PetscScalar constants[NUM_CONSTANTS];
194 
195     PetscCall(PetscBagGetData(user->bag, (void **) &param));
196 
197     constants[EPSILON] = param->epsilon;
198     PetscCall(PetscDSSetConstants(ds, NUM_CONSTANTS, constants));
199   }
200   PetscFunctionReturn(0);
201 }
202 
203 static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
204 {
205   DM             cdm = dm;
206   PetscFE        fe;
207   PetscBool      simplex;
208   PetscInt       dim;
209   char           prefix[PETSC_MAX_PATH_LEN];
210 
211   PetscFunctionBeginUser;
212   PetscCall(DMGetDimension(dm, &dim));
213   PetscCall(DMPlexIsSimplex(dm, &simplex));
214   /* Create finite element */
215   PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
216   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe));
217   PetscCall(PetscObjectSetName((PetscObject) fe, name));
218   /* Set discretization and boundary conditions for each mesh */
219   PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe));
220   PetscCall(DMCreateDS(dm));
221   PetscCall((*setup)(dm, user));
222   while (cdm) {
223     PetscCall(DMCopyDisc(dm,cdm));
224     PetscCall(DMGetCoarseDM(cdm, &cdm));
225   }
226   PetscCall(PetscFEDestroy(&fe));
227   PetscFunctionReturn(0);
228 }
229 
230 static PetscErrorCode CompareView(Vec u)
231 {
232   DM                dm;
233   Vec               v[2], lv[2], exact;
234   PetscOptions      options;
235   PetscViewer       viewer;
236   PetscViewerFormat format;
237   PetscBool         flg;
238   PetscInt          i;
239   const char       *name, *prefix;
240 
241   PetscFunctionBeginUser;
242   PetscCall(VecGetDM(u, &dm));
243   PetscCall(PetscObjectGetOptions((PetscObject) dm, &options));
244   PetscCall(PetscObjectGetOptionsPrefix((PetscObject) dm, &prefix));
245   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject) dm), options, prefix, "-compare_view", &viewer, &format, &flg));
246   if (flg) {
247     PetscCall(DMGetGlobalVector(dm, &exact));
248     PetscCall(DMComputeExactSolution(dm, 0.0, exact, NULL));
249     v[0] = u;
250     v[1] = exact;
251     for (i = 0; i < 2; ++i) {
252       PetscCall(DMGetLocalVector(dm, &lv[i]));
253       PetscCall(PetscObjectGetName((PetscObject) v[i], &name));
254       PetscCall(PetscObjectSetName((PetscObject) lv[i], name));
255       PetscCall(DMGlobalToLocalBegin(dm, v[i], INSERT_VALUES, lv[i]));
256       PetscCall(DMGlobalToLocalEnd(dm, v[i], INSERT_VALUES, lv[i]));
257       PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, lv[i], 0., NULL, NULL, NULL));
258     }
259     PetscCall(DMPlexVecView1D(dm, 2, lv, viewer));
260     for (i = 0; i < 2; ++i) PetscCall(DMRestoreLocalVector(dm, &lv[i]));
261     PetscCall(DMRestoreGlobalVector(dm, &exact));
262     PetscCall(PetscViewerDestroy(&viewer));
263   }
264   PetscFunctionReturn(0);
265 }
266 
267 typedef struct
268 {
269   Mat Mcoarse;   /* Mass matrix on the coarse space */
270   Mat Mfine;     /* Mass matrix on the fine space */
271   Mat Ifine;     /* Interpolator from coarse to fine */
272   Vec Iscale;    /* Scaling vector for restriction */
273   KSP kspCoarse; /* Solver for the coarse mass matrix */
274   Vec tmpfine;   /* Temporary vector in the fine space */
275   Vec tmpcoarse; /* Temporary vector in the coarse space */
276 } ProjStruct;
277 
278 static PetscErrorCode DestroyCoarseProjection(Mat Pi)
279 {
280   ProjStruct    *ctx;
281 
282   PetscFunctionBegin;
283   PetscCall(MatShellGetContext(Pi, (void **) &ctx));
284   PetscCall(MatDestroy(&ctx->Mcoarse));
285   PetscCall(MatDestroy(&ctx->Mfine));
286   PetscCall(MatDestroy(&ctx->Ifine));
287   PetscCall(VecDestroy(&ctx->Iscale));
288   PetscCall(KSPDestroy(&ctx->kspCoarse));
289   PetscCall(VecDestroy(&ctx->tmpcoarse));
290   PetscCall(VecDestroy(&ctx->tmpfine));
291   PetscCall(PetscFree(ctx));
292   PetscCall(MatShellSetContext(Pi, NULL));
293   PetscFunctionReturn(0);
294 }
295 
296 static PetscErrorCode CoarseProjection(Mat Pi, Vec x, Vec y)
297 {
298   ProjStruct    *ctx;
299 
300   PetscFunctionBegin;
301   PetscCall(MatShellGetContext(Pi, (void **) &ctx));
302   PetscCall(MatMult(ctx->Mfine, x, ctx->tmpfine));
303   PetscCall(PetscObjectSetName((PetscObject) ctx->tmpfine, "Fine DG RHS"));
304   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) ctx->tmpfine, "fine_dg_"));
305   PetscCall(VecViewFromOptions(ctx->tmpfine, NULL, "-rhs_view"));
306   PetscCall(MatMultTranspose(ctx->Ifine, ctx->tmpfine, ctx->tmpcoarse));
307   PetscCall(PetscObjectSetName((PetscObject) ctx->tmpcoarse, "Coarse DG RHS"));
308   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) ctx->tmpcoarse, "coarse_dg_"));
309   PetscCall(VecViewFromOptions(ctx->tmpcoarse, NULL, "-rhs_view"));
310   PetscCall(VecPointwiseMult(ctx->tmpcoarse, ctx->Iscale, ctx->tmpcoarse));
311   PetscCall(VecViewFromOptions(ctx->tmpcoarse, NULL, "-rhs_view"));
312   PetscCall(KSPSolve(ctx->kspCoarse, ctx->tmpcoarse, y));
313   PetscFunctionReturn(0);
314 }
315 
316 static PetscErrorCode CreateCoarseProjection(DM dmc, DM dmf, Mat *Pi)
317 {
318   ProjStruct    *ctx;
319   PetscInt       m, n, M, N;
320 
321   PetscFunctionBegin;
322   PetscCall(PetscMalloc1(1, &ctx));
323   PetscCall(DMCreateGlobalVector(dmc, &ctx->tmpcoarse));
324   PetscCall(DMCreateGlobalVector(dmf, &ctx->tmpfine));
325   PetscCall(VecGetLocalSize(ctx->tmpcoarse, &m));
326   PetscCall(VecGetSize(ctx->tmpcoarse, &M));
327   PetscCall(VecGetLocalSize(ctx->tmpfine, &n));
328   PetscCall(VecGetSize(ctx->tmpfine, &N));
329   PetscCall(DMCreateMassMatrix(dmc, dmc, &ctx->Mcoarse));
330   PetscCall(DMCreateMassMatrix(dmf, dmf, &ctx->Mfine));
331   PetscCall(DMCreateInterpolation(dmc, dmf, &ctx->Ifine, &ctx->Iscale));
332   PetscCall(KSPCreate(PetscObjectComm((PetscObject) dmc), &ctx->kspCoarse));
333   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) ctx->kspCoarse, "coarse_"));
334   PetscCall(KSPSetOperators(ctx->kspCoarse, ctx->Mcoarse, ctx->Mcoarse));
335   PetscCall(KSPSetFromOptions(ctx->kspCoarse));
336   PetscCall(MatCreateShell(PetscObjectComm((PetscObject) dmc), m, n, M, N, ctx, Pi));
337   PetscCall(MatShellSetOperation(*Pi, MATOP_DESTROY, (void (*)(void)) DestroyCoarseProjection));
338   PetscCall(MatShellSetOperation(*Pi, MATOP_MULT, (void (*)(void)) CoarseProjection));
339   PetscFunctionReturn(0);
340 }
341 
342 typedef struct
343 {
344   Mat Ifdg; /* Embed the fine space into the DG version */
345   Mat Pi;   /* The L_2 stable projection to the DG coarse space */
346   Vec tmpc; /* A temporary vector in the DG coarse space */
347   Vec tmpf; /* A temporary vector in the DG fine space */
348 } QuasiInterp;
349 
350 static PetscErrorCode DestroyQuasiInterpolator(Mat P)
351 {
352   QuasiInterp   *ctx;
353 
354   PetscFunctionBegin;
355   PetscCall(MatShellGetContext(P, (void **) &ctx));
356   PetscCall(MatDestroy(&ctx->Ifdg));
357   PetscCall(MatDestroy(&ctx->Pi));
358   PetscCall(VecDestroy(&ctx->tmpc));
359   PetscCall(VecDestroy(&ctx->tmpf));
360   PetscCall(PetscFree(ctx));
361   PetscCall(MatShellSetContext(P, NULL));
362   PetscFunctionReturn(0);
363 }
364 
365 static PetscErrorCode QuasiInterpolate(Mat P, Vec x, Vec y)
366 {
367   QuasiInterp   *ctx;
368   DM             dmcdg, dmc;
369   Vec            ly;
370 
371   PetscFunctionBegin;
372   PetscCall(MatShellGetContext(P, (void **) &ctx));
373   PetscCall(MatMult(ctx->Ifdg, x, ctx->tmpf));
374 
375   PetscCall(PetscObjectSetName((PetscObject) ctx->tmpf, "Fine DG Potential"));
376   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) ctx->tmpf, "fine_dg_"));
377   PetscCall(VecViewFromOptions(ctx->tmpf, NULL, "-vec_view"));
378   PetscCall(MatMult(ctx->Pi, x, ctx->tmpc));
379 
380   PetscCall(PetscObjectSetName((PetscObject) ctx->tmpc, "Coarse DG Potential"));
381   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) ctx->tmpc, "coarse_dg_"));
382   PetscCall(VecViewFromOptions(ctx->tmpc, NULL, "-vec_view"));
383   PetscCall(VecGetDM(ctx->tmpc, &dmcdg));
384 
385   PetscCall(VecGetDM(y, &dmc));
386   PetscCall(DMGetLocalVector(dmc, &ly));
387   PetscCall(DMPlexComputeClementInterpolant(dmcdg, ctx->tmpc, ly));
388   PetscCall(DMLocalToGlobalBegin(dmc, ly, INSERT_VALUES, y));
389   PetscCall(DMLocalToGlobalEnd(dmc, ly, INSERT_VALUES, y));
390   PetscCall(DMRestoreLocalVector(dmc, &ly));
391   PetscFunctionReturn(0);
392 }
393 
394 static PetscErrorCode CreateQuasiInterpolator(DM dmc, DM dmf, Mat *P)
395 {
396   QuasiInterp   *ctx;
397   DM             dmcdg, dmfdg;
398   PetscFE        fe;
399   Vec            x, y;
400   DMPolytopeType ct;
401   PetscInt       dim, cStart, m, n, M, N;
402 
403   PetscFunctionBegin;
404   PetscCall(PetscCalloc1(1, &ctx));
405   PetscCall(DMGetGlobalVector(dmc, &x));
406   PetscCall(DMGetGlobalVector(dmf, &y));
407   PetscCall(VecGetLocalSize(x, &m));
408   PetscCall(VecGetSize(x, &M));
409   PetscCall(VecGetLocalSize(y, &n));
410   PetscCall(VecGetSize(y, &N));
411   PetscCall(DMRestoreGlobalVector(dmc, &x));
412   PetscCall(DMRestoreGlobalVector(dmf, &y));
413 
414   PetscCall(DMClone(dmf, &dmfdg));
415   PetscCall(DMGetDimension(dmfdg, &dim));
416   PetscCall(DMPlexGetHeightStratum(dmfdg, 0, &cStart, NULL));
417   PetscCall(DMPlexGetCellType(dmfdg, cStart, &ct));
418   PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, "fine_dg_", PETSC_DETERMINE, &fe));
419   PetscCall(DMSetField(dmfdg, 0, NULL, (PetscObject) fe));
420   PetscCall(PetscFEDestroy(&fe));
421   PetscCall(DMCreateDS(dmfdg));
422   PetscCall(DMCreateInterpolation(dmf, dmfdg, &ctx->Ifdg, NULL));
423   PetscCall(DMCreateGlobalVector(dmfdg, &ctx->tmpf));
424   PetscCall(DMDestroy(&dmfdg));
425 
426   PetscCall(DMClone(dmc, &dmcdg));
427   PetscCall(DMGetDimension(dmcdg, &dim));
428   PetscCall(DMPlexGetHeightStratum(dmcdg, 0, &cStart, NULL));
429   PetscCall(DMPlexGetCellType(dmcdg, cStart, &ct));
430   PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, "coarse_dg_", PETSC_DETERMINE, &fe));
431   PetscCall(DMSetField(dmcdg, 0, NULL, (PetscObject) fe));
432   PetscCall(PetscFEDestroy(&fe));
433   PetscCall(DMCreateDS(dmcdg));
434 
435   PetscCall(CreateCoarseProjection(dmcdg, dmf, &ctx->Pi));
436   PetscCall(DMCreateGlobalVector(dmcdg, &ctx->tmpc));
437   PetscCall(DMDestroy(&dmcdg));
438 
439   PetscCall(MatCreateShell(PetscObjectComm((PetscObject) dmc), m, n, M, N, ctx, P));
440   PetscCall(MatShellSetOperation(*P, MATOP_DESTROY, (void (*)(void)) DestroyQuasiInterpolator));
441   PetscCall(MatShellSetOperation(*P, MATOP_MULT, (void (*)(void)) QuasiInterpolate));
442   PetscFunctionReturn(0);
443 }
444 
445 static PetscErrorCode CoarseTest(DM dm, Vec u, AppCtx *user)
446 {
447   DM             dmc;
448   Mat            P;    /* The quasi-interpolator to the coarse space */
449   Vec            uc;
450 
451   PetscFunctionBegin;
452   if (user->modType == MOD_CONSTANT) PetscFunctionReturn(0);
453   PetscCall(DMCreate(PetscObjectComm((PetscObject) dm), &dmc));
454   PetscCall(DMSetType(dmc, DMPLEX));
455   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) dmc, "coarse_"));
456   PetscCall(DMSetApplicationContext(dmc, user));
457   PetscCall(DMSetFromOptions(dmc));
458   PetscCall(DMViewFromOptions(dmc, NULL, "-dm_view"));
459 
460   PetscCall(SetupDiscretization(dmc, "potential", SetupPrimalProblem, user));
461   PetscCall(DMCreateGlobalVector(dmc, &uc));
462   PetscCall(PetscObjectSetName((PetscObject) uc, "potential"));
463   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) uc, "coarse_"));
464 
465   PetscCall(CreateQuasiInterpolator(dmc, dm, &P));
466 #if 1
467   PetscCall(MatMult(P, u, uc));
468 #else
469   {
470     Mat In;
471     Vec sc;
472 
473     PetscCall(DMCreateInterpolation(dmc, dm, &In, &sc));
474     PetscCall(MatMultTranspose(In, u, uc));
475     PetscCall(VecPointwiseMult(uc, sc, uc));
476     PetscCall(MatDestroy(&In));
477     PetscCall(VecDestroy(&sc));
478   }
479 #endif
480   PetscCall(CompareView(uc));
481 
482   PetscCall(MatDestroy(&P));
483   PetscCall(VecDestroy(&uc));
484   PetscCall(DMDestroy(&dmc));
485   PetscFunctionReturn(0);
486 }
487 
488 int main(int argc, char **argv)
489 {
490   DM             dm;   /* Problem specification */
491   SNES           snes; /* Nonlinear solver */
492   Vec            u;    /* Solutions */
493   AppCtx         user; /* User-defined work context */
494 
495   PetscFunctionBeginUser;
496   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
497   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
498   PetscCall(SetupParameters(PETSC_COMM_WORLD, &user));
499   /* Primal system */
500   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
501   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
502   PetscCall(SNESSetDM(snes, dm));
503   PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user));
504   PetscCall(DMCreateGlobalVector(dm, &u));
505   PetscCall(VecSet(u, 0.0));
506   PetscCall(PetscObjectSetName((PetscObject) u, "potential"));
507   PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user));
508   PetscCall(SNESSetFromOptions(snes));
509   PetscCall(DMSNESCheckFromOptions(snes, u));
510   PetscCall(SNESSolve(snes, NULL, u));
511   PetscCall(SNESGetSolution(snes, &u));
512   PetscCall(VecViewFromOptions(u, NULL, "-potential_view"));
513   PetscCall(CompareView(u));
514   /* Looking at a coarse problem */
515   PetscCall(CoarseTest(dm, u, &user));
516   /* Cleanup */
517   PetscCall(VecDestroy(&u));
518   PetscCall(SNESDestroy(&snes));
519   PetscCall(DMDestroy(&dm));
520   PetscCall(PetscBagDestroy(&user.bag));
521   PetscCall(PetscFinalize());
522   return 0;
523 }
524 
525 /*TEST
526 
527   test:
528     suffix: 1d_p1_constant
529     args: -dm_plex_dim 1 -dm_plex_box_faces 4 -potential_petscspace_degree 1 -dmsnes_check
530 
531   test:
532     suffix: 1d_p1_constant_conv
533     args: -dm_plex_dim 1 -dm_plex_box_faces 4 -potential_petscspace_degree 1 \
534           -snes_convergence_estimate -convest_num_refine 2
535 
536   test:
537     suffix: 1d_p1_oscillatory
538     args: -mod_type oscillatory -epsilon 0.03125 \
539           -dm_plex_dim 1 -dm_plex_box_faces 4 -potential_petscspace_degree 1 -dm_refine 2 -dmsnes_check \
540           -coarse_dm_plex_dim 1 -coarse_dm_plex_box_faces 4 -coarse_dm_plex_hash_location \
541           -fine_dg_petscspace_degree 1 -fine_dg_petscdualspace_lagrange_continuity 0 \
542           -coarse_dg_petscspace_degree 1 -coarse_dg_petscdualspace_lagrange_continuity 0
543 
544 TEST*/
545