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