xref: /petsc/src/snes/tutorials/ex13.c (revision c78f73d15fb0051bc87be86ac52a7c5747d72d5c)
1c4762a1bSJed Brown static char help[] = "Poisson Problem in 2d and 3d with finite elements.\n\
2c4762a1bSJed Brown We solve the Poisson problem in a rectangular\n\
3c4762a1bSJed Brown domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\
4c4762a1bSJed Brown This example supports automatic convergence estimation\n\
5c4762a1bSJed Brown and eventually adaptivity.\n\n\n";
6c4762a1bSJed Brown 
7c4762a1bSJed Brown #include <petscdmplex.h>
8c4762a1bSJed Brown #include <petscsnes.h>
9c4762a1bSJed Brown #include <petscds.h>
10c4762a1bSJed Brown #include <petscconvest.h>
11c4762a1bSJed Brown 
12c4762a1bSJed Brown typedef struct {
13c4762a1bSJed Brown   /* Domain and mesh definition */
14c4762a1bSJed Brown   PetscBool spectral;    /* Look at the spectrum along planes in the solution */
15c4762a1bSJed Brown   PetscBool shear;       /* Shear the domain */
16c4762a1bSJed Brown   PetscBool adjoint;     /* Solve the adjoint problem */
17ad1a7c93SMatthew Knepley   PetscBool homogeneous; /* Use homogeneous boudnary conditions */
18ad1a7c93SMatthew Knepley   PetscBool viewError;   /* Output the solution error */
19c4762a1bSJed Brown } AppCtx;
20c4762a1bSJed Brown 
21d71ae5a4SJacob Faibussowitsch static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
22d71ae5a4SJacob Faibussowitsch {
23c4762a1bSJed Brown   *u = 0.0;
24c4762a1bSJed Brown   return 0;
25c4762a1bSJed Brown }
26c4762a1bSJed Brown 
27d71ae5a4SJacob Faibussowitsch static PetscErrorCode trig_inhomogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
28d71ae5a4SJacob Faibussowitsch {
29c4762a1bSJed Brown   PetscInt d;
30c4762a1bSJed Brown   *u = 0.0;
31c4762a1bSJed Brown   for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0 * PETSC_PI * x[d]);
32c4762a1bSJed Brown   return 0;
33c4762a1bSJed Brown }
34c4762a1bSJed Brown 
35d71ae5a4SJacob Faibussowitsch static PetscErrorCode trig_homogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
36d71ae5a4SJacob Faibussowitsch {
37b724ec41SToby Isaac   PetscInt d;
38b724ec41SToby Isaac   *u = 1.0;
39b724ec41SToby Isaac   for (d = 0; d < dim; ++d) *u *= PetscSinReal(2.0 * PETSC_PI * x[d]);
40b724ec41SToby Isaac   return 0;
41b724ec41SToby Isaac }
42b724ec41SToby Isaac 
43c4762a1bSJed Brown /* Compute integral of (residual of solution)*(adjoint solution - projection of adjoint solution) */
44d71ae5a4SJacob Faibussowitsch 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[])
45d71ae5a4SJacob Faibussowitsch {
46c4762a1bSJed Brown   obj[0] = a[aOff[0]] * (u[0] - a[aOff[1]]);
47c4762a1bSJed Brown }
48c4762a1bSJed Brown 
49d71ae5a4SJacob Faibussowitsch 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[])
50d71ae5a4SJacob Faibussowitsch {
51c4762a1bSJed Brown   PetscInt d;
52c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f0[0] += -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]);
53c4762a1bSJed Brown }
54c4762a1bSJed Brown 
55d71ae5a4SJacob Faibussowitsch 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[])
56d71ae5a4SJacob Faibussowitsch {
57b724ec41SToby Isaac   PetscInt d;
58b724ec41SToby Isaac   for (d = 0; d < dim; ++d) {
59b724ec41SToby Isaac     PetscScalar v = 1.;
60b724ec41SToby Isaac     for (PetscInt e = 0; e < dim; e++) {
61b724ec41SToby Isaac       if (e == d) {
62b724ec41SToby Isaac         v *= -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]);
63b724ec41SToby Isaac       } else {
64b724ec41SToby Isaac         v *= PetscSinReal(2.0 * PETSC_PI * x[d]);
65b724ec41SToby Isaac       }
66b724ec41SToby Isaac     }
67b724ec41SToby Isaac     f0[0] += v;
68b724ec41SToby Isaac   }
69b724ec41SToby Isaac }
70b724ec41SToby Isaac 
71d71ae5a4SJacob Faibussowitsch 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[])
72d71ae5a4SJacob Faibussowitsch {
73c4762a1bSJed Brown   f0[0] = 1.0;
74c4762a1bSJed Brown }
75c4762a1bSJed Brown 
76d71ae5a4SJacob Faibussowitsch 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[])
77d71ae5a4SJacob Faibussowitsch {
78c4762a1bSJed Brown   f0[0] = a[0];
79c4762a1bSJed Brown }
80c4762a1bSJed Brown 
81d71ae5a4SJacob Faibussowitsch 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[])
82d71ae5a4SJacob Faibussowitsch {
83c4762a1bSJed Brown   PetscInt d;
84c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
85c4762a1bSJed Brown }
86c4762a1bSJed Brown 
87d71ae5a4SJacob Faibussowitsch 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[])
88d71ae5a4SJacob Faibussowitsch {
89c4762a1bSJed Brown   PetscInt d;
90c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
91c4762a1bSJed Brown }
92c4762a1bSJed Brown 
93d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
94d71ae5a4SJacob Faibussowitsch {
95c4762a1bSJed Brown   PetscFunctionBeginUser;
96c4762a1bSJed Brown   options->shear       = PETSC_FALSE;
97c4762a1bSJed Brown   options->spectral    = PETSC_FALSE;
98c4762a1bSJed Brown   options->adjoint     = PETSC_FALSE;
99b724ec41SToby Isaac   options->homogeneous = PETSC_FALSE;
100ad1a7c93SMatthew Knepley   options->viewError   = PETSC_FALSE;
101c4762a1bSJed Brown 
102d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");
1039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-shear", "Shear the domain", "ex13.c", options->shear, &options->shear, NULL));
1049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-spectral", "Look at the spectrum along planes of the solution", "ex13.c", options->spectral, &options->spectral, NULL));
1059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-adjoint", "Solve the adjoint problem", "ex13.c", options->adjoint, &options->adjoint, NULL));
1069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-homogeneous", "Use homogeneous boundary conditions", "ex13.c", options->homogeneous, &options->homogeneous, NULL));
1079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-error_view", "Output the solution error", "ex13.c", options->viewError, &options->viewError, NULL));
108d0609cedSBarry Smith   PetscOptionsEnd();
109c4762a1bSJed Brown   PetscFunctionReturn(0);
110c4762a1bSJed Brown }
111c4762a1bSJed Brown 
112d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateSpectralPlanes(DM dm, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user)
113d71ae5a4SJacob Faibussowitsch {
114c4762a1bSJed Brown   PetscSection       coordSection;
115c4762a1bSJed Brown   Vec                coordinates;
116c4762a1bSJed Brown   const PetscScalar *coords;
117c4762a1bSJed Brown   PetscInt           dim, p, vStart, vEnd, v;
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   PetscFunctionBeginUser;
1209566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &dim));
1219566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
1229566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
1239566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(dm, &coordSection));
1249566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(coordinates, &coords));
125c4762a1bSJed Brown   for (p = 0; p < numPlanes; ++p) {
126c4762a1bSJed Brown     DMLabel label;
127c4762a1bSJed Brown     char    name[PETSC_MAX_PATH_LEN];
128c4762a1bSJed Brown 
12963a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%" PetscInt_FMT, p));
1309566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(dm, name));
1319566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dm, name, &label));
1329566063dSJacob Faibussowitsch     PetscCall(DMLabelAddStratum(label, 1));
133c4762a1bSJed Brown     for (v = vStart; v < vEnd; ++v) {
134c4762a1bSJed Brown       PetscInt off;
135c4762a1bSJed Brown 
1369566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(coordSection, v, &off));
13748a46eb9SPierre Jolivet       if (PetscAbsReal(planeCoord[p] - PetscRealPart(coords[off + planeDir[p]])) < PETSC_SMALL) PetscCall(DMLabelSetValue(label, v, 1));
138c4762a1bSJed Brown     }
139c4762a1bSJed Brown   }
1409566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(coordinates, &coords));
141c4762a1bSJed Brown   PetscFunctionReturn(0);
142c4762a1bSJed Brown }
143c4762a1bSJed Brown 
144d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
145d71ae5a4SJacob Faibussowitsch {
146c4762a1bSJed Brown   PetscFunctionBeginUser;
1479566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
1489566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
1499566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
1509566063dSJacob Faibussowitsch   if (user->shear) PetscCall(DMPlexShearGeometry(*dm, DM_X, NULL));
1519566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(*dm, user));
1529566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
153c4762a1bSJed Brown   if (user->spectral) {
154c4762a1bSJed Brown     PetscInt  planeDir[2]   = {0, 1};
155c4762a1bSJed Brown     PetscReal planeCoord[2] = {0., 1.};
156c4762a1bSJed Brown 
1579566063dSJacob Faibussowitsch     PetscCall(CreateSpectralPlanes(*dm, 2, planeDir, planeCoord, user));
158c4762a1bSJed Brown   }
159c4762a1bSJed Brown   PetscFunctionReturn(0);
160c4762a1bSJed Brown }
161c4762a1bSJed Brown 
162d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
163d71ae5a4SJacob Faibussowitsch {
16445480ffeSMatthew G. Knepley   PetscDS        ds;
16545480ffeSMatthew G. Knepley   DMLabel        label;
166c4762a1bSJed Brown   const PetscInt id                                                                             = 1;
167b724ec41SToby Isaac   PetscPointFunc f0                                                                             = user->homogeneous ? f0_trig_homogeneous_u : f0_trig_inhomogeneous_u;
168b724ec41SToby Isaac   PetscErrorCode (*ex)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = user->homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u;
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   PetscFunctionBeginUser;
1719566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
1729566063dSJacob Faibussowitsch   PetscCall(PetscDSSetResidual(ds, 0, f0, f1_u));
1739566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu));
1749566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(ds, 0, ex, user));
1759566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
176*c78f73d1SMatthew G. Knepley   if (label) PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))ex, NULL, user, NULL));
177c4762a1bSJed Brown   PetscFunctionReturn(0);
178c4762a1bSJed Brown }
179c4762a1bSJed Brown 
180d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupAdjointProblem(DM dm, AppCtx *user)
181d71ae5a4SJacob Faibussowitsch {
18245480ffeSMatthew G. Knepley   PetscDS        ds;
18345480ffeSMatthew G. Knepley   DMLabel        label;
184c4762a1bSJed Brown   const PetscInt id = 1;
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   PetscFunctionBeginUser;
1879566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
1889566063dSJacob Faibussowitsch   PetscCall(PetscDSSetResidual(ds, 0, f0_unity_u, f1_u));
1899566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu));
1909566063dSJacob Faibussowitsch   PetscCall(PetscDSSetObjective(ds, 0, obj_error_u));
1919566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
1929566063dSJacob Faibussowitsch   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))zero, NULL, user, NULL));
193c4762a1bSJed Brown   PetscFunctionReturn(0);
194c4762a1bSJed Brown }
195c4762a1bSJed Brown 
196d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupErrorProblem(DM dm, AppCtx *user)
197d71ae5a4SJacob Faibussowitsch {
198c4762a1bSJed Brown   PetscDS prob;
199c4762a1bSJed Brown 
200c4762a1bSJed Brown   PetscFunctionBeginUser;
2019566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &prob));
202c4762a1bSJed Brown   PetscFunctionReturn(0);
203c4762a1bSJed Brown }
204c4762a1bSJed Brown 
205d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
206d71ae5a4SJacob Faibussowitsch {
207c4762a1bSJed Brown   DM             cdm = dm;
208c4762a1bSJed Brown   PetscFE        fe;
2095be41b18SMatthew G. Knepley   DMPolytopeType ct;
2105be41b18SMatthew G. Knepley   PetscBool      simplex;
2115be41b18SMatthew G. Knepley   PetscInt       dim, cStart;
212c4762a1bSJed Brown   char           prefix[PETSC_MAX_PATH_LEN];
213c4762a1bSJed Brown 
214c4762a1bSJed Brown   PetscFunctionBeginUser;
2159566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
2169566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
2179566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
2185be41b18SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
219c4762a1bSJed Brown   /* Create finite element */
2209566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
2219566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe));
2229566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe, name));
223c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
2249566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
2259566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
2269566063dSJacob Faibussowitsch   PetscCall((*setup)(dm, user));
227c4762a1bSJed Brown   while (cdm) {
2289566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm, cdm));
229c4762a1bSJed Brown     /* TODO: Check whether the boundary of coarse meshes is marked */
2309566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
231c4762a1bSJed Brown   }
2329566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
233c4762a1bSJed Brown   PetscFunctionReturn(0);
234c4762a1bSJed Brown }
235c4762a1bSJed Brown 
236d71ae5a4SJacob Faibussowitsch static PetscErrorCode ComputeSpectral(DM dm, Vec u, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user)
237d71ae5a4SJacob Faibussowitsch {
238c4762a1bSJed Brown   MPI_Comm           comm;
239c4762a1bSJed Brown   PetscSection       coordSection, section;
240c4762a1bSJed Brown   Vec                coordinates, uloc;
241c4762a1bSJed Brown   const PetscScalar *coords, *array;
242c4762a1bSJed Brown   PetscInt           p;
243c4762a1bSJed Brown   PetscMPIInt        size, rank;
244c4762a1bSJed Brown 
245c4762a1bSJed Brown   PetscFunctionBeginUser;
2469566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2479566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
2489566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2499566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &uloc));
2509566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, u, INSERT_VALUES, uloc));
2519566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, u, INSERT_VALUES, uloc));
2529566063dSJacob Faibussowitsch   PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, uloc, 0.0, NULL, NULL, NULL));
2539566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(uloc, NULL, "-sol_view"));
2549566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
2559566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(uloc, &array));
2569566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
2579566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(dm, &coordSection));
2589566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(coordinates, &coords));
259c4762a1bSJed Brown   for (p = 0; p < numPlanes; ++p) {
260c4762a1bSJed Brown     DMLabel         label;
261c4762a1bSJed Brown     char            name[PETSC_MAX_PATH_LEN];
262c4762a1bSJed Brown     Mat             F;
263c4762a1bSJed Brown     Vec             x, y;
264c4762a1bSJed Brown     IS              stratum;
265c4762a1bSJed Brown     PetscReal      *ray, *gray;
266c4762a1bSJed Brown     PetscScalar    *rvals, *svals, *gsvals;
267c4762a1bSJed Brown     PetscInt       *perm, *nperm;
268c4762a1bSJed Brown     PetscInt        n, N, i, j, off, offu;
269c4762a1bSJed Brown     const PetscInt *points;
270c4762a1bSJed Brown 
27163a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%" PetscInt_FMT, p));
2729566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dm, name, &label));
2739566063dSJacob Faibussowitsch     PetscCall(DMLabelGetStratumIS(label, 1, &stratum));
2749566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(stratum, &n));
2759566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(stratum, &points));
2769566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(n, &ray, n, &svals));
277c4762a1bSJed Brown     for (i = 0; i < n; ++i) {
2789566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(coordSection, points[i], &off));
2799566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(section, points[i], &offu));
280c4762a1bSJed Brown       ray[i]   = PetscRealPart(coords[off + ((planeDir[p] + 1) % 2)]);
281c4762a1bSJed Brown       svals[i] = array[offu];
282c4762a1bSJed Brown     }
283c4762a1bSJed Brown     /* Gather the ray data to proc 0 */
284c4762a1bSJed Brown     if (size > 1) {
285c4762a1bSJed Brown       PetscMPIInt *cnt, *displs, p;
286c4762a1bSJed Brown 
2879566063dSJacob Faibussowitsch       PetscCall(PetscCalloc2(size, &cnt, size, &displs));
2889566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Gather(&n, 1, MPIU_INT, cnt, 1, MPIU_INT, 0, comm));
289c4762a1bSJed Brown       for (p = 1; p < size; ++p) displs[p] = displs[p - 1] + cnt[p - 1];
290c4762a1bSJed Brown       N = displs[size - 1] + cnt[size - 1];
2919566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(N, &gray, N, &gsvals));
2929566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Gatherv(ray, n, MPIU_REAL, gray, cnt, displs, MPIU_REAL, 0, comm));
2939566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Gatherv(svals, n, MPIU_SCALAR, gsvals, cnt, displs, MPIU_SCALAR, 0, comm));
2949566063dSJacob Faibussowitsch       PetscCall(PetscFree2(cnt, displs));
295c4762a1bSJed Brown     } else {
296c4762a1bSJed Brown       N      = n;
297c4762a1bSJed Brown       gray   = ray;
298c4762a1bSJed Brown       gsvals = svals;
299c4762a1bSJed Brown     }
300dd400576SPatrick Sanan     if (rank == 0) {
301c4762a1bSJed Brown       /* Sort point along ray */
3029566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(N, &perm, N, &nperm));
303ad540459SPierre Jolivet       for (i = 0; i < N; ++i) perm[i] = i;
3049566063dSJacob Faibussowitsch       PetscCall(PetscSortRealWithPermutation(N, gray, perm));
305c4762a1bSJed Brown       /* Count duplicates and squish mapping */
306c4762a1bSJed Brown       nperm[0] = perm[0];
307c4762a1bSJed Brown       for (i = 1, j = 1; i < N; ++i) {
308c4762a1bSJed Brown         if (PetscAbsReal(gray[perm[i]] - gray[perm[i - 1]]) > PETSC_SMALL) nperm[j++] = perm[i];
309c4762a1bSJed Brown       }
310c4762a1bSJed Brown       /* Create FFT structs */
3119566063dSJacob Faibussowitsch       PetscCall(MatCreateFFT(PETSC_COMM_SELF, 1, &j, MATFFTW, &F));
3129566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(F, &x, &y));
3139566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)y, name));
3149566063dSJacob Faibussowitsch       PetscCall(VecGetArray(x, &rvals));
315c4762a1bSJed Brown       for (i = 0, j = 0; i < N; ++i) {
316c4762a1bSJed Brown         if (i > 0 && PetscAbsReal(gray[perm[i]] - gray[perm[i - 1]]) < PETSC_SMALL) continue;
317c4762a1bSJed Brown         rvals[j] = gsvals[nperm[j]];
318c4762a1bSJed Brown         ++j;
319c4762a1bSJed Brown       }
3209566063dSJacob Faibussowitsch       PetscCall(PetscFree2(perm, nperm));
3219566063dSJacob Faibussowitsch       if (size > 1) PetscCall(PetscFree2(gray, gsvals));
3229566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(x, &rvals));
323c4762a1bSJed Brown       /* Do FFT along the ray */
3249566063dSJacob Faibussowitsch       PetscCall(MatMult(F, x, y));
325c4762a1bSJed Brown       /* Chop FFT */
3269566063dSJacob Faibussowitsch       PetscCall(VecChop(y, PETSC_SMALL));
3279566063dSJacob Faibussowitsch       PetscCall(VecViewFromOptions(x, NULL, "-real_view"));
3289566063dSJacob Faibussowitsch       PetscCall(VecViewFromOptions(y, NULL, "-fft_view"));
3299566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&x));
3309566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&y));
3319566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&F));
332c4762a1bSJed Brown     }
3339566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(stratum, &points));
3349566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&stratum));
3359566063dSJacob Faibussowitsch     PetscCall(PetscFree2(ray, svals));
336c4762a1bSJed Brown   }
3379566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(coordinates, &coords));
3389566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(uloc, &array));
3399566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &uloc));
340c4762a1bSJed Brown   PetscFunctionReturn(0);
341c4762a1bSJed Brown }
342c4762a1bSJed Brown 
343d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
344d71ae5a4SJacob Faibussowitsch {
345c4762a1bSJed Brown   DM     dm;   /* Problem specification */
346c4762a1bSJed Brown   SNES   snes; /* Nonlinear solver */
347c4762a1bSJed Brown   Vec    u;    /* Solutions */
348c4762a1bSJed Brown   AppCtx user; /* User-defined work context */
349c4762a1bSJed Brown 
350327415f7SBarry Smith   PetscFunctionBeginUser;
3519566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3529566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
353c4762a1bSJed Brown   /* Primal system */
3549566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
3559566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
3569566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, dm));
3579566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user));
3589566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &u));
3599566063dSJacob Faibussowitsch   PetscCall(VecSet(u, 0.0));
3609566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)u, "potential"));
3619566063dSJacob Faibussowitsch   PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user));
3629566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
3639566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, u));
3649566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes, &u));
3659566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-potential_view"));
366ad1a7c93SMatthew Knepley   if (user.viewError) {
367ad1a7c93SMatthew Knepley     PetscErrorCode (*sol)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
368ad1a7c93SMatthew Knepley     void     *ctx;
369ad1a7c93SMatthew Knepley     PetscDS   ds;
370ad1a7c93SMatthew Knepley     PetscReal error;
371ad1a7c93SMatthew Knepley     PetscInt  N;
372ad1a7c93SMatthew Knepley 
3739566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dm, &ds));
3749566063dSJacob Faibussowitsch     PetscCall(PetscDSGetExactSolution(ds, 0, &sol, &ctx));
3759566063dSJacob Faibussowitsch     PetscCall(VecGetSize(u, &N));
3769566063dSJacob Faibussowitsch     PetscCall(DMComputeL2Diff(dm, 0.0, &sol, &ctx, u, &error));
37763a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "N: %" PetscInt_FMT " L2 error: %g\n", N, (double)error));
378ad1a7c93SMatthew Knepley   }
379c4762a1bSJed Brown   if (user.spectral) {
380c4762a1bSJed Brown     PetscInt  planeDir[2]   = {0, 1};
381c4762a1bSJed Brown     PetscReal planeCoord[2] = {0., 1.};
382c4762a1bSJed Brown 
3839566063dSJacob Faibussowitsch     PetscCall(ComputeSpectral(dm, u, 2, planeDir, planeCoord, &user));
384c4762a1bSJed Brown   }
385c4762a1bSJed Brown   /* Adjoint system */
386c4762a1bSJed Brown   if (user.adjoint) {
387c4762a1bSJed Brown     DM   dmAdj;
388c4762a1bSJed Brown     SNES snesAdj;
389c4762a1bSJed Brown     Vec  uAdj;
390c4762a1bSJed Brown 
3919566063dSJacob Faibussowitsch     PetscCall(SNESCreate(PETSC_COMM_WORLD, &snesAdj));
3929566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)snesAdj, "adjoint_"));
3939566063dSJacob Faibussowitsch     PetscCall(DMClone(dm, &dmAdj));
3949566063dSJacob Faibussowitsch     PetscCall(SNESSetDM(snesAdj, dmAdj));
3959566063dSJacob Faibussowitsch     PetscCall(SetupDiscretization(dmAdj, "adjoint", SetupAdjointProblem, &user));
3969566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(dmAdj, &uAdj));
3979566063dSJacob Faibussowitsch     PetscCall(VecSet(uAdj, 0.0));
3989566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)uAdj, "adjoint"));
3999566063dSJacob Faibussowitsch     PetscCall(DMPlexSetSNESLocalFEM(dmAdj, &user, &user, &user));
4009566063dSJacob Faibussowitsch     PetscCall(SNESSetFromOptions(snesAdj));
4019566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snesAdj, NULL, uAdj));
4029566063dSJacob Faibussowitsch     PetscCall(SNESGetSolution(snesAdj, &uAdj));
4039566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(uAdj, NULL, "-adjoint_view"));
404c4762a1bSJed Brown     /* Error representation */
405c4762a1bSJed Brown     {
406c4762a1bSJed Brown       DM        dmErr, dmErrAux, dms[2];
407c4762a1bSJed Brown       Vec       errorEst, errorL2, uErr, uErrLoc, uAdjLoc, uAdjProj;
408c4762a1bSJed Brown       IS       *subis;
409c4762a1bSJed Brown       PetscReal errorEstTot, errorL2Norm, errorL2Tot;
410c4762a1bSJed Brown       PetscInt  N, i;
411b724ec41SToby Isaac       PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {user.homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u};
4129371c9d4SSatish Balay       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};
413c4762a1bSJed Brown       void *ctxs[1] = {0};
414c4762a1bSJed Brown 
415c4762a1bSJed Brown       ctxs[0] = &user;
4169566063dSJacob Faibussowitsch       PetscCall(DMClone(dm, &dmErr));
4179566063dSJacob Faibussowitsch       PetscCall(SetupDiscretization(dmErr, "error", SetupErrorProblem, &user));
4189566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(dmErr, &errorEst));
4199566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(dmErr, &errorL2));
420c4762a1bSJed Brown       /*   Compute auxiliary data (solution and projection of adjoint solution) */
4219566063dSJacob Faibussowitsch       PetscCall(DMGetLocalVector(dmAdj, &uAdjLoc));
4229566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalBegin(dmAdj, uAdj, INSERT_VALUES, uAdjLoc));
4239566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalEnd(dmAdj, uAdj, INSERT_VALUES, uAdjLoc));
4249566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(dm, &uAdjProj));
4259566063dSJacob Faibussowitsch       PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, uAdjLoc));
4269566063dSJacob Faibussowitsch       PetscCall(DMProjectField(dm, 0.0, u, identity, INSERT_VALUES, uAdjProj));
4279566063dSJacob Faibussowitsch       PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, NULL));
4289566063dSJacob Faibussowitsch       PetscCall(DMRestoreLocalVector(dmAdj, &uAdjLoc));
429c4762a1bSJed Brown       /*   Attach auxiliary data */
4309371c9d4SSatish Balay       dms[0] = dm;
4319371c9d4SSatish Balay       dms[1] = dm;
4329566063dSJacob Faibussowitsch       PetscCall(DMCreateSuperDM(dms, 2, &subis, &dmErrAux));
433c4762a1bSJed Brown       if (0) {
434c4762a1bSJed Brown         PetscSection sec;
435c4762a1bSJed Brown 
4369566063dSJacob Faibussowitsch         PetscCall(DMGetLocalSection(dms[0], &sec));
4379566063dSJacob Faibussowitsch         PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD));
4389566063dSJacob Faibussowitsch         PetscCall(DMGetLocalSection(dms[1], &sec));
4399566063dSJacob Faibussowitsch         PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD));
4409566063dSJacob Faibussowitsch         PetscCall(DMGetLocalSection(dmErrAux, &sec));
4419566063dSJacob Faibussowitsch         PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD));
442c4762a1bSJed Brown       }
4439566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(dmErrAux, NULL, "-dm_err_view"));
4449566063dSJacob Faibussowitsch       PetscCall(ISViewFromOptions(subis[0], NULL, "-super_is_view"));
4459566063dSJacob Faibussowitsch       PetscCall(ISViewFromOptions(subis[1], NULL, "-super_is_view"));
4469566063dSJacob Faibussowitsch       PetscCall(DMGetGlobalVector(dmErrAux, &uErr));
4479566063dSJacob Faibussowitsch       PetscCall(VecViewFromOptions(u, NULL, "-map_vec_view"));
4489566063dSJacob Faibussowitsch       PetscCall(VecViewFromOptions(uAdjProj, NULL, "-map_vec_view"));
4499566063dSJacob Faibussowitsch       PetscCall(VecViewFromOptions(uErr, NULL, "-map_vec_view"));
4509566063dSJacob Faibussowitsch       PetscCall(VecISCopy(uErr, subis[0], SCATTER_FORWARD, u));
4519566063dSJacob Faibussowitsch       PetscCall(VecISCopy(uErr, subis[1], SCATTER_FORWARD, uAdjProj));
4529566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(dm, &uAdjProj));
4539566063dSJacob Faibussowitsch       for (i = 0; i < 2; ++i) PetscCall(ISDestroy(&subis[i]));
4549566063dSJacob Faibussowitsch       PetscCall(PetscFree(subis));
4559566063dSJacob Faibussowitsch       PetscCall(DMGetLocalVector(dmErrAux, &uErrLoc));
4569566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalBegin(dm, uErr, INSERT_VALUES, uErrLoc));
4579566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalEnd(dm, uErr, INSERT_VALUES, uErrLoc));
4589566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(dmErrAux, &uErr));
4599566063dSJacob Faibussowitsch       PetscCall(DMSetAuxiliaryVec(dmAdj, NULL, 0, 0, uErrLoc));
460c4762a1bSJed Brown       /*   Compute cellwise error estimate */
4619566063dSJacob Faibussowitsch       PetscCall(VecSet(errorEst, 0.0));
4629566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellwiseIntegralFEM(dmAdj, uAdj, errorEst, &user));
4639566063dSJacob Faibussowitsch       PetscCall(DMSetAuxiliaryVec(dmAdj, NULL, 0, 0, NULL));
4649566063dSJacob Faibussowitsch       PetscCall(DMRestoreLocalVector(dmErrAux, &uErrLoc));
4659566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&dmErrAux));
466c4762a1bSJed Brown       /*   Plot cellwise error vector */
4679566063dSJacob Faibussowitsch       PetscCall(VecViewFromOptions(errorEst, NULL, "-error_view"));
468c4762a1bSJed Brown       /*   Compute ratio of estimate (sum over cells) with actual L_2 error */
4699566063dSJacob Faibussowitsch       PetscCall(DMComputeL2Diff(dm, 0.0, funcs, ctxs, u, &errorL2Norm));
4709566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeL2DiffVec(dm, 0.0, funcs, ctxs, u, errorL2));
4719566063dSJacob Faibussowitsch       PetscCall(VecViewFromOptions(errorL2, NULL, "-l2_error_view"));
4729566063dSJacob Faibussowitsch       PetscCall(VecNorm(errorL2, NORM_INFINITY, &errorL2Tot));
4739566063dSJacob Faibussowitsch       PetscCall(VecNorm(errorEst, NORM_INFINITY, &errorEstTot));
4749566063dSJacob Faibussowitsch       PetscCall(VecGetSize(errorEst, &N));
4759566063dSJacob Faibussowitsch       PetscCall(VecPointwiseDivide(errorEst, errorEst, errorL2));
4769566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)errorEst, "Error ratio"));
4779566063dSJacob Faibussowitsch       PetscCall(VecViewFromOptions(errorEst, NULL, "-error_ratio_view"));
47863a3b9bcSJacob Faibussowitsch       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))));
4799566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(dmErr, &errorEst));
4809566063dSJacob Faibussowitsch       PetscCall(DMRestoreGlobalVector(dmErr, &errorL2));
4819566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&dmErr));
482c4762a1bSJed Brown     }
4839566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dmAdj));
4849566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&uAdj));
4859566063dSJacob Faibussowitsch     PetscCall(SNESDestroy(&snesAdj));
486c4762a1bSJed Brown   }
487c4762a1bSJed Brown   /* Cleanup */
4889566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
4899566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
4909566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
4919566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
492b122ec5aSJacob Faibussowitsch   return 0;
493c4762a1bSJed Brown }
494c4762a1bSJed Brown 
495c4762a1bSJed Brown /*TEST
496c4762a1bSJed Brown 
497c4762a1bSJed Brown   test:
4985be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9
4995be41b18SMatthew G. Knepley     suffix: 2d_p1_conv
500c4762a1bSJed Brown     requires: triangle
5015be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2
5025be41b18SMatthew G. Knepley   test:
5035be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9
5045be41b18SMatthew G. Knepley     suffix: 2d_p2_conv
5055be41b18SMatthew G. Knepley     requires: triangle
5065be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2
5075be41b18SMatthew G. Knepley   test:
5085be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9
5095be41b18SMatthew G. Knepley     suffix: 2d_p3_conv
5105be41b18SMatthew G. Knepley     requires: triangle
5115be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2
5125be41b18SMatthew G. Knepley   test:
5135be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9
5145be41b18SMatthew G. Knepley     suffix: 2d_q1_conv
51530602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2
5165be41b18SMatthew G. Knepley   test:
5175be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9
5185be41b18SMatthew G. Knepley     suffix: 2d_q2_conv
51930602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2
5205be41b18SMatthew G. Knepley   test:
5215be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9
5225be41b18SMatthew G. Knepley     suffix: 2d_q3_conv
52330602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2
5245be41b18SMatthew G. Knepley   test:
5255be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9
5265be41b18SMatthew G. Knepley     suffix: 2d_q1_shear_conv
52730602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2
5285be41b18SMatthew G. Knepley   test:
5295be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9
5305be41b18SMatthew G. Knepley     suffix: 2d_q2_shear_conv
53130602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2
5325be41b18SMatthew G. Knepley   test:
5335be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9
5345be41b18SMatthew G. Knepley     suffix: 2d_q3_shear_conv
53530602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2
5365be41b18SMatthew G. Knepley   test:
5370fdc7489SMatthew Knepley     # Using -convest_num_refine 3 we get L_2 convergence rate: 1.7
5385be41b18SMatthew G. Knepley     suffix: 3d_p1_conv
5395be41b18SMatthew G. Knepley     requires: ctetgen
54030602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_refine 1 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1
5415be41b18SMatthew G. Knepley   test:
5425be41b18SMatthew G. Knepley     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 2.8
5435be41b18SMatthew G. Knepley     suffix: 3d_p2_conv
5445be41b18SMatthew G. Knepley     requires: ctetgen
54530602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 1
5465be41b18SMatthew G. Knepley   test:
5475be41b18SMatthew G. Knepley     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 4.0
5485be41b18SMatthew G. Knepley     suffix: 3d_p3_conv
5495be41b18SMatthew G. Knepley     requires: ctetgen
55030602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 1
5515be41b18SMatthew G. Knepley   test:
5525be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.8
5535be41b18SMatthew G. Knepley     suffix: 3d_q1_conv
55430602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_refine 1 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1
5555be41b18SMatthew G. Knepley   test:
5565be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.8
5575be41b18SMatthew G. Knepley     suffix: 3d_q2_conv
55830602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 1
5595be41b18SMatthew G. Knepley   test:
5605be41b18SMatthew G. Knepley     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 3.8
5615be41b18SMatthew G. Knepley     suffix: 3d_q3_conv
56230602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 1
5639139b27eSToby Isaac   test:
5649139b27eSToby Isaac     suffix: 2d_p1_fas_full
5659139b27eSToby Isaac     requires: triangle
566e600fa54SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_refine_hierarchy 5 \
5679139b27eSToby Isaac       -snes_max_it 1 -snes_type fas -snes_fas_levels 5 -snes_fas_type full -snes_fas_full_total \
5689139b27eSToby Isaac         -fas_coarse_snes_monitor -fas_coarse_snes_max_it 1 -fas_coarse_ksp_atol 1.e-13 \
5699139b27eSToby Isaac         -fas_levels_snes_monitor -fas_levels_snes_max_it 1 -fas_levels_snes_type newtonls \
5709139b27eSToby Isaac           -fas_levels_pc_type none -fas_levels_ksp_max_it 2 -fas_levels_ksp_converged_maxits -fas_levels_ksp_type chebyshev \
57173f7197eSJed Brown             -fas_levels_esteig_ksp_type cg -fas_levels_ksp_chebyshev_esteig 0,0.25,0,1.1 -fas_levels_esteig_ksp_max_it 10
5729139b27eSToby Isaac   test:
5739139b27eSToby Isaac     suffix: 2d_p1_fas_full_homogeneous
5749139b27eSToby Isaac     requires: triangle
575e600fa54SMatthew G. Knepley     args: -homogeneous -potential_petscspace_degree 1 -dm_refine_hierarchy 5 \
5769139b27eSToby Isaac       -snes_max_it 1 -snes_type fas -snes_fas_levels 5 -snes_fas_type full \
5779139b27eSToby Isaac         -fas_coarse_snes_monitor -fas_coarse_snes_max_it 1 -fas_coarse_ksp_atol 1.e-13 \
5789139b27eSToby Isaac         -fas_levels_snes_monitor -fas_levels_snes_max_it 1 -fas_levels_snes_type newtonls \
5799139b27eSToby Isaac           -fas_levels_pc_type none -fas_levels_ksp_max_it 2 -fas_levels_ksp_converged_maxits -fas_levels_ksp_type chebyshev \
58073f7197eSJed Brown             -fas_levels_esteig_ksp_type cg -fas_levels_ksp_chebyshev_esteig 0,0.25,0,1.1 -fas_levels_esteig_ksp_max_it 10
5815be41b18SMatthew G. Knepley 
582c4762a1bSJed Brown   test:
583c4762a1bSJed Brown     suffix: 2d_p1_scalable
5845be41b18SMatthew G. Knepley     requires: triangle
5855be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_refine 3 \
586c4762a1bSJed Brown       -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned \
58773f7197eSJed Brown       -pc_type gamg -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 \
588c4762a1bSJed Brown         -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \
589c4762a1bSJed Brown         -pc_gamg_coarse_eq_limit 1000 \
590c4762a1bSJed Brown         -pc_gamg_threshold 0.05 \
591c4762a1bSJed Brown         -pc_gamg_threshold_scale .0 \
592c4762a1bSJed Brown         -mg_levels_ksp_type chebyshev \
593c4762a1bSJed Brown         -mg_levels_ksp_max_it 1 \
594c4762a1bSJed Brown         -mg_levels_pc_type jacobi \
595c4762a1bSJed Brown       -matptap_via scalable
596c4762a1bSJed Brown   test:
597c4762a1bSJed Brown     suffix: 2d_p1_gmg_vcycle
598c4762a1bSJed Brown     requires: triangle
5995be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
600c4762a1bSJed Brown           -ksp_rtol 5e-10 -pc_type mg \
601c4762a1bSJed Brown             -mg_levels_ksp_max_it 1 \
602c4762a1bSJed Brown             -mg_levels_esteig_ksp_type cg \
603c4762a1bSJed Brown             -mg_levels_esteig_ksp_max_it 10 \
60473f7197eSJed Brown             -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \
605c4762a1bSJed Brown             -mg_levels_pc_type jacobi
606c4762a1bSJed Brown   test:
607c4762a1bSJed Brown     suffix: 2d_p1_gmg_fcycle
608c4762a1bSJed Brown     requires: triangle
6095be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
610c4762a1bSJed Brown           -ksp_rtol 5e-10 -pc_type mg -pc_mg_type full \
611c4762a1bSJed Brown             -mg_levels_ksp_max_it 2 \
612c4762a1bSJed Brown             -mg_levels_esteig_ksp_type cg \
613c4762a1bSJed Brown             -mg_levels_esteig_ksp_max_it 10 \
61473f7197eSJed Brown             -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \
615c4762a1bSJed Brown             -mg_levels_pc_type jacobi
616c4762a1bSJed Brown   test:
617d6837840SMatthew G. Knepley     suffix: 2d_p1_gmg_vcycle_adapt
6182b3cbbdaSStefano Zampini     requires: triangle
619908b9b43SStefano Zampini     args: -petscpartitioner_type simple -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
6202b3cbbdaSStefano Zampini           -ksp_rtol 5e-10 -pc_type mg -pc_mg_galerkin -pc_mg_adapt_interp_coarse_space harmonic -pc_mg_adapt_interp_n 8 \
621d6837840SMatthew G. Knepley             -mg_levels_ksp_max_it 1 \
622d6837840SMatthew G. Knepley             -mg_levels_esteig_ksp_type cg \
623d6837840SMatthew G. Knepley             -mg_levels_esteig_ksp_max_it 10 \
62473f7197eSJed Brown             -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \
625d6837840SMatthew G. Knepley             -mg_levels_pc_type jacobi
626d6837840SMatthew G. Knepley   test:
627c4762a1bSJed Brown     suffix: 2d_p1_spectral_0
628c4762a1bSJed Brown     requires: triangle fftw !complex
6295be41b18SMatthew G. Knepley     args: -dm_plex_box_faces 1,1 -potential_petscspace_degree 1 -dm_refine 6 -spectral -fft_view
630c4762a1bSJed Brown   test:
631c4762a1bSJed Brown     suffix: 2d_p1_spectral_1
632c4762a1bSJed Brown     requires: triangle fftw !complex
633c4762a1bSJed Brown     nsize: 2
634e600fa54SMatthew G. Knepley     args: -dm_plex_box_faces 4,4 -potential_petscspace_degree 1 -spectral -fft_view
635c4762a1bSJed Brown   test:
636c4762a1bSJed Brown     suffix: 2d_p1_adj_0
637c4762a1bSJed Brown     requires: triangle
6385be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_refine 1 -adjoint -adjoint_petscspace_degree 1 -error_petscspace_degree 0
639b1ad28e3SJunchao Zhang   test:
640b1ad28e3SJunchao Zhang     nsize: 2
641dcfd994dSJunchao Zhang     requires: kokkos_kernels
642b1ad28e3SJunchao Zhang     suffix: kokkos
643e600fa54SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_box_faces 2,3,6 -petscpartitioner_type simple -dm_plex_simplex 0 -potential_petscspace_degree 1 \
644b1ad28e3SJunchao Zhang          -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 \
64573f7197eSJed Brown          -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 \
64673f7197eSJed Brown          -ksp_monitor -snes_monitor -dm_view -dm_mat_type aijkokkos -dm_vec_type kokkos
647c4762a1bSJed Brown 
648c4762a1bSJed Brown TEST*/
649