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