xref: /petsc/src/snes/tutorials/ex13.c (revision 4e278199b78715991f5c71ebbd945c1489263e6c)
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   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
173   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
174   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
175   if (user->shear) {ierr = DMPlexShearGeometry(*dm, DM_X, NULL);CHKERRQ(ierr);}
176   ierr = DMSetApplicationContext(*dm, user);CHKERRQ(ierr);
177   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
178   if (user->spectral) {
179     PetscInt  planeDir[2]   = {0,  1};
180     PetscReal planeCoord[2] = {0., 1.};
181 
182     ierr = CreateSpectralPlanes(*dm, 2, planeDir, planeCoord, user);CHKERRQ(ierr);
183   }
184   PetscFunctionReturn(0);
185 }
186 
187 static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
188 {
189   PetscDS        ds;
190   DMLabel        label;
191   const PetscInt id = 1;
192   PetscPointFunc f0 = user->homogeneous ? f0_trig_homogeneous_u : f0_trig_inhomogeneous_u;
193   PetscErrorCode (*ex)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *) = user->homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u;
194   PetscErrorCode ierr;
195 
196   PetscFunctionBeginUser;
197   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
198   ierr = PetscDSSetResidual(ds, 0, f0, f1_u);CHKERRQ(ierr);
199   ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr);
200   ierr = PetscDSSetExactSolution(ds, 0, ex, user);CHKERRQ(ierr);
201   ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr);
202   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) ex, NULL, user, NULL);CHKERRQ(ierr);
203   PetscFunctionReturn(0);
204 }
205 
206 static PetscErrorCode SetupAdjointProblem(DM dm, AppCtx *user)
207 {
208   PetscDS        ds;
209   DMLabel        label;
210   const PetscInt id = 1;
211   PetscErrorCode ierr;
212 
213   PetscFunctionBeginUser;
214   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
215   ierr = PetscDSSetResidual(ds, 0, f0_unity_u, f1_u);CHKERRQ(ierr);
216   ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr);
217   ierr = PetscDSSetObjective(ds, 0, obj_error_u);CHKERRQ(ierr);
218   ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr);
219   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) zero, NULL, user, NULL);CHKERRQ(ierr);
220   PetscFunctionReturn(0);
221 }
222 
223 static PetscErrorCode SetupErrorProblem(DM dm, AppCtx *user)
224 {
225   PetscDS        prob;
226   PetscErrorCode ierr;
227 
228   PetscFunctionBeginUser;
229   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
230   PetscFunctionReturn(0);
231 }
232 
233 static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
234 {
235   DM             cdm = dm;
236   PetscFE        fe;
237   DMPolytopeType ct;
238   PetscBool      simplex;
239   PetscInt       dim, cStart;
240   char           prefix[PETSC_MAX_PATH_LEN];
241   PetscErrorCode ierr;
242 
243   PetscFunctionBeginUser;
244   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
245   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr);
246   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
247   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
248   /* Create finite element */
249   ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name);CHKERRQ(ierr);
250   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe);CHKERRQ(ierr);
251   ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr);
252   /* Set discretization and boundary conditions for each mesh */
253   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr);
254   ierr = DMCreateDS(dm);CHKERRQ(ierr);
255   ierr = (*setup)(dm, user);CHKERRQ(ierr);
256   while (cdm) {
257     ierr = DMCopyDisc(dm,cdm);CHKERRQ(ierr);
258     /* TODO: Check whether the boundary of coarse meshes is marked */
259     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
260   }
261   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
262   PetscFunctionReturn(0);
263 }
264 
265 static PetscErrorCode ComputeSpectral(DM dm, Vec u, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user)
266 {
267   MPI_Comm           comm;
268   PetscSection       coordSection, section;
269   Vec                coordinates, uloc;
270   const PetscScalar *coords, *array;
271   PetscInt           p;
272   PetscMPIInt        size, rank;
273   PetscErrorCode     ierr;
274 
275   PetscFunctionBeginUser;
276   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
277   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
278   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
279   ierr = DMGetLocalVector(dm, &uloc);CHKERRQ(ierr);
280   ierr = DMGlobalToLocalBegin(dm, u, INSERT_VALUES, uloc);CHKERRQ(ierr);
281   ierr = DMGlobalToLocalEnd(dm, u, INSERT_VALUES, uloc);CHKERRQ(ierr);
282   ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, uloc, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
283   ierr = VecViewFromOptions(uloc, NULL, "-sol_view");CHKERRQ(ierr);
284   ierr = DMGetLocalSection(dm, &section);CHKERRQ(ierr);
285   ierr = VecGetArrayRead(uloc, &array);CHKERRQ(ierr);
286   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
287   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
288   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
289   for (p = 0; p < numPlanes; ++p) {
290     DMLabel         label;
291     char            name[PETSC_MAX_PATH_LEN];
292     Mat             F;
293     Vec             x, y;
294     IS              stratum;
295     PetscReal      *ray, *gray;
296     PetscScalar    *rvals, *svals, *gsvals;
297     PetscInt       *perm, *nperm;
298     PetscInt        n, N, i, j, off, offu;
299     const PetscInt *points;
300 
301     ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%D", p);CHKERRQ(ierr);
302     ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr);
303     ierr = DMLabelGetStratumIS(label, 1, &stratum);CHKERRQ(ierr);
304     ierr = ISGetLocalSize(stratum, &n);CHKERRQ(ierr);
305     ierr = ISGetIndices(stratum, &points);CHKERRQ(ierr);
306     ierr = PetscMalloc2(n, &ray, n, &svals);CHKERRQ(ierr);
307     for (i = 0; i < n; ++i) {
308       ierr = PetscSectionGetOffset(coordSection, points[i], &off);CHKERRQ(ierr);
309       ierr = PetscSectionGetOffset(section, points[i], &offu);CHKERRQ(ierr);
310       ray[i]   = PetscRealPart(coords[off+((planeDir[p]+1)%2)]);
311       svals[i] = array[offu];
312     }
313     /* Gather the ray data to proc 0 */
314     if (size > 1) {
315       PetscMPIInt *cnt, *displs, p;
316 
317       ierr = PetscCalloc2(size, &cnt, size, &displs);CHKERRQ(ierr);
318       ierr = MPI_Gather(&n, 1, MPIU_INT, cnt, 1, MPIU_INT, 0, comm);CHKERRMPI(ierr);
319       for (p = 1; p < size; ++p) displs[p] = displs[p-1] + cnt[p-1];
320       N = displs[size-1] + cnt[size-1];
321       ierr = PetscMalloc2(N, &gray, N, &gsvals);CHKERRQ(ierr);
322       ierr = MPI_Gatherv(ray, n, MPIU_REAL, gray, cnt, displs, MPIU_REAL, 0, comm);CHKERRMPI(ierr);
323       ierr = MPI_Gatherv(svals, n, MPIU_SCALAR, gsvals, cnt, displs, MPIU_SCALAR, 0, comm);CHKERRMPI(ierr);
324       ierr = PetscFree2(cnt, displs);CHKERRQ(ierr);
325     } else {
326       N      = n;
327       gray   = ray;
328       gsvals = svals;
329     }
330     if (!rank) {
331       /* Sort point along ray */
332       ierr = PetscMalloc2(N, &perm, N, &nperm);CHKERRQ(ierr);
333       for (i = 0; i < N; ++i) {perm[i] = i;}
334       ierr = PetscSortRealWithPermutation(N, gray, perm);CHKERRQ(ierr);
335       /* Count duplicates and squish mapping */
336       nperm[0] = perm[0];
337       for (i = 1, j = 1; i < N; ++i) {
338         if (PetscAbsReal(gray[perm[i]] - gray[perm[i-1]]) > PETSC_SMALL) nperm[j++] = perm[i];
339       }
340       /* Create FFT structs */
341       ierr = MatCreateFFT(PETSC_COMM_SELF, 1, &j, MATFFTW, &F);CHKERRQ(ierr);
342       ierr = MatCreateVecs(F, &x, &y);CHKERRQ(ierr);
343       ierr = PetscObjectSetName((PetscObject) y, name);CHKERRQ(ierr);
344       ierr = VecGetArray(x, &rvals);CHKERRQ(ierr);
345       for (i = 0, j = 0; i < N; ++i) {
346         if (i > 0 && PetscAbsReal(gray[perm[i]] - gray[perm[i-1]]) < PETSC_SMALL) continue;
347         rvals[j] = gsvals[nperm[j]];
348         ++j;
349       }
350       ierr = PetscFree2(perm, nperm);CHKERRQ(ierr);
351       if (size > 1) {ierr = PetscFree2(gray, gsvals);CHKERRQ(ierr);}
352       ierr = VecRestoreArray(x, &rvals);CHKERRQ(ierr);
353       /* Do FFT along the ray */
354       ierr = MatMult(F, x, y);CHKERRQ(ierr);
355       /* Chop FFT */
356       ierr = VecChop(y, PETSC_SMALL);CHKERRQ(ierr);
357       ierr = VecViewFromOptions(x, NULL, "-real_view");CHKERRQ(ierr);
358       ierr = VecViewFromOptions(y, NULL, "-fft_view");CHKERRQ(ierr);
359       ierr = VecDestroy(&x);CHKERRQ(ierr);
360       ierr = VecDestroy(&y);CHKERRQ(ierr);
361       ierr = MatDestroy(&F);CHKERRQ(ierr);
362     }
363     ierr = ISRestoreIndices(stratum, &points);CHKERRQ(ierr);
364     ierr = ISDestroy(&stratum);CHKERRQ(ierr);
365     ierr = PetscFree2(ray, svals);CHKERRQ(ierr);
366   }
367   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
368   ierr = VecRestoreArrayRead(uloc, &array);CHKERRQ(ierr);
369   ierr = DMRestoreLocalVector(dm, &uloc);CHKERRQ(ierr);
370   PetscFunctionReturn(0);
371 }
372 
373 int main(int argc, char **argv)
374 {
375   DM             dm;   /* Problem specification */
376   SNES           snes; /* Nonlinear solver */
377   Vec            u;    /* Solutions */
378   AppCtx         user; /* User-defined work context */
379   PetscErrorCode ierr;
380 
381   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
382   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
383   /* Primal system */
384   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
385   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
386   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
387   ierr = SetupDiscretization(dm, "potential", SetupPrimalProblem, &user);CHKERRQ(ierr);
388   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
389   ierr = VecSet(u, 0.0);CHKERRQ(ierr);
390   ierr = PetscObjectSetName((PetscObject) u, "potential");CHKERRQ(ierr);
391   ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr);
392   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
393   ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
394   ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
395   ierr = VecViewFromOptions(u, NULL, "-potential_view");CHKERRQ(ierr);
396   if (user.spectral) {
397     PetscInt  planeDir[2]   = {0,  1};
398     PetscReal planeCoord[2] = {0., 1.};
399 
400     ierr = ComputeSpectral(dm, u, 2, planeDir, planeCoord, &user);CHKERRQ(ierr);
401   }
402   /* Adjoint system */
403   if (user.adjoint) {
404     DM   dmAdj;
405     SNES snesAdj;
406     Vec  uAdj;
407 
408     ierr = SNESCreate(PETSC_COMM_WORLD, &snesAdj);CHKERRQ(ierr);
409     ierr = PetscObjectSetOptionsPrefix((PetscObject) snesAdj, "adjoint_");CHKERRQ(ierr);
410     ierr = DMClone(dm, &dmAdj);CHKERRQ(ierr);
411     ierr = SNESSetDM(snesAdj, dmAdj);CHKERRQ(ierr);
412     ierr = SetupDiscretization(dmAdj, "adjoint", SetupAdjointProblem, &user);CHKERRQ(ierr);
413     ierr = DMCreateGlobalVector(dmAdj, &uAdj);CHKERRQ(ierr);
414     ierr = VecSet(uAdj, 0.0);CHKERRQ(ierr);
415     ierr = PetscObjectSetName((PetscObject) uAdj, "adjoint");CHKERRQ(ierr);
416     ierr = DMPlexSetSNESLocalFEM(dmAdj, &user, &user, &user);CHKERRQ(ierr);
417     ierr = SNESSetFromOptions(snesAdj);CHKERRQ(ierr);
418     ierr = SNESSolve(snesAdj, NULL, uAdj);CHKERRQ(ierr);
419     ierr = SNESGetSolution(snesAdj, &uAdj);CHKERRQ(ierr);
420     ierr = VecViewFromOptions(uAdj, NULL, "-adjoint_view");CHKERRQ(ierr);
421     /* Error representation */
422     {
423       DM        dmErr, dmErrAux, dms[2];
424       Vec       errorEst, errorL2, uErr, uErrLoc, uAdjLoc, uAdjProj;
425       IS       *subis;
426       PetscReal errorEstTot, errorL2Norm, errorL2Tot;
427       PetscInt  N, i;
428       PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *) = {user.homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u};
429       void (*identity[1])(PetscInt, PetscInt, PetscInt,
430                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
431                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
432                           PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]) = {f0_identityaux_u};
433       void            *ctxs[1] = {0};
434 
435       ctxs[0] = &user;
436       ierr = DMClone(dm, &dmErr);CHKERRQ(ierr);
437       ierr = SetupDiscretization(dmErr, "error", SetupErrorProblem, &user);CHKERRQ(ierr);
438       ierr = DMGetGlobalVector(dmErr, &errorEst);CHKERRQ(ierr);
439       ierr = DMGetGlobalVector(dmErr, &errorL2);CHKERRQ(ierr);
440       /*   Compute auxiliary data (solution and projection of adjoint solution) */
441       ierr = DMGetLocalVector(dmAdj, &uAdjLoc);CHKERRQ(ierr);
442       ierr = DMGlobalToLocalBegin(dmAdj, uAdj, INSERT_VALUES, uAdjLoc);CHKERRQ(ierr);
443       ierr = DMGlobalToLocalEnd(dmAdj, uAdj, INSERT_VALUES, uAdjLoc);CHKERRQ(ierr);
444       ierr = DMGetGlobalVector(dm, &uAdjProj);CHKERRQ(ierr);
445       ierr = DMSetAuxiliaryVec(dm, NULL, 0, uAdjLoc);CHKERRQ(ierr);
446       ierr = DMProjectField(dm, 0.0, u, identity, INSERT_VALUES, uAdjProj);CHKERRQ(ierr);
447       ierr = DMSetAuxiliaryVec(dm, NULL, 0, NULL);CHKERRQ(ierr);
448       ierr = DMRestoreLocalVector(dmAdj, &uAdjLoc);CHKERRQ(ierr);
449       /*   Attach auxiliary data */
450       dms[0] = dm; dms[1] = dm;
451       ierr = DMCreateSuperDM(dms, 2, &subis, &dmErrAux);CHKERRQ(ierr);
452       if (0) {
453         PetscSection sec;
454 
455         ierr = DMGetLocalSection(dms[0], &sec);CHKERRQ(ierr);
456         ierr = PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
457         ierr = DMGetLocalSection(dms[1], &sec);CHKERRQ(ierr);
458         ierr = PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
459         ierr = DMGetLocalSection(dmErrAux, &sec);CHKERRQ(ierr);
460         ierr = PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
461       }
462       ierr = DMViewFromOptions(dmErrAux, NULL, "-dm_err_view");CHKERRQ(ierr);
463       ierr = ISViewFromOptions(subis[0], NULL, "-super_is_view");CHKERRQ(ierr);
464       ierr = ISViewFromOptions(subis[1], NULL, "-super_is_view");CHKERRQ(ierr);
465       ierr = DMGetGlobalVector(dmErrAux, &uErr);CHKERRQ(ierr);
466       ierr = VecViewFromOptions(u, NULL, "-map_vec_view");CHKERRQ(ierr);
467       ierr = VecViewFromOptions(uAdjProj, NULL, "-map_vec_view");CHKERRQ(ierr);
468       ierr = VecViewFromOptions(uErr, NULL, "-map_vec_view");CHKERRQ(ierr);
469       ierr = VecISCopy(uErr, subis[0], SCATTER_FORWARD, u);CHKERRQ(ierr);
470       ierr = VecISCopy(uErr, subis[1], SCATTER_FORWARD, uAdjProj);CHKERRQ(ierr);
471       ierr = DMRestoreGlobalVector(dm, &uAdjProj);CHKERRQ(ierr);
472       for (i = 0; i < 2; ++i) {ierr = ISDestroy(&subis[i]);CHKERRQ(ierr);}
473       ierr = PetscFree(subis);CHKERRQ(ierr);
474       ierr = DMGetLocalVector(dmErrAux, &uErrLoc);CHKERRQ(ierr);
475       ierr = DMGlobalToLocalBegin(dm, uErr, INSERT_VALUES, uErrLoc);CHKERRQ(ierr);
476       ierr = DMGlobalToLocalEnd(dm, uErr, INSERT_VALUES, uErrLoc);CHKERRQ(ierr);
477       ierr = DMRestoreGlobalVector(dmErrAux, &uErr);CHKERRQ(ierr);
478       ierr = DMSetAuxiliaryVec(dmAdj, NULL, 0, uErrLoc);CHKERRQ(ierr);
479       /*   Compute cellwise error estimate */
480       ierr = VecSet(errorEst, 0.0);CHKERRQ(ierr);
481       ierr = DMPlexComputeCellwiseIntegralFEM(dmAdj, uAdj, errorEst, &user);CHKERRQ(ierr);
482       ierr = DMSetAuxiliaryVec(dmAdj, NULL, 0, NULL);CHKERRQ(ierr);
483       ierr = DMRestoreLocalVector(dmErrAux, &uErrLoc);CHKERRQ(ierr);
484       ierr = DMDestroy(&dmErrAux);CHKERRQ(ierr);
485       /*   Plot cellwise error vector */
486       ierr = VecViewFromOptions(errorEst, NULL, "-error_view");CHKERRQ(ierr);
487       /*   Compute ratio of estimate (sum over cells) with actual L_2 error */
488       ierr = DMComputeL2Diff(dm, 0.0, funcs, ctxs, u, &errorL2Norm);CHKERRQ(ierr);
489       ierr = DMPlexComputeL2DiffVec(dm, 0.0, funcs, ctxs, u, errorL2);CHKERRQ(ierr);
490       ierr = VecViewFromOptions(errorL2, NULL, "-l2_error_view");CHKERRQ(ierr);
491       ierr = VecNorm(errorL2,  NORM_INFINITY, &errorL2Tot);CHKERRQ(ierr);
492       ierr = VecNorm(errorEst, NORM_INFINITY, &errorEstTot);CHKERRQ(ierr);
493       ierr = VecGetSize(errorEst, &N);CHKERRQ(ierr);
494       ierr = VecPointwiseDivide(errorEst, errorEst, errorL2);CHKERRQ(ierr);
495       ierr = PetscObjectSetName((PetscObject) errorEst, "Error ratio");CHKERRQ(ierr);
496       ierr = VecViewFromOptions(errorEst, NULL, "-error_ratio_view");CHKERRQ(ierr);
497       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);
498       ierr = DMRestoreGlobalVector(dmErr, &errorEst);CHKERRQ(ierr);
499       ierr = DMRestoreGlobalVector(dmErr, &errorL2);CHKERRQ(ierr);
500       ierr = DMDestroy(&dmErr);CHKERRQ(ierr);
501     }
502     ierr = DMDestroy(&dmAdj);CHKERRQ(ierr);
503     ierr = VecDestroy(&uAdj);CHKERRQ(ierr);
504     ierr = SNESDestroy(&snesAdj);CHKERRQ(ierr);
505   }
506   /* Cleanup */
507   ierr = VecDestroy(&u);CHKERRQ(ierr);
508   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
509   ierr = DMDestroy(&dm);CHKERRQ(ierr);
510   ierr = PetscFinalize();
511   return ierr;
512 }
513 
514 /*TEST
515 
516   test:
517     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9
518     suffix: 2d_p1_conv
519     requires: triangle
520     args: -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2
521   test:
522     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9
523     suffix: 2d_p2_conv
524     requires: triangle
525     args: -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2
526   test:
527     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9
528     suffix: 2d_p3_conv
529     requires: triangle
530     args: -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2
531   test:
532     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9
533     suffix: 2d_q1_conv
534     args: -dm_plex_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2
535   test:
536     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9
537     suffix: 2d_q2_conv
538     args: -dm_plex_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2
539   test:
540     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9
541     suffix: 2d_q3_conv
542     args: -dm_plex_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2
543   test:
544     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9
545     suffix: 2d_q1_shear_conv
546     args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2
547   test:
548     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9
549     suffix: 2d_q2_shear_conv
550     args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2
551   test:
552     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9
553     suffix: 2d_q3_shear_conv
554     args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2
555   test:
556     # Using -convest_num_refine 3 we get L_2 convergence rate: 1.7
557     suffix: 3d_p1_conv
558     requires: ctetgen
559     args: -dm_plex_dim 3 -dm_refine 1 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1
560   test:
561     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 2.8
562     suffix: 3d_p2_conv
563     requires: ctetgen
564     args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 1
565   test:
566     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 4.0
567     suffix: 3d_p3_conv
568     requires: ctetgen
569     args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 1
570   test:
571     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.8
572     suffix: 3d_q1_conv
573     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_refine 1 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1
574   test:
575     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.8
576     suffix: 3d_q2_conv
577     args: -dm_plex_dim 3 -dm_plex_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 1
578   test:
579     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 3.8
580     suffix: 3d_q3_conv
581     args: -dm_plex_dim 3 -dm_plex_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 1
582   test:
583     suffix: 2d_p1_fas_full
584     requires: triangle
585     args: -potential_petscspace_degree 1 -dm_refine_hierarchy 5 -dm_distribute \
586       -snes_max_it 1 -snes_type fas -snes_fas_levels 5 -snes_fas_type full -snes_fas_full_total \
587         -fas_coarse_snes_monitor -fas_coarse_snes_max_it 1 -fas_coarse_ksp_atol 1.e-13 \
588         -fas_levels_snes_monitor -fas_levels_snes_max_it 1 -fas_levels_snes_type newtonls \
589           -fas_levels_pc_type none -fas_levels_ksp_max_it 2 -fas_levels_ksp_converged_maxits -fas_levels_ksp_type chebyshev \
590             -fas_levels_esteig_ksp_type cg -fas_levels_ksp_chebyshev_esteig 0,0.25,0,1.05 -fas_levels_esteig_ksp_max_it 10
591   test:
592     suffix: 2d_p1_fas_full_homogeneous
593     requires: triangle
594     args: -homogeneous -potential_petscspace_degree 1 -dm_refine_hierarchy 5 -dm_distribute \
595       -snes_max_it 1 -snes_type fas -snes_fas_levels 5 -snes_fas_type full \
596         -fas_coarse_snes_monitor -fas_coarse_snes_max_it 1 -fas_coarse_ksp_atol 1.e-13 \
597         -fas_levels_snes_monitor -fas_levels_snes_max_it 1 -fas_levels_snes_type newtonls \
598           -fas_levels_pc_type none -fas_levels_ksp_max_it 2 -fas_levels_ksp_converged_maxits -fas_levels_ksp_type chebyshev \
599             -fas_levels_esteig_ksp_type cg -fas_levels_ksp_chebyshev_esteig 0,0.25,0,1.05 -fas_levels_esteig_ksp_max_it 10
600 
601   test:
602     suffix: 2d_p1_scalable
603     requires: triangle
604     args: -potential_petscspace_degree 1 -dm_refine 3 \
605       -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned \
606       -pc_type gamg \
607         -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \
608         -pc_gamg_coarse_eq_limit 1000 \
609         -pc_gamg_square_graph 1 \
610         -pc_gamg_threshold 0.05 \
611         -pc_gamg_threshold_scale .0 \
612       -mg_levels_ksp_type chebyshev \
613         -mg_levels_ksp_max_it 1 \
614         -mg_levels_esteig_ksp_type cg \
615         -mg_levels_esteig_ksp_max_it 10 \
616         -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \
617         -mg_levels_pc_type jacobi \
618       -matptap_via scalable
619   test:
620     suffix: 2d_p1_gmg_vcycle
621     requires: triangle
622     args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
623           -ksp_rtol 5e-10 -pc_type mg \
624             -mg_levels_ksp_max_it 1 \
625             -mg_levels_esteig_ksp_type cg \
626             -mg_levels_esteig_ksp_max_it 10 \
627             -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \
628             -mg_levels_pc_type jacobi
629   test:
630     suffix: 2d_p1_gmg_fcycle
631     requires: triangle
632     args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
633           -ksp_rtol 5e-10 -pc_type mg -pc_mg_type full \
634             -mg_levels_ksp_max_it 2 \
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   test:
640     suffix: 2d_p1_gmg_vcycle_adapt
641     requires: triangle bamg
642     args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
643           -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 \
644             -mg_levels_ksp_max_it 1 \
645             -mg_levels_esteig_ksp_type cg \
646             -mg_levels_esteig_ksp_max_it 10 \
647             -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \
648             -mg_levels_pc_type jacobi
649   test:
650     suffix: 2d_p1_spectral_0
651     requires: triangle fftw !complex
652     args: -dm_plex_box_faces 1,1 -potential_petscspace_degree 1 -dm_refine 6 -spectral -fft_view
653   test:
654     suffix: 2d_p1_spectral_1
655     requires: triangle fftw !complex
656     nsize: 2
657     args: -dm_plex_box_faces 4,4 -dm_distribute -potential_petscspace_degree 1 -spectral -fft_view
658   test:
659     suffix: 2d_p1_adj_0
660     requires: triangle
661     args: -potential_petscspace_degree 1 -dm_refine 1 -adjoint -adjoint_petscspace_degree 1 -error_petscspace_degree 0
662   test:
663     nsize: 2
664     requires: kokkos_kernels
665     suffix: kokkos
666     args: -dm_plex_dim 3 -dm_plex_box_faces 2,3,6 -dm_distribute -petscpartitioner_type simple -dm_plex_simplex 0 -potential_petscspace_degree 1 \
667          -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 \
668          -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 \
669          -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
670 
671 TEST*/
672