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