xref: /petsc/src/snes/tutorials/ex13.c (revision 1fb7b255396822dd64a4ea9c656967637c985cf5)
1 static char help[] = "Poisson Problem in 2d and 3d with finite elements.\n\
2 We solve the Poisson problem in a rectangular\n\
3 domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\
4 This example supports automatic convergence estimation\n\
5 and eventually adaptivity.\n\n\n";
6 
7 #include <petscdmplex.h>
8 #include <petscsnes.h>
9 #include <petscds.h>
10 #include <petscconvest.h>
11 
12 typedef struct {
13   /* Domain and mesh definition */
14   PetscBool spectral; /* Look at the spectrum along planes in the solution */
15   PetscBool shear;    /* Shear the domain */
16   PetscBool adjoint;  /* Solve the adjoint problem */
17   PetscBool homogeneous;
18 } AppCtx;
19 
20 static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
21 {
22   *u = 0.0;
23   return 0;
24 }
25 
26 static PetscErrorCode trig_inhomogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
27 {
28   PetscInt d;
29   *u = 0.0;
30   for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0*PETSC_PI*x[d]);
31   return 0;
32 }
33 
34 static PetscErrorCode trig_homogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
35 {
36   PetscInt d;
37   *u = 1.0;
38   for (d = 0; d < dim; ++d) *u *= PetscSinReal(2.0*PETSC_PI*x[d]);
39   return 0;
40 }
41 
42 /* Compute integral of (residual of solution)*(adjoint solution - projection of adjoint solution) */
43 static void obj_error_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
44                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
45                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
46                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar obj[])
47 {
48   obj[0] = a[aOff[0]]*(u[0] - a[aOff[1]]);
49 }
50 
51 static void f0_trig_inhomogeneous_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
52                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
53                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
54                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
55 {
56   PetscInt d;
57   for (d = 0; d < dim; ++d) f0[0] += -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]);
58 }
59 
60 static void f0_trig_homogeneous_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
61                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
62                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
63                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
64 {
65   PetscInt d;
66   for (d = 0; d < dim; ++d) {
67     PetscScalar v = 1.;
68     for (PetscInt e = 0; e < dim; e++) {
69       if (e == d) {
70         v *= -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]);
71       } else {
72         v *= PetscSinReal(2.0*PETSC_PI*x[d]);
73       }
74     }
75     f0[0] += v;
76   }
77 }
78 
79 static void f0_unity_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
80                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
81                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
82                        PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
83 {
84   f0[0] = 1.0;
85 }
86 
87 static void f0_identityaux_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
88                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
89                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
90                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
91 {
92   f0[0] = a[0];
93 }
94 
95 static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
96                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
97                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
98                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
99 {
100   PetscInt d;
101   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
102 }
103 
104 static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
105                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
106                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
107                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
108 {
109   PetscInt d;
110   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
111 }
112 
113 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
114 {
115   PetscErrorCode ierr;
116 
117   PetscFunctionBeginUser;
118   options->shear    = PETSC_FALSE;
119   options->spectral = PETSC_FALSE;
120   options->adjoint  = PETSC_FALSE;
121   options->homogeneous = PETSC_FALSE;
122 
123   ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");CHKERRQ(ierr);
124   ierr = PetscOptionsBool("-shear", "Shear the domain", "ex13.c", options->shear, &options->shear, NULL);CHKERRQ(ierr);
125   ierr = PetscOptionsBool("-spectral", "Look at the spectrum along planes of the solution", "ex13.c", options->spectral, &options->spectral, NULL);CHKERRQ(ierr);
126   ierr = PetscOptionsBool("-adjoint", "Solve the adjoint problem", "ex13.c", options->adjoint, &options->adjoint, NULL);CHKERRQ(ierr);
127   ierr = PetscOptionsBool("-homogeneous", "Use homogeneous boundary conditions", "ex13.c", options->homogeneous, &options->homogeneous, NULL);CHKERRQ(ierr);
128   ierr = PetscOptionsEnd();
129   PetscFunctionReturn(0);
130 }
131 
132 static PetscErrorCode CreateSpectralPlanes(DM dm, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user)
133 {
134   PetscSection       coordSection;
135   Vec                coordinates;
136   const PetscScalar *coords;
137   PetscInt           dim, p, vStart, vEnd, v;
138   PetscErrorCode     ierr;
139 
140   PetscFunctionBeginUser;
141   ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr);
142   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
143   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
144   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
145   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
146   for (p = 0; p < numPlanes; ++p) {
147     DMLabel label;
148     char    name[PETSC_MAX_PATH_LEN];
149 
150     ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%D", p);CHKERRQ(ierr);
151     ierr = DMCreateLabel(dm, name);CHKERRQ(ierr);
152     ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr);
153     ierr = DMLabelAddStratum(label, 1);CHKERRQ(ierr);
154     for (v = vStart; v < vEnd; ++v) {
155       PetscInt off;
156 
157       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
158       if (PetscAbsReal(planeCoord[p] - PetscRealPart(coords[off+planeDir[p]])) < PETSC_SMALL) {
159         ierr = DMLabelSetValue(label, v, 1);CHKERRQ(ierr);
160       }
161     }
162   }
163   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
164   PetscFunctionReturn(0);
165 }
166 
167 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
168 {
169   PetscErrorCode ierr;
170 
171   PetscFunctionBeginUser;
172   /* Create box mesh */
173   ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
174   /* TODO: This should be pulled into the library */
175   {
176     char      convType[256];
177     PetscBool flg;
178 
179     ierr = PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");CHKERRQ(ierr);
180     ierr = PetscOptionsFList("-dm_plex_convert_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);CHKERRQ(ierr);
181     ierr = PetscOptionsEnd();
182     if (flg) {
183       DM dmConv;
184 
185       ierr = DMConvert(*dm,convType,&dmConv);CHKERRQ(ierr);
186       if (dmConv) {
187         ierr = DMDestroy(dm);CHKERRQ(ierr);
188         *dm  = dmConv;
189       }
190     }
191   }
192   if (user->shear) {ierr = DMPlexShearGeometry(*dm, DM_X, NULL);CHKERRQ(ierr);}
193   ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr);
194 
195   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
196   ierr = DMSetApplicationContext(*dm, user);CHKERRQ(ierr);
197   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
198   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
199   if (user->spectral) {
200     PetscInt  planeDir[2]   = {0,  1};
201     PetscReal planeCoord[2] = {0., 1.};
202 
203     ierr = CreateSpectralPlanes(*dm, 2, planeDir, planeCoord, user);CHKERRQ(ierr);
204   }
205   PetscFunctionReturn(0);
206 }
207 
208 static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
209 {
210   PetscDS        prob;
211   const PetscInt id = 1;
212   PetscPointFunc f0 = user->homogeneous ? f0_trig_homogeneous_u : f0_trig_inhomogeneous_u;
213   PetscErrorCode (*ex)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *) = user->homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u;
214   PetscErrorCode ierr;
215 
216   PetscFunctionBeginUser;
217   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
218   ierr = PetscDSSetResidual(prob, 0, f0, f1_u);CHKERRQ(ierr);
219   ierr = PetscDSSetExactSolution(prob, 0, ex, user);CHKERRQ(ierr);
220   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) ex, NULL, 1, &id, user);CHKERRQ(ierr);
221   ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr);
222   PetscFunctionReturn(0);
223 }
224 
225 static PetscErrorCode SetupAdjointProblem(DM dm, AppCtx *user)
226 {
227   PetscDS        prob;
228   const PetscInt id = 1;
229   PetscErrorCode ierr;
230 
231   PetscFunctionBeginUser;
232   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
233   ierr = PetscDSSetResidual(prob, 0, f0_unity_u, f1_u);CHKERRQ(ierr);
234   ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr);
235   ierr = PetscDSSetObjective(prob, 0, obj_error_u);CHKERRQ(ierr);
236   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) zero, NULL, 1, &id, user);CHKERRQ(ierr);
237   PetscFunctionReturn(0);
238 }
239 
240 static PetscErrorCode SetupErrorProblem(DM dm, AppCtx *user)
241 {
242   PetscDS        prob;
243   PetscErrorCode ierr;
244 
245   PetscFunctionBeginUser;
246   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
247   PetscFunctionReturn(0);
248 }
249 
250 static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
251 {
252   DM             cdm = dm;
253   PetscFE        fe;
254   DMPolytopeType ct;
255   PetscBool      simplex;
256   PetscInt       dim, cStart;
257   char           prefix[PETSC_MAX_PATH_LEN];
258   PetscErrorCode ierr;
259 
260   PetscFunctionBeginUser;
261   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
262   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr);
263   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
264   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
265   /* Create finite element */
266   ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name);CHKERRQ(ierr);
267   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe);CHKERRQ(ierr);
268   ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr);
269   /* Set discretization and boundary conditions for each mesh */
270   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr);
271   ierr = DMCreateDS(dm);CHKERRQ(ierr);
272   ierr = (*setup)(dm, user);CHKERRQ(ierr);
273   while (cdm) {
274     ierr = DMCopyDisc(dm,cdm);CHKERRQ(ierr);
275     /* TODO: Check whether the boundary of coarse meshes is marked */
276     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
277   }
278   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
279   PetscFunctionReturn(0);
280 }
281 
282 static PetscErrorCode ComputeSpectral(DM dm, Vec u, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user)
283 {
284   MPI_Comm           comm;
285   PetscSection       coordSection, section;
286   Vec                coordinates, uloc;
287   const PetscScalar *coords, *array;
288   PetscInt           p;
289   PetscMPIInt        size, rank;
290   PetscErrorCode     ierr;
291 
292   PetscFunctionBeginUser;
293   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
294   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
295   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
296   ierr = DMGetLocalVector(dm, &uloc);CHKERRQ(ierr);
297   ierr = DMGlobalToLocalBegin(dm, u, INSERT_VALUES, uloc);CHKERRQ(ierr);
298   ierr = DMGlobalToLocalEnd(dm, u, INSERT_VALUES, uloc);CHKERRQ(ierr);
299   ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, uloc, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
300   ierr = VecViewFromOptions(uloc, NULL, "-sol_view");CHKERRQ(ierr);
301   ierr = DMGetLocalSection(dm, &section);CHKERRQ(ierr);
302   ierr = VecGetArrayRead(uloc, &array);CHKERRQ(ierr);
303   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
304   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
305   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
306   for (p = 0; p < numPlanes; ++p) {
307     DMLabel         label;
308     char            name[PETSC_MAX_PATH_LEN];
309     Mat             F;
310     Vec             x, y;
311     IS              stratum;
312     PetscReal      *ray, *gray;
313     PetscScalar    *rvals, *svals, *gsvals;
314     PetscInt       *perm, *nperm;
315     PetscInt        n, N, i, j, off, offu;
316     const PetscInt *points;
317 
318     ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%D", p);CHKERRQ(ierr);
319     ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr);
320     ierr = DMLabelGetStratumIS(label, 1, &stratum);CHKERRQ(ierr);
321     ierr = ISGetLocalSize(stratum, &n);CHKERRQ(ierr);
322     ierr = ISGetIndices(stratum, &points);CHKERRQ(ierr);
323     ierr = PetscMalloc2(n, &ray, n, &svals);CHKERRQ(ierr);
324     for (i = 0; i < n; ++i) {
325       ierr = PetscSectionGetOffset(coordSection, points[i], &off);CHKERRQ(ierr);
326       ierr = PetscSectionGetOffset(section, points[i], &offu);CHKERRQ(ierr);
327       ray[i]   = PetscRealPart(coords[off+((planeDir[p]+1)%2)]);
328       svals[i] = array[offu];
329     }
330     /* Gather the ray data to proc 0 */
331     if (size > 1) {
332       PetscMPIInt *cnt, *displs, p;
333 
334       ierr = PetscCalloc2(size, &cnt, size, &displs);CHKERRQ(ierr);
335       ierr = MPI_Gather(&n, 1, MPIU_INT, cnt, 1, MPIU_INT, 0, comm);CHKERRMPI(ierr);
336       for (p = 1; p < size; ++p) displs[p] = displs[p-1] + cnt[p-1];
337       N = displs[size-1] + cnt[size-1];
338       ierr = PetscMalloc2(N, &gray, N, &gsvals);CHKERRQ(ierr);
339       ierr = MPI_Gatherv(ray, n, MPIU_REAL, gray, cnt, displs, MPIU_REAL, 0, comm);CHKERRMPI(ierr);
340       ierr = MPI_Gatherv(svals, n, MPIU_SCALAR, gsvals, cnt, displs, MPIU_SCALAR, 0, comm);CHKERRMPI(ierr);
341       ierr = PetscFree2(cnt, displs);CHKERRQ(ierr);
342     } else {
343       N      = n;
344       gray   = ray;
345       gsvals = svals;
346     }
347     if (!rank) {
348       /* Sort point along ray */
349       ierr = PetscMalloc2(N, &perm, N, &nperm);CHKERRQ(ierr);
350       for (i = 0; i < N; ++i) {perm[i] = i;}
351       ierr = PetscSortRealWithPermutation(N, gray, perm);CHKERRQ(ierr);
352       /* Count duplicates and squish mapping */
353       nperm[0] = perm[0];
354       for (i = 1, j = 1; i < N; ++i) {
355         if (PetscAbsReal(gray[perm[i]] - gray[perm[i-1]]) > PETSC_SMALL) nperm[j++] = perm[i];
356       }
357       /* Create FFT structs */
358       ierr = MatCreateFFT(PETSC_COMM_SELF, 1, &j, MATFFTW, &F);CHKERRQ(ierr);
359       ierr = MatCreateVecs(F, &x, &y);CHKERRQ(ierr);
360       ierr = PetscObjectSetName((PetscObject) y, name);CHKERRQ(ierr);
361       ierr = VecGetArray(x, &rvals);CHKERRQ(ierr);
362       for (i = 0, j = 0; i < N; ++i) {
363         if (i > 0 && PetscAbsReal(gray[perm[i]] - gray[perm[i-1]]) < PETSC_SMALL) continue;
364         rvals[j] = gsvals[nperm[j]];
365         ++j;
366       }
367       ierr = PetscFree2(perm, nperm);CHKERRQ(ierr);
368       if (size > 1) {ierr = PetscFree2(gray, gsvals);CHKERRQ(ierr);}
369       ierr = VecRestoreArray(x, &rvals);CHKERRQ(ierr);
370       /* Do FFT along the ray */
371       ierr = MatMult(F, x, y);CHKERRQ(ierr);
372       /* Chop FFT */
373       ierr = VecChop(y, PETSC_SMALL);CHKERRQ(ierr);
374       ierr = VecViewFromOptions(x, NULL, "-real_view");CHKERRQ(ierr);
375       ierr = VecViewFromOptions(y, NULL, "-fft_view");CHKERRQ(ierr);
376       ierr = VecDestroy(&x);CHKERRQ(ierr);
377       ierr = VecDestroy(&y);CHKERRQ(ierr);
378       ierr = MatDestroy(&F);CHKERRQ(ierr);
379     }
380     ierr = ISRestoreIndices(stratum, &points);CHKERRQ(ierr);
381     ierr = ISDestroy(&stratum);CHKERRQ(ierr);
382     ierr = PetscFree2(ray, svals);CHKERRQ(ierr);
383   }
384   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
385   ierr = VecRestoreArrayRead(uloc, &array);CHKERRQ(ierr);
386   ierr = DMRestoreLocalVector(dm, &uloc);CHKERRQ(ierr);
387   PetscFunctionReturn(0);
388 }
389 
390 int main(int argc, char **argv)
391 {
392   DM             dm;   /* Problem specification */
393   SNES           snes; /* Nonlinear solver */
394   Vec            u;    /* Solutions */
395   AppCtx         user; /* User-defined work context */
396   PetscErrorCode ierr;
397 
398   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
399   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
400   /* Primal system */
401   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
402   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
403   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
404   ierr = SetupDiscretization(dm, "potential", SetupPrimalProblem, &user);CHKERRQ(ierr);
405   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
406   ierr = VecSet(u, 0.0);CHKERRQ(ierr);
407   ierr = PetscObjectSetName((PetscObject) u, "potential");CHKERRQ(ierr);
408   ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr);
409   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
410   ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
411   ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
412   ierr = VecViewFromOptions(u, NULL, "-potential_view");CHKERRQ(ierr);
413   if (user.spectral) {
414     PetscInt  planeDir[2]   = {0,  1};
415     PetscReal planeCoord[2] = {0., 1.};
416 
417     ierr = ComputeSpectral(dm, u, 2, planeDir, planeCoord, &user);CHKERRQ(ierr);
418   }
419   /* Adjoint system */
420   if (user.adjoint) {
421     DM   dmAdj;
422     SNES snesAdj;
423     Vec  uAdj;
424 
425     ierr = SNESCreate(PETSC_COMM_WORLD, &snesAdj);CHKERRQ(ierr);
426     ierr = PetscObjectSetOptionsPrefix((PetscObject) snesAdj, "adjoint_");CHKERRQ(ierr);
427     ierr = DMClone(dm, &dmAdj);CHKERRQ(ierr);
428     ierr = SNESSetDM(snesAdj, dmAdj);CHKERRQ(ierr);
429     ierr = SetupDiscretization(dmAdj, "adjoint", SetupAdjointProblem, &user);CHKERRQ(ierr);
430     ierr = DMCreateGlobalVector(dmAdj, &uAdj);CHKERRQ(ierr);
431     ierr = VecSet(uAdj, 0.0);CHKERRQ(ierr);
432     ierr = PetscObjectSetName((PetscObject) uAdj, "adjoint");CHKERRQ(ierr);
433     ierr = DMPlexSetSNESLocalFEM(dmAdj, &user, &user, &user);CHKERRQ(ierr);
434     ierr = SNESSetFromOptions(snesAdj);CHKERRQ(ierr);
435     ierr = SNESSolve(snesAdj, NULL, uAdj);CHKERRQ(ierr);
436     ierr = SNESGetSolution(snesAdj, &uAdj);CHKERRQ(ierr);
437     ierr = VecViewFromOptions(uAdj, NULL, "-adjoint_view");CHKERRQ(ierr);
438     /* Error representation */
439     {
440       DM        dmErr, dmErrAux, dms[2];
441       Vec       errorEst, errorL2, uErr, uErrLoc, uAdjLoc, uAdjProj;
442       IS       *subis;
443       PetscReal errorEstTot, errorL2Norm, errorL2Tot;
444       PetscInt  N, i;
445       PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *) = {user.homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u};
446       void (*identity[1])(PetscInt, PetscInt, PetscInt,
447                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
448                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
449                           PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]) = {f0_identityaux_u};
450       void            *ctxs[1] = {0};
451 
452       ctxs[0] = &user;
453       ierr = DMClone(dm, &dmErr);CHKERRQ(ierr);
454       ierr = SetupDiscretization(dmErr, "error", SetupErrorProblem, &user);CHKERRQ(ierr);
455       ierr = DMGetGlobalVector(dmErr, &errorEst);CHKERRQ(ierr);
456       ierr = DMGetGlobalVector(dmErr, &errorL2);CHKERRQ(ierr);
457       /*   Compute auxiliary data (solution and projection of adjoint solution) */
458       ierr = DMGetLocalVector(dmAdj, &uAdjLoc);CHKERRQ(ierr);
459       ierr = DMGlobalToLocalBegin(dmAdj, uAdj, INSERT_VALUES, uAdjLoc);CHKERRQ(ierr);
460       ierr = DMGlobalToLocalEnd(dmAdj, uAdj, INSERT_VALUES, uAdjLoc);CHKERRQ(ierr);
461       ierr = DMGetGlobalVector(dm, &uAdjProj);CHKERRQ(ierr);
462       ierr = PetscObjectCompose((PetscObject) dm, "dmAux", (PetscObject) dmAdj);CHKERRQ(ierr);
463       ierr = PetscObjectCompose((PetscObject) dm, "A", (PetscObject) uAdjLoc);CHKERRQ(ierr);
464       ierr = DMProjectField(dm, 0.0, u, identity, INSERT_VALUES, uAdjProj);CHKERRQ(ierr);
465       ierr = PetscObjectCompose((PetscObject) dm, "dmAux", NULL);CHKERRQ(ierr);
466       ierr = PetscObjectCompose((PetscObject) dm, "A", NULL);CHKERRQ(ierr);
467       ierr = DMRestoreLocalVector(dmAdj, &uAdjLoc);CHKERRQ(ierr);
468       /*   Attach auxiliary data */
469       dms[0] = dm; dms[1] = dm;
470       ierr = DMCreateSuperDM(dms, 2, &subis, &dmErrAux);CHKERRQ(ierr);
471       if (0) {
472         PetscSection sec;
473 
474         ierr = DMGetLocalSection(dms[0], &sec);CHKERRQ(ierr);
475         ierr = PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
476         ierr = DMGetLocalSection(dms[1], &sec);CHKERRQ(ierr);
477         ierr = PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
478         ierr = DMGetLocalSection(dmErrAux, &sec);CHKERRQ(ierr);
479         ierr = PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
480       }
481       ierr = DMViewFromOptions(dmErrAux, NULL, "-dm_err_view");CHKERRQ(ierr);
482       ierr = ISViewFromOptions(subis[0], NULL, "-super_is_view");CHKERRQ(ierr);
483       ierr = ISViewFromOptions(subis[1], NULL, "-super_is_view");CHKERRQ(ierr);
484       ierr = DMGetGlobalVector(dmErrAux, &uErr);CHKERRQ(ierr);
485       ierr = VecViewFromOptions(u, NULL, "-map_vec_view");CHKERRQ(ierr);
486       ierr = VecViewFromOptions(uAdjProj, NULL, "-map_vec_view");CHKERRQ(ierr);
487       ierr = VecViewFromOptions(uErr, NULL, "-map_vec_view");CHKERRQ(ierr);
488       ierr = VecISCopy(uErr, subis[0], SCATTER_FORWARD, u);CHKERRQ(ierr);
489       ierr = VecISCopy(uErr, subis[1], SCATTER_FORWARD, uAdjProj);CHKERRQ(ierr);
490       ierr = DMRestoreGlobalVector(dm, &uAdjProj);CHKERRQ(ierr);
491       for (i = 0; i < 2; ++i) {ierr = ISDestroy(&subis[i]);CHKERRQ(ierr);}
492       ierr = PetscFree(subis);CHKERRQ(ierr);
493       ierr = DMGetLocalVector(dmErrAux, &uErrLoc);CHKERRQ(ierr);
494       ierr = DMGlobalToLocalBegin(dm, uErr, INSERT_VALUES, uErrLoc);CHKERRQ(ierr);
495       ierr = DMGlobalToLocalEnd(dm, uErr, INSERT_VALUES, uErrLoc);CHKERRQ(ierr);
496       ierr = DMRestoreGlobalVector(dmErrAux, &uErr);CHKERRQ(ierr);
497       ierr = PetscObjectCompose((PetscObject) dmAdj, "dmAux", (PetscObject) dmErrAux);CHKERRQ(ierr);
498       ierr = PetscObjectCompose((PetscObject) dmAdj, "A", (PetscObject) uErrLoc);CHKERRQ(ierr);
499       /*   Compute cellwise error estimate */
500       ierr = VecSet(errorEst, 0.0);CHKERRQ(ierr);
501       ierr = DMPlexComputeCellwiseIntegralFEM(dmAdj, uAdj, errorEst, &user);CHKERRQ(ierr);
502       ierr = PetscObjectCompose((PetscObject) dmAdj, "dmAux", NULL);CHKERRQ(ierr);
503       ierr = PetscObjectCompose((PetscObject) dmAdj, "A", NULL);CHKERRQ(ierr);
504       ierr = DMRestoreLocalVector(dmErrAux, &uErrLoc);CHKERRQ(ierr);
505       ierr = DMDestroy(&dmErrAux);CHKERRQ(ierr);
506       /*   Plot cellwise error vector */
507       ierr = VecViewFromOptions(errorEst, NULL, "-error_view");CHKERRQ(ierr);
508       /*   Compute ratio of estimate (sum over cells) with actual L_2 error */
509       ierr = DMComputeL2Diff(dm, 0.0, funcs, ctxs, u, &errorL2Norm);CHKERRQ(ierr);
510       ierr = DMPlexComputeL2DiffVec(dm, 0.0, funcs, ctxs, u, errorL2);CHKERRQ(ierr);
511       ierr = VecViewFromOptions(errorL2, NULL, "-l2_error_view");CHKERRQ(ierr);
512       ierr = VecNorm(errorL2,  NORM_INFINITY, &errorL2Tot);CHKERRQ(ierr);
513       ierr = VecNorm(errorEst, NORM_INFINITY, &errorEstTot);CHKERRQ(ierr);
514       ierr = VecGetSize(errorEst, &N);CHKERRQ(ierr);
515       ierr = VecPointwiseDivide(errorEst, errorEst, errorL2);CHKERRQ(ierr);
516       ierr = PetscObjectSetName((PetscObject) errorEst, "Error ratio");CHKERRQ(ierr);
517       ierr = VecViewFromOptions(errorEst, NULL, "-error_ratio_view");CHKERRQ(ierr);
518       ierr = PetscPrintf(PETSC_COMM_WORLD, "N: %D L2 error: %g Error Ratio: %g/%g = %g\n", N, (double) errorL2Norm, (double) errorEstTot, (double) PetscSqrtReal(errorL2Tot), (double) errorEstTot/PetscSqrtReal(errorL2Tot));CHKERRQ(ierr);
519       ierr = DMRestoreGlobalVector(dmErr, &errorEst);CHKERRQ(ierr);
520       ierr = DMRestoreGlobalVector(dmErr, &errorL2);CHKERRQ(ierr);
521       ierr = DMDestroy(&dmErr);CHKERRQ(ierr);
522     }
523     ierr = DMDestroy(&dmAdj);CHKERRQ(ierr);
524     ierr = VecDestroy(&uAdj);CHKERRQ(ierr);
525     ierr = SNESDestroy(&snesAdj);CHKERRQ(ierr);
526   }
527   /* Cleanup */
528   ierr = VecDestroy(&u);CHKERRQ(ierr);
529   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
530   ierr = DMDestroy(&dm);CHKERRQ(ierr);
531   ierr = PetscFinalize();
532   return ierr;
533 }
534 
535 /*TEST
536 
537   test:
538     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9
539     suffix: 2d_p1_conv
540     requires: triangle
541     args: -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2
542   test:
543     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9
544     suffix: 2d_p2_conv
545     requires: triangle
546     args: -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2
547   test:
548     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9
549     suffix: 2d_p3_conv
550     requires: triangle
551     args: -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2
552   test:
553     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9
554     suffix: 2d_q1_conv
555     args: -dm_plex_box_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2
556   test:
557     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9
558     suffix: 2d_q2_conv
559     args: -dm_plex_box_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2
560   test:
561     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9
562     suffix: 2d_q3_conv
563     args: -dm_plex_box_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2
564   test:
565     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9
566     suffix: 2d_q1_shear_conv
567     args: -dm_plex_box_simplex 0 -shear -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2
568   test:
569     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9
570     suffix: 2d_q2_shear_conv
571     args: -dm_plex_box_simplex 0 -shear -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2
572   test:
573     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9
574     suffix: 2d_q3_shear_conv
575     args: -dm_plex_box_simplex 0 -shear -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2
576   test:
577     # Using -convest_num_refine 3 we get L_2 convergence rate: 1.7
578     suffix: 3d_p1_conv
579     requires: ctetgen
580     args: -dm_plex_box_dim 3 -dm_refine 1 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1
581   test:
582     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 2.8
583     suffix: 3d_p2_conv
584     requires: ctetgen
585     args: -dm_plex_box_dim 3 -dm_plex_box_faces 2,2,2 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 1
586   test:
587     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 4.0
588     suffix: 3d_p3_conv
589     requires: ctetgen
590     args: -dm_plex_box_dim 3 -dm_plex_box_faces 2,2,2 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 1
591   test:
592     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.8
593     suffix: 3d_q1_conv
594     args: -dm_plex_box_dim 3 -dm_plex_box_simplex 0 -dm_refine 1 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1
595   test:
596     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.8
597     suffix: 3d_q2_conv
598     args: -dm_plex_box_dim 3 -dm_plex_box_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 1
599   test:
600     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 3.8
601     suffix: 3d_q3_conv
602     args: -dm_plex_box_dim 3 -dm_plex_box_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 1
603   test:
604     suffix: 2d_p1_fas_full
605     requires: triangle
606     args: -potential_petscspace_degree 1 -dm_refine_hierarchy 5 -dm_distribute \
607       -snes_max_it 1 -snes_type fas -snes_fas_levels 5 -snes_fas_type full -snes_fas_full_total \
608         -fas_coarse_snes_monitor -fas_coarse_snes_max_it 1 -fas_coarse_ksp_atol 1.e-13 \
609         -fas_levels_snes_monitor -fas_levels_snes_max_it 1 -fas_levels_snes_type newtonls \
610           -fas_levels_pc_type none -fas_levels_ksp_max_it 2 -fas_levels_ksp_converged_maxits -fas_levels_ksp_type chebyshev \
611             -fas_levels_esteig_ksp_type cg -fas_levels_ksp_chebyshev_esteig 0,0.25,0,1.05 -fas_levels_esteig_ksp_max_it 10
612   test:
613     suffix: 2d_p1_fas_full_homogeneous
614     requires: triangle
615     args: -homogeneous -potential_petscspace_degree 1 -dm_refine_hierarchy 5 -dm_distribute \
616       -snes_max_it 1 -snes_type fas -snes_fas_levels 5 -snes_fas_type full \
617         -fas_coarse_snes_monitor -fas_coarse_snes_max_it 1 -fas_coarse_ksp_atol 1.e-13 \
618         -fas_levels_snes_monitor -fas_levels_snes_max_it 1 -fas_levels_snes_type newtonls \
619           -fas_levels_pc_type none -fas_levels_ksp_max_it 2 -fas_levels_ksp_converged_maxits -fas_levels_ksp_type chebyshev \
620             -fas_levels_esteig_ksp_type cg -fas_levels_ksp_chebyshev_esteig 0,0.25,0,1.05 -fas_levels_esteig_ksp_max_it 10
621 
622   test:
623     suffix: 2d_p1_scalable
624     requires: triangle
625     args: -potential_petscspace_degree 1 -dm_refine 3 \
626       -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned \
627       -pc_type gamg \
628         -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \
629         -pc_gamg_coarse_eq_limit 1000 \
630         -pc_gamg_square_graph 1 \
631         -pc_gamg_threshold 0.05 \
632         -pc_gamg_threshold_scale .0 \
633       -mg_levels_ksp_type chebyshev \
634         -mg_levels_ksp_max_it 1 \
635         -mg_levels_esteig_ksp_type cg \
636         -mg_levels_esteig_ksp_max_it 10 \
637         -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \
638         -mg_levels_pc_type jacobi \
639       -matptap_via scalable
640   test:
641     suffix: 2d_p1_gmg_vcycle
642     requires: triangle
643     args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
644           -ksp_rtol 5e-10 -pc_type mg \
645             -mg_levels_ksp_max_it 1 \
646             -mg_levels_esteig_ksp_type cg \
647             -mg_levels_esteig_ksp_max_it 10 \
648             -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \
649             -mg_levels_pc_type jacobi
650   test:
651     suffix: 2d_p1_gmg_fcycle
652     requires: triangle
653     args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
654           -ksp_rtol 5e-10 -pc_type mg -pc_mg_type full \
655             -mg_levels_ksp_max_it 2 \
656             -mg_levels_esteig_ksp_type cg \
657             -mg_levels_esteig_ksp_max_it 10 \
658             -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \
659             -mg_levels_pc_type jacobi
660   test:
661     suffix: 2d_p1_gmg_vcycle_adapt
662     requires: triangle bamg
663     args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
664           -ksp_rtol 5e-10 -pc_type mg -pc_mg_galerkin -pc_mg_adapt_interp -pc_mg_adapt_interp_coarse_space harmonic -pc_mg_adapt_interp_n 8 \
665             -mg_levels_ksp_max_it 1 \
666             -mg_levels_esteig_ksp_type cg \
667             -mg_levels_esteig_ksp_max_it 10 \
668             -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \
669             -mg_levels_pc_type jacobi
670   test:
671     suffix: 2d_p1_spectral_0
672     requires: triangle fftw !complex
673     args: -dm_plex_box_faces 1,1 -potential_petscspace_degree 1 -dm_refine 6 -spectral -fft_view
674   test:
675     suffix: 2d_p1_spectral_1
676     requires: triangle fftw !complex
677     nsize: 2
678     args: -dm_plex_box_faces 4,4 -dm_distribute -potential_petscspace_degree 1 -spectral -fft_view
679   test:
680     suffix: 2d_p1_adj_0
681     requires: triangle
682     args: -potential_petscspace_degree 1 -dm_refine 1 -adjoint -adjoint_petscspace_degree 1 -error_petscspace_degree 0
683   test:
684     nsize: 2
685     requires: kokkos_kernels
686     suffix: kokkos
687     args: -dm_plex_box_dim 3 -dm_plex_box_faces 2,3,6 -dm_distribute -petscpartitioner_type simple -dm_plex_box_simplex 0 -potential_petscspace_degree 1 \
688          -dm_refine 0 -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned -pc_type gamg -pc_gamg_coarse_eq_limit 1000 -pc_gamg_threshold 0.0 \
689          -pc_gamg_threshold_scale .5 -mg_levels_ksp_type chebyshev -mg_levels_ksp_max_it 2 -mg_levels_esteig_ksp_type cg -mg_levels_esteig_ksp_max_it 10 \
690          -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 -mg_levels_pc_type jacobi -ksp_monitor -snes_monitor -dm_view -dm_mat_type aijkokkos -dm_vec_type kokkos
691 
692 TEST*/
693