Actual source code: ex13.c

  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";

  7: #include <petscdmplex.h>
  8: #include <petscdmceed.h>
  9: #include <petscsnes.h>
 10: #include <petscds.h>
 11: #include <petscconvest.h>

 13: typedef struct {
 14:   /* Domain and mesh definition */
 15:   PetscBool spectral;    /* Look at the spectrum along planes in the solution */
 16:   PetscBool shear;       /* Shear the domain */
 17:   PetscBool adjoint;     /* Solve the adjoint problem */
 18:   PetscBool homogeneous; /* Use homogeneous boundary conditions */
 19:   PetscBool viewError;   /* Output the solution error */
 20: } AppCtx;

 22: static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
 23: {
 24:   *u = 0.0;
 25:   return PETSC_SUCCESS;
 26: }

 28: static PetscErrorCode trig_inhomogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
 29: {
 30:   PetscInt d;
 31:   *u = 0.0;
 32:   for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0 * PETSC_PI * x[d]);
 33:   return PETSC_SUCCESS;
 34: }

 36: static PetscErrorCode trig_homogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
 37: {
 38:   PetscInt d;
 39:   *u = 1.0;
 40:   for (d = 0; d < dim; ++d) *u *= PetscSinReal(2.0 * PETSC_PI * x[d]);
 41:   return PETSC_SUCCESS;
 42: }

 44: /* Compute integral of (residual of solution)*(adjoint solution - projection of adjoint solution) */
 45: 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[])
 46: {
 47:   obj[0] = a[aOff[0]] * (u[0] - a[aOff[1]]);
 48: }

 50: 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[])
 51: {
 52:   PetscInt d;
 53:   for (d = 0; d < dim; ++d) f0[0] += -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]);
 54: }

 56: 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[])
 57: {
 58:   for (PetscInt d = 0; d < dim; ++d) {
 59:     PetscScalar v = 1.;
 60:     for (PetscInt e = 0; e < dim; e++) {
 61:       if (e == d) {
 62:         v *= -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]);
 63:       } else {
 64:         v *= PetscSinReal(2.0 * PETSC_PI * x[d]);
 65:       }
 66:     }
 67:     f0[0] += v;
 68:   }
 69: }

 71: 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[])
 72: {
 73:   f0[0] = 1.0;
 74: }

 76: 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[])
 77: {
 78:   f0[0] = a[0];
 79: }

 81: 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[])
 82: {
 83:   for (PetscInt d = 0; d < dim; ++d) f1[d] = u_x[d];
 84: }

 86: 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[])
 87: {
 88:   for (PetscInt d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
 89: }

 91: PLEXFE_QFUNCTION(Laplace, f0_trig_inhomogeneous_u, f1_u)

 93: static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
 94: {
 95:   PetscFunctionBeginUser;
 96:   options->shear       = PETSC_FALSE;
 97:   options->spectral    = PETSC_FALSE;
 98:   options->adjoint     = PETSC_FALSE;
 99:   options->homogeneous = PETSC_FALSE;
100:   options->viewError   = PETSC_FALSE;

102:   PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");
103:   PetscCall(PetscOptionsBool("-shear", "Shear the domain", "ex13.c", options->shear, &options->shear, NULL));
104:   PetscCall(PetscOptionsBool("-spectral", "Look at the spectrum along planes of the solution", "ex13.c", options->spectral, &options->spectral, NULL));
105:   PetscCall(PetscOptionsBool("-adjoint", "Solve the adjoint problem", "ex13.c", options->adjoint, &options->adjoint, NULL));
106:   PetscCall(PetscOptionsBool("-homogeneous", "Use homogeneous boundary conditions", "ex13.c", options->homogeneous, &options->homogeneous, NULL));
107:   PetscCall(PetscOptionsBool("-error_view", "Output the solution error", "ex13.c", options->viewError, &options->viewError, NULL));
108:   PetscOptionsEnd();
109:   PetscFunctionReturn(PETSC_SUCCESS);
110: }

112: static PetscErrorCode CreateSpectralPlanes(DM dm, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user)
113: {
114:   PetscSection       coordSection;
115:   Vec                coordinates;
116:   const PetscScalar *coords;
117:   PetscInt           dim, p, vStart, vEnd, v;

119:   PetscFunctionBeginUser;
120:   PetscCall(DMGetCoordinateDim(dm, &dim));
121:   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
122:   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
123:   PetscCall(DMGetCoordinateSection(dm, &coordSection));
124:   PetscCall(VecGetArrayRead(coordinates, &coords));
125:   for (p = 0; p < numPlanes; ++p) {
126:     DMLabel label;
127:     char    name[PETSC_MAX_PATH_LEN];

129:     PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%" PetscInt_FMT, p));
130:     PetscCall(DMCreateLabel(dm, name));
131:     PetscCall(DMGetLabel(dm, name, &label));
132:     PetscCall(DMLabelAddStratum(label, 1));
133:     for (v = vStart; v < vEnd; ++v) {
134:       PetscInt off;

136:       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
137:       if (PetscAbsReal(planeCoord[p] - PetscRealPart(coords[off + planeDir[p]])) < PETSC_SMALL) PetscCall(DMLabelSetValue(label, v, 1));
138:     }
139:   }
140:   PetscCall(VecRestoreArrayRead(coordinates, &coords));
141:   PetscFunctionReturn(PETSC_SUCCESS);
142: }

144: static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
145: {
146:   PetscFunctionBeginUser;
147:   PetscCall(DMCreate(comm, dm));
148:   PetscCall(DMSetType(*dm, DMPLEX));
149:   PetscCall(DMSetFromOptions(*dm));
150:   if (user->shear) PetscCall(DMPlexShearGeometry(*dm, DM_X, NULL));
151:   PetscCall(DMSetApplicationContext(*dm, user));
152:   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
153:   if (user->spectral) {
154:     PetscInt  planeDir[2]   = {0, 1};
155:     PetscReal planeCoord[2] = {0., 1.};

157:     PetscCall(CreateSpectralPlanes(*dm, 2, planeDir, planeCoord, user));
158:   }
159:   PetscFunctionReturn(PETSC_SUCCESS);
160: }

162: static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
163: {
164:   PetscDS        ds;
165:   DMLabel        label;
166:   const PetscInt id                                                                             = 1;
167:   PetscPointFn  *f0                                                                             = user->homogeneous ? f0_trig_homogeneous_u : f0_trig_inhomogeneous_u;
168:   PetscErrorCode (*ex)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = user->homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u;

170:   PetscFunctionBeginUser;
171:   PetscCall(DMGetDS(dm, &ds));
172:   PetscCall(PetscDSSetResidual(ds, 0, f0, f1_u));
173:   PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu));
174:   PetscCall(PetscDSSetExactSolution(ds, 0, ex, user));
175:   PetscCall(DMGetLabel(dm, "marker", &label));
176:   if (label) PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)ex, NULL, user, NULL));
177:   PetscFunctionReturn(PETSC_SUCCESS);
178: }

180: static PetscErrorCode SetupAdjointProblem(DM dm, AppCtx *user)
181: {
182:   PetscDS        ds;
183:   DMLabel        label;
184:   const PetscInt id = 1;

186:   PetscFunctionBeginUser;
187:   PetscCall(DMGetDS(dm, &ds));
188:   PetscCall(PetscDSSetResidual(ds, 0, f0_unity_u, f1_u));
189:   PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu));
190:   PetscCall(PetscDSSetObjective(ds, 0, obj_error_u));
191:   PetscCall(DMGetLabel(dm, "marker", &label));
192:   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)zero, NULL, user, NULL));
193:   PetscFunctionReturn(PETSC_SUCCESS);
194: }

196: static PetscErrorCode SetupErrorProblem(DM dm, AppCtx *user)
197: {
198:   PetscDS prob;

200:   PetscFunctionBeginUser;
201:   PetscCall(DMGetDS(dm, &prob));
202:   PetscFunctionReturn(PETSC_SUCCESS);
203: }

205: static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
206: {
207:   DM             cdm = dm;
208:   PetscFE        fe;
209:   DMPolytopeType ct;
210:   PetscBool      simplex;
211:   PetscInt       dim, cStart;
212:   char           prefix[PETSC_MAX_PATH_LEN];

214:   PetscFunctionBeginUser;
215:   PetscCall(DMGetDimension(dm, &dim));
216:   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
217:   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
218:   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
219:   /* Create finite element */
220:   PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
221:   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe));
222:   PetscCall(PetscObjectSetName((PetscObject)fe, name));
223:   /* Set discretization and boundary conditions for each mesh */
224:   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
225:   PetscCall(DMCreateDS(dm));
226:   PetscCall((*setup)(dm, user));
227:   while (cdm) {
228:     PetscCall(DMCopyDisc(dm, cdm));
229:     /* TODO: Check whether the boundary of coarse meshes is marked */
230:     PetscCall(DMGetCoarseDM(cdm, &cdm));
231:   }
232:   PetscCall(PetscFEDestroy(&fe));
233: #ifdef PETSC_HAVE_LIBCEED
234:   PetscBool useCeed;
235:   PetscCall(DMPlexGetUseCeed(dm, &useCeed));
236:   if (useCeed) PetscCall(DMCeedCreate(dm, PETSC_TRUE, PlexQFunctionLaplace, PlexQFunctionLaplace_loc));
237: #endif
238:   PetscFunctionReturn(PETSC_SUCCESS);
239: }

241: static PetscErrorCode ComputeSpectral(Vec u, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user)
242: {
243:   MPI_Comm           comm;
244:   DM                 dm;
245:   PetscSection       coordSection, section;
246:   Vec                coordinates, uloc;
247:   const PetscScalar *coords, *array;
248:   PetscMPIInt        size, rank;

250:   PetscFunctionBeginUser;
251:   if (!user->spectral) PetscFunctionReturn(PETSC_SUCCESS);
252:   PetscCall(VecGetDM(u, &dm));
253:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
254:   PetscCallMPI(MPI_Comm_size(comm, &size));
255:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
256:   PetscCall(DMGetLocalVector(dm, &uloc));
257:   PetscCall(DMGlobalToLocalBegin(dm, u, INSERT_VALUES, uloc));
258:   PetscCall(DMGlobalToLocalEnd(dm, u, INSERT_VALUES, uloc));
259:   PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, uloc, 0.0, NULL, NULL, NULL));
260:   PetscCall(VecViewFromOptions(uloc, NULL, "-sol_view"));
261:   PetscCall(DMGetLocalSection(dm, &section));
262:   PetscCall(VecGetArrayRead(uloc, &array));
263:   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
264:   PetscCall(DMGetCoordinateSection(dm, &coordSection));
265:   PetscCall(VecGetArrayRead(coordinates, &coords));
266:   for (PetscInt p = 0; p < numPlanes; ++p) {
267:     DMLabel         label;
268:     char            name[PETSC_MAX_PATH_LEN];
269:     Mat             F;
270:     Vec             x, y;
271:     IS              stratum;
272:     PetscReal      *ray, *gray;
273:     PetscScalar    *rvals, *svals, *gsvals;
274:     PetscInt       *perm, *nperm;
275:     PetscInt        n, N, i, j, off, offu;
276:     PetscMPIInt     in;
277:     const PetscInt *points;

279:     PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%" PetscInt_FMT, p));
280:     PetscCall(DMGetLabel(dm, name, &label));
281:     PetscCall(DMLabelGetStratumIS(label, 1, &stratum));
282:     PetscCall(ISGetLocalSize(stratum, &n));
283:     PetscCall(PetscMPIIntCast(n, &in));
284:     PetscCall(ISGetIndices(stratum, &points));
285:     PetscCall(PetscMalloc2(n, &ray, n, &svals));
286:     for (i = 0; i < n; ++i) {
287:       PetscCall(PetscSectionGetOffset(coordSection, points[i], &off));
288:       PetscCall(PetscSectionGetOffset(section, points[i], &offu));
289:       ray[i]   = PetscRealPart(coords[off + ((planeDir[p] + 1) % 2)]);
290:       svals[i] = array[offu];
291:     }
292:     /* Gather the ray data to proc 0 */
293:     if (size > 1) {
294:       PetscMPIInt *cnt, *displs;

296:       PetscCall(PetscCalloc2(size, &cnt, size, &displs));
297:       PetscCallMPI(MPI_Gather(&n, 1, MPIU_INT, cnt, 1, MPIU_INT, 0, comm));
298:       for (PetscInt p = 1; p < size; ++p) displs[p] = displs[p - 1] + cnt[p - 1];
299:       N = displs[size - 1] + cnt[size - 1];
300:       PetscCall(PetscMalloc2(N, &gray, N, &gsvals));
301:       PetscCallMPI(MPI_Gatherv(ray, in, MPIU_REAL, gray, cnt, displs, MPIU_REAL, 0, comm));
302:       PetscCallMPI(MPI_Gatherv(svals, in, MPIU_SCALAR, gsvals, cnt, displs, MPIU_SCALAR, 0, comm));
303:       PetscCall(PetscFree2(cnt, displs));
304:     } else {
305:       N      = n;
306:       gray   = ray;
307:       gsvals = svals;
308:     }
309:     if (rank == 0) {
310:       /* Sort point along ray */
311:       PetscCall(PetscMalloc2(N, &perm, N, &nperm));
312:       for (i = 0; i < N; ++i) perm[i] = i;
313:       PetscCall(PetscSortRealWithPermutation(N, gray, perm));
314:       /* Count duplicates and squish mapping */
315:       nperm[0] = perm[0];
316:       for (i = 1, j = 1; i < N; ++i) {
317:         if (PetscAbsReal(gray[perm[i]] - gray[perm[i - 1]]) > PETSC_SMALL) nperm[j++] = perm[i];
318:       }
319:       /* Create FFT structs */
320:       PetscCall(MatCreateFFT(PETSC_COMM_SELF, 1, &j, MATFFTW, &F));
321:       PetscCall(MatCreateVecs(F, &x, &y));
322:       PetscCall(PetscObjectSetName((PetscObject)y, name));
323:       PetscCall(VecGetArray(x, &rvals));
324:       for (i = 0, j = 0; i < N; ++i) {
325:         if (i > 0 && PetscAbsReal(gray[perm[i]] - gray[perm[i - 1]]) < PETSC_SMALL) continue;
326:         rvals[j] = gsvals[nperm[j]];
327:         ++j;
328:       }
329:       PetscCall(PetscFree2(perm, nperm));
330:       if (size > 1) PetscCall(PetscFree2(gray, gsvals));
331:       PetscCall(VecRestoreArray(x, &rvals));
332:       /* Do FFT along the ray */
333:       PetscCall(MatMult(F, x, y));
334:       /* Chop FFT */
335:       PetscCall(VecFilter(y, PETSC_SMALL));
336:       PetscCall(VecViewFromOptions(x, NULL, "-real_view"));
337:       PetscCall(VecViewFromOptions(y, NULL, "-fft_view"));
338:       PetscCall(VecDestroy(&x));
339:       PetscCall(VecDestroy(&y));
340:       PetscCall(MatDestroy(&F));
341:     }
342:     PetscCall(ISRestoreIndices(stratum, &points));
343:     PetscCall(ISDestroy(&stratum));
344:     PetscCall(PetscFree2(ray, svals));
345:   }
346:   PetscCall(VecRestoreArrayRead(coordinates, &coords));
347:   PetscCall(VecRestoreArrayRead(uloc, &array));
348:   PetscCall(DMRestoreLocalVector(dm, &uloc));
349:   PetscFunctionReturn(PETSC_SUCCESS);
350: }

352: static PetscErrorCode ComputeAdjoint(Vec u, AppCtx *user)
353: {
354:   PetscFunctionBegin;
355:   if (!user->adjoint) PetscFunctionReturn(PETSC_SUCCESS);
356:   DM   dm, dmAdj;
357:   SNES snesAdj;
358:   Vec  uAdj;

360:   PetscCall(VecGetDM(u, &dm));
361:   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snesAdj));
362:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)snesAdj, "adjoint_"));
363:   PetscCall(DMClone(dm, &dmAdj));
364:   PetscCall(SNESSetDM(snesAdj, dmAdj));
365:   PetscCall(SetupDiscretization(dmAdj, "adjoint", SetupAdjointProblem, user));
366:   PetscCall(DMCreateGlobalVector(dmAdj, &uAdj));
367:   PetscCall(PetscObjectSetName((PetscObject)uAdj, "adjoint"));
368:   PetscCall(DMPlexSetSNESLocalFEM(dmAdj, PETSC_FALSE, &user));
369:   PetscCall(SNESSetFromOptions(snesAdj));
370:   PetscCall(SNESSolve(snesAdj, NULL, uAdj));
371:   PetscCall(SNESGetSolution(snesAdj, &uAdj));
372:   PetscCall(VecViewFromOptions(uAdj, NULL, "-adjoint_view"));
373:   /* Error representation */
374:   {
375:     DM        dmErr, dmErrAux, dms[2];
376:     Vec       errorEst, errorL2, uErr, uErrLoc, uAdjLoc, uAdjProj;
377:     IS       *subis;
378:     PetscReal errorEstTot, errorL2Norm, errorL2Tot;
379:     PetscInt  N;
380:     PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {user->homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u};
381:     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};
382:     PetscCtx ctxs[1] = {0};

384:     ctxs[0] = user;
385:     PetscCall(DMClone(dm, &dmErr));
386:     PetscCall(SetupDiscretization(dmErr, "error", SetupErrorProblem, user));
387:     PetscCall(DMGetGlobalVector(dmErr, &errorEst));
388:     PetscCall(DMGetGlobalVector(dmErr, &errorL2));
389:     /*   Compute auxiliary data (solution and projection of adjoint solution) */
390:     PetscCall(DMGetLocalVector(dmAdj, &uAdjLoc));
391:     PetscCall(DMGlobalToLocalBegin(dmAdj, uAdj, INSERT_VALUES, uAdjLoc));
392:     PetscCall(DMGlobalToLocalEnd(dmAdj, uAdj, INSERT_VALUES, uAdjLoc));
393:     PetscCall(DMGetGlobalVector(dm, &uAdjProj));
394:     PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, uAdjLoc));
395:     PetscCall(DMProjectField(dm, 0.0, u, identity, INSERT_VALUES, uAdjProj));
396:     PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, NULL));
397:     PetscCall(DMRestoreLocalVector(dmAdj, &uAdjLoc));
398:     /*   Attach auxiliary data */
399:     dms[0] = dm;
400:     dms[1] = dm;
401:     PetscCall(DMCreateSuperDM(dms, 2, &subis, &dmErrAux));
402:     if (0) {
403:       PetscSection sec;

405:       PetscCall(DMGetLocalSection(dms[0], &sec));
406:       PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD));
407:       PetscCall(DMGetLocalSection(dms[1], &sec));
408:       PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD));
409:       PetscCall(DMGetLocalSection(dmErrAux, &sec));
410:       PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD));
411:     }
412:     PetscCall(DMViewFromOptions(dmErrAux, NULL, "-dm_err_view"));
413:     PetscCall(ISViewFromOptions(subis[0], NULL, "-super_is_view"));
414:     PetscCall(ISViewFromOptions(subis[1], NULL, "-super_is_view"));
415:     PetscCall(DMGetGlobalVector(dmErrAux, &uErr));
416:     PetscCall(VecViewFromOptions(u, NULL, "-map_vec_view"));
417:     PetscCall(VecViewFromOptions(uAdjProj, NULL, "-map_vec_view"));
418:     PetscCall(VecViewFromOptions(uErr, NULL, "-map_vec_view"));
419:     PetscCall(VecISCopy(uErr, subis[0], SCATTER_FORWARD, u));
420:     PetscCall(VecISCopy(uErr, subis[1], SCATTER_FORWARD, uAdjProj));
421:     PetscCall(DMRestoreGlobalVector(dm, &uAdjProj));
422:     for (PetscInt i = 0; i < 2; ++i) PetscCall(ISDestroy(&subis[i]));
423:     PetscCall(PetscFree(subis));
424:     PetscCall(DMGetLocalVector(dmErrAux, &uErrLoc));
425:     PetscCall(DMGlobalToLocalBegin(dm, uErr, INSERT_VALUES, uErrLoc));
426:     PetscCall(DMGlobalToLocalEnd(dm, uErr, INSERT_VALUES, uErrLoc));
427:     PetscCall(DMRestoreGlobalVector(dmErrAux, &uErr));
428:     PetscCall(DMSetAuxiliaryVec(dmAdj, NULL, 0, 0, uErrLoc));
429:     /*   Compute cellwise error estimate */
430:     PetscCall(VecSet(errorEst, 0.0));
431:     PetscCall(DMPlexComputeCellwiseIntegralFEM(dmAdj, uAdj, errorEst, user));
432:     PetscCall(DMSetAuxiliaryVec(dmAdj, NULL, 0, 0, NULL));
433:     PetscCall(DMRestoreLocalVector(dmErrAux, &uErrLoc));
434:     PetscCall(DMDestroy(&dmErrAux));
435:     /*   Plot cellwise error vector */
436:     PetscCall(VecViewFromOptions(errorEst, NULL, "-error_view"));
437:     /*   Compute ratio of estimate (sum over cells) with actual L_2 error */
438:     PetscCall(DMComputeL2Diff(dm, 0.0, funcs, ctxs, u, &errorL2Norm));
439:     PetscCall(DMPlexComputeL2DiffVec(dm, 0.0, funcs, ctxs, u, errorL2));
440:     PetscCall(VecViewFromOptions(errorL2, NULL, "-l2_error_view"));
441:     PetscCall(VecNorm(errorL2, NORM_INFINITY, &errorL2Tot));
442:     PetscCall(VecNorm(errorEst, NORM_INFINITY, &errorEstTot));
443:     PetscCall(VecGetSize(errorEst, &N));
444:     PetscCall(VecPointwiseDivide(errorEst, errorEst, errorL2));
445:     PetscCall(PetscObjectSetName((PetscObject)errorEst, "Error ratio"));
446:     PetscCall(VecViewFromOptions(errorEst, NULL, "-error_ratio_view"));
447:     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))));
448:     PetscCall(DMRestoreGlobalVector(dmErr, &errorEst));
449:     PetscCall(DMRestoreGlobalVector(dmErr, &errorL2));
450:     PetscCall(DMDestroy(&dmErr));
451:   }
452:   PetscCall(DMDestroy(&dmAdj));
453:   PetscCall(VecDestroy(&uAdj));
454:   PetscCall(SNESDestroy(&snesAdj));
455:   PetscFunctionReturn(PETSC_SUCCESS);
456: }

458: static PetscErrorCode ErrorView(Vec u, AppCtx *user)
459: {
460:   PetscErrorCode (*sol)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
461:   void     *ctx;
462:   DM        dm;
463:   PetscDS   ds;
464:   PetscReal error;
465:   PetscInt  N;

467:   PetscFunctionBegin;
468:   if (!user->viewError) PetscFunctionReturn(PETSC_SUCCESS);
469:   PetscCall(VecGetDM(u, &dm));
470:   PetscCall(DMGetDS(dm, &ds));
471:   PetscCall(PetscDSGetExactSolution(ds, 0, &sol, &ctx));
472:   PetscCall(VecGetSize(u, &N));
473:   PetscCall(DMComputeL2Diff(dm, 0.0, &sol, &ctx, u, &error));
474:   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "N: %" PetscInt_FMT " L2 error: %g\n", N, (double)error));
475:   PetscFunctionReturn(PETSC_SUCCESS);
476: }

478: int main(int argc, char **argv)
479: {
480:   DM        dm;   /* Problem specification */
481:   SNES      snes; /* Nonlinear solver */
482:   Vec       u;    /* Solutions */
483:   AppCtx    user; /* User-defined work context */
484:   PetscInt  planeDir[2]   = {0, 1};
485:   PetscReal planeCoord[2] = {0., 1.};

487:   PetscFunctionBeginUser;
488:   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
489:   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
490:   /* Primal system */
491:   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
492:   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
493:   PetscCall(SNESSetDM(snes, dm));
494:   PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user));
495:   PetscCall(DMCreateGlobalVector(dm, &u));
496:   PetscCall(PetscObjectSetName((PetscObject)u, "potential"));
497:   PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user));
498:   PetscCall(SNESSetFromOptions(snes));
499:   PetscCall(SNESSolve(snes, NULL, u));
500:   PetscCall(SNESGetSolution(snes, &u));
501:   PetscCall(VecViewFromOptions(u, NULL, "-potential_view"));
502:   PetscCall(ErrorView(u, &user));
503:   PetscCall(ComputeSpectral(u, 2, planeDir, planeCoord, &user));
504:   PetscCall(ComputeAdjoint(u, &user));
505:   /* Cleanup */
506:   PetscCall(VecDestroy(&u));
507:   PetscCall(SNESDestroy(&snes));
508:   PetscCall(DMDestroy(&dm));
509:   PetscCall(PetscFinalize());
510:   return 0;
511: }

513: /*TEST

515:   test:
516:     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9
517:     suffix: 2d_p1_conv
518:     requires: triangle
519:     args: -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2
520:   test:
521:     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9
522:     suffix: 2d_p2_conv
523:     requires: triangle
524:     args: -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2
525:   test:
526:     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9
527:     suffix: 2d_p3_conv
528:     requires: triangle
529:     args: -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2
530:   test:
531:     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9
532:     suffix: 2d_q1_conv
533:     args: -dm_plex_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2
534:   test:
535:     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9
536:     suffix: 2d_q2_conv
537:     args: -dm_plex_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2
538:   test:
539:     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9
540:     suffix: 2d_q3_conv
541:     args: -dm_plex_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2
542:   test:
543:     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9
544:     suffix: 2d_q1_ceed_conv
545:     requires: libceed
546:     args: -dm_plex_use_ceed -dm_plex_simplex 0 -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_ceed_conv
550:     requires: libceed
551:     args: -dm_plex_use_ceed -dm_plex_simplex 0 -potential_petscspace_degree 2 -cdm_default_quadrature_order 2 \
552:           -snes_convergence_estimate -convest_num_refine 2
553:   test:
554:     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9
555:     suffix: 2d_q3_ceed_conv
556:     requires: libceed
557:     args: -dm_plex_use_ceed -dm_plex_simplex 0 -potential_petscspace_degree 3 -cdm_default_quadrature_order 3 \
558:           -snes_convergence_estimate -convest_num_refine 2
559:   test:
560:     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9
561:     suffix: 2d_q1_shear_conv
562:     args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2
563:   test:
564:     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9
565:     suffix: 2d_q2_shear_conv
566:     args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2
567:   test:
568:     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9
569:     suffix: 2d_q3_shear_conv
570:     args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2
571:   test:
572:     # Using -convest_num_refine 3 we get L_2 convergence rate: 1.7
573:     suffix: 3d_p1_conv
574:     requires: ctetgen
575:     args: -dm_plex_dim 3 -dm_refine 1 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1
576:   test:
577:     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 2.8
578:     suffix: 3d_p2_conv
579:     requires: ctetgen
580:     args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 1
581:   test:
582:     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 4.0
583:     suffix: 3d_p3_conv
584:     requires: ctetgen
585:     args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 1
586:   test:
587:     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.8
588:     suffix: 3d_q1_conv
589:     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_refine 1 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1
590:   test:
591:     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.8
592:     suffix: 3d_q2_conv
593:     args: -dm_plex_dim 3 -dm_plex_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 1
594:   test:
595:     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 3.8
596:     suffix: 3d_q3_conv
597:     args: -dm_plex_dim 3 -dm_plex_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 1
598:   test:
599:     suffix: 2d_p1_fas_full
600:     requires: triangle
601:     args: -potential_petscspace_degree 1 -dm_refine_hierarchy 5 \
602:       -snes_max_it 1 -snes_type fas -snes_fas_levels 5 -snes_fas_type full -snes_fas_full_total \
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:   test:
608:     suffix: 2d_p1_fas_full_homogeneous
609:     requires: triangle
610:     args: -homogeneous -potential_petscspace_degree 1 -dm_refine_hierarchy 5 \
611:       -snes_max_it 1 -snes_type fas -snes_fas_levels 5 -snes_fas_type full \
612:         -fas_coarse_snes_monitor -fas_coarse_snes_max_it 1 -fas_coarse_ksp_atol 1.e-13 \
613:         -fas_levels_snes_monitor -fas_levels_snes_max_it 1 -fas_levels_snes_type newtonls \
614:           -fas_levels_pc_type none -fas_levels_ksp_max_it 2 -fas_levels_ksp_converged_maxits -fas_levels_ksp_type chebyshev \
615:             -fas_levels_esteig_ksp_type cg -fas_levels_ksp_chebyshev_esteig 0,0.25,0,1.1 -fas_levels_esteig_ksp_max_it 10

617:   test:
618:     suffix: 2d_p1_scalable
619:     requires: triangle
620:     args: -potential_petscspace_degree 1 -dm_refine 3 \
621:       -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned \
622:       -pc_type gamg -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 \
623:         -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \
624:         -pc_gamg_coarse_eq_limit 1000 \
625:         -pc_gamg_threshold 0.05 \
626:         -pc_gamg_threshold_scale .0 \
627:         -mg_levels_ksp_type chebyshev \
628:         -mg_levels_ksp_max_it 1 \
629:         -mg_levels_pc_type jacobi \
630:       -matptap_via scalable
631:     output_file: output/empty.out
632:   test:
633:     suffix: 2d_p1_gmg_vcycle
634:     requires: triangle
635:     output_file: output/empty.out
636:     args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
637:           -ksp_rtol 5e-10 -pc_type mg \
638:             -mg_levels_ksp_max_it 1 \
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:   # Run with -dm_refine_hierarchy 3 to get a better idea of the solver
644:   testset:
645:     args: -potential_petscspace_degree 1 -dm_refine_hierarchy 2 \
646:           -ksp_rtol 5e-10 -pc_type mg -pc_mg_type full \
647:             -mg_levels_ksp_max_it 2 \
648:             -mg_levels_esteig_ksp_type cg \
649:             -mg_levels_esteig_ksp_max_it 10 \
650:             -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \
651:             -mg_levels_pc_type jacobi
652:     output_file: output/empty.out
653:     test:
654:       suffix: 2d_p1_gmg_fcycle
655:       requires: triangle
656:       args: -dm_plex_box_faces 2,2
657:     test:
658:       suffix: 2d_q1_gmg_fcycle
659:       args: -dm_plex_simplex 0 -dm_plex_box_faces 2,2
660:     test:
661:       suffix: 3d_p1_gmg_fcycle
662:       requires: ctetgen
663:       args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,1
664:     test:
665:       suffix: 3d_q1_gmg_fcycle
666:       args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 2,2,1
667:   test:
668:     suffix: 2d_p1_gmg_vcycle_adapt
669:     requires: triangle
670:     output_file: output/empty.out
671:     args: -petscpartitioner_type simple -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
672:           -ksp_rtol 5e-10 -pc_type mg -pc_mg_galerkin -pc_mg_adapt_interp_coarse_space harmonic -pc_mg_adapt_interp_n 8 \
673:             -mg_levels_ksp_max_it 1 \
674:             -mg_levels_esteig_ksp_type cg \
675:             -mg_levels_esteig_ksp_max_it 10 \
676:             -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \
677:             -mg_levels_pc_type jacobi
678:   test:
679:     suffix: 2d_p1_spectral_0
680:     requires: triangle fftw !complex
681:     args: -dm_plex_box_faces 1,1 -potential_petscspace_degree 1 -dm_refine 6 -spectral -fft_view
682:   test:
683:     suffix: 2d_p1_spectral_1
684:     requires: triangle fftw !complex
685:     nsize: 2
686:     args: -dm_plex_box_faces 4,4 -potential_petscspace_degree 1 -spectral -fft_view
687:   test:
688:     suffix: 2d_p1_adj_0
689:     requires: triangle
690:     args: -potential_petscspace_degree 1 -dm_refine 1 -adjoint -adjoint_petscspace_degree 1 -error_petscspace_degree 0
691:   test:
692:     nsize: 2
693:     requires: kokkos_kernels
694:     suffix: kokkos
695:     args: -dm_plex_dim 3 -dm_plex_box_faces 2,3,6 -petscpartitioner_type simple -dm_plex_simplex 0 -potential_petscspace_degree 1 \
696:          -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 \
697:          -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 \
698:          -ksp_monitor -snes_monitor -dm_view -dm_mat_type aijkokkos -dm_vec_type kokkos

700: TEST*/