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