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