xref: /petsc/src/snes/tutorials/ex13.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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 
21c4762a1bSJed Brown static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
22c4762a1bSJed Brown {
23c4762a1bSJed Brown   *u = 0.0;
24c4762a1bSJed Brown   return 0;
25c4762a1bSJed Brown }
26c4762a1bSJed Brown 
27b724ec41SToby Isaac static PetscErrorCode trig_inhomogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
28c4762a1bSJed Brown {
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 
35b724ec41SToby Isaac static PetscErrorCode trig_homogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
36b724ec41SToby Isaac {
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) */
44c4762a1bSJed Brown static void obj_error_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
45c4762a1bSJed Brown                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
46c4762a1bSJed Brown                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
47c4762a1bSJed Brown                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar obj[])
48c4762a1bSJed Brown {
49c4762a1bSJed Brown   obj[0] = a[aOff[0]]*(u[0] - a[aOff[1]]);
50c4762a1bSJed Brown }
51c4762a1bSJed Brown 
52b724ec41SToby Isaac static void f0_trig_inhomogeneous_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
53c4762a1bSJed Brown                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
54c4762a1bSJed Brown                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
55c4762a1bSJed Brown                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
56c4762a1bSJed Brown {
57c4762a1bSJed Brown   PetscInt d;
58c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f0[0] += -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]);
59c4762a1bSJed Brown }
60c4762a1bSJed Brown 
61b724ec41SToby Isaac static void f0_trig_homogeneous_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
62b724ec41SToby Isaac                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
63b724ec41SToby Isaac                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
64b724ec41SToby Isaac                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
65b724ec41SToby Isaac {
66b724ec41SToby Isaac   PetscInt d;
67b724ec41SToby Isaac   for (d = 0; d < dim; ++d) {
68b724ec41SToby Isaac     PetscScalar v = 1.;
69b724ec41SToby Isaac     for (PetscInt e = 0; e < dim; e++) {
70b724ec41SToby Isaac       if (e == d) {
71b724ec41SToby Isaac         v *= -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]);
72b724ec41SToby Isaac       } else {
73b724ec41SToby Isaac         v *= PetscSinReal(2.0*PETSC_PI*x[d]);
74b724ec41SToby Isaac       }
75b724ec41SToby Isaac     }
76b724ec41SToby Isaac     f0[0] += v;
77b724ec41SToby Isaac   }
78b724ec41SToby Isaac }
79b724ec41SToby Isaac 
80c4762a1bSJed Brown static void f0_unity_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
81c4762a1bSJed Brown                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
82c4762a1bSJed Brown                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
83c4762a1bSJed Brown                        PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
84c4762a1bSJed Brown {
85c4762a1bSJed Brown   f0[0] = 1.0;
86c4762a1bSJed Brown }
87c4762a1bSJed Brown 
88c4762a1bSJed Brown static void f0_identityaux_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
89c4762a1bSJed Brown                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
90c4762a1bSJed Brown                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
91c4762a1bSJed Brown                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
92c4762a1bSJed Brown {
93c4762a1bSJed Brown   f0[0] = a[0];
94c4762a1bSJed Brown }
95c4762a1bSJed Brown 
96c4762a1bSJed Brown static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
97c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
98c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
99c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
100c4762a1bSJed Brown {
101c4762a1bSJed Brown   PetscInt d;
102c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
103c4762a1bSJed Brown }
104c4762a1bSJed Brown 
105c4762a1bSJed Brown static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
106c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
107c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
108c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
109c4762a1bSJed Brown {
110c4762a1bSJed Brown   PetscInt d;
111c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
112c4762a1bSJed Brown }
113c4762a1bSJed Brown 
114c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
115c4762a1bSJed Brown {
116c4762a1bSJed Brown   PetscErrorCode ierr;
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   PetscFunctionBeginUser;
119c4762a1bSJed Brown   options->shear       = PETSC_FALSE;
120c4762a1bSJed Brown   options->spectral    = PETSC_FALSE;
121c4762a1bSJed Brown   options->adjoint     = PETSC_FALSE;
122b724ec41SToby Isaac   options->homogeneous = PETSC_FALSE;
123ad1a7c93SMatthew Knepley   options->viewError   = PETSC_FALSE;
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");CHKERRQ(ierr);
126*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-shear", "Shear the domain", "ex13.c", options->shear, &options->shear, NULL));
127*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-spectral", "Look at the spectrum along planes of the solution", "ex13.c", options->spectral, &options->spectral, NULL));
128*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-adjoint", "Solve the adjoint problem", "ex13.c", options->adjoint, &options->adjoint, NULL));
129*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-homogeneous", "Use homogeneous boundary conditions", "ex13.c", options->homogeneous, &options->homogeneous, NULL));
130*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-error_view", "Output the solution error", "ex13.c", options->viewError, &options->viewError, NULL));
131ad1a7c93SMatthew Knepley   ierr = PetscOptionsEnd();
132c4762a1bSJed Brown   PetscFunctionReturn(0);
133c4762a1bSJed Brown }
134c4762a1bSJed Brown 
135c4762a1bSJed Brown static PetscErrorCode CreateSpectralPlanes(DM dm, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user)
136c4762a1bSJed Brown {
137c4762a1bSJed Brown   PetscSection       coordSection;
138c4762a1bSJed Brown   Vec                coordinates;
139c4762a1bSJed Brown   const PetscScalar *coords;
140c4762a1bSJed Brown   PetscInt           dim, p, vStart, vEnd, v;
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   PetscFunctionBeginUser;
143*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDim(dm, &dim));
144*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
145*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(dm, &coordinates));
146*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateSection(dm, &coordSection));
147*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(coordinates, &coords));
148c4762a1bSJed Brown   for (p = 0; p < numPlanes; ++p) {
149c4762a1bSJed Brown     DMLabel label;
150c4762a1bSJed Brown     char    name[PETSC_MAX_PATH_LEN];
151c4762a1bSJed Brown 
152*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%D", p));
153*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCreateLabel(dm, name));
154*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetLabel(dm, name, &label));
155*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLabelAddStratum(label, 1));
156c4762a1bSJed Brown     for (v = vStart; v < vEnd; ++v) {
157c4762a1bSJed Brown       PetscInt off;
158c4762a1bSJed Brown 
159*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetOffset(coordSection, v, &off));
160c4762a1bSJed Brown       if (PetscAbsReal(planeCoord[p] - PetscRealPart(coords[off+planeDir[p]])) < PETSC_SMALL) {
161*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMLabelSetValue(label, v, 1));
162c4762a1bSJed Brown       }
163c4762a1bSJed Brown     }
164c4762a1bSJed Brown   }
165*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(coordinates, &coords));
166c4762a1bSJed Brown   PetscFunctionReturn(0);
167c4762a1bSJed Brown }
168c4762a1bSJed Brown 
169c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
170c4762a1bSJed Brown {
171c4762a1bSJed Brown   PetscFunctionBeginUser;
172*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(comm, dm));
173*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*dm, DMPLEX));
174*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*dm));
175*5f80ce2aSJacob Faibussowitsch   if (user->shear) CHKERRQ(DMPlexShearGeometry(*dm, DM_X, NULL));
176*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetApplicationContext(*dm, user));
177*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view"));
178c4762a1bSJed Brown   if (user->spectral) {
179c4762a1bSJed Brown     PetscInt  planeDir[2]   = {0,  1};
180c4762a1bSJed Brown     PetscReal planeCoord[2] = {0., 1.};
181c4762a1bSJed Brown 
182*5f80ce2aSJacob Faibussowitsch     CHKERRQ(CreateSpectralPlanes(*dm, 2, planeDir, planeCoord, user));
183c4762a1bSJed Brown   }
184c4762a1bSJed Brown   PetscFunctionReturn(0);
185c4762a1bSJed Brown }
186c4762a1bSJed Brown 
187c4762a1bSJed Brown static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
188c4762a1bSJed Brown {
18945480ffeSMatthew G. Knepley   PetscDS        ds;
19045480ffeSMatthew G. Knepley   DMLabel        label;
191c4762a1bSJed Brown   const PetscInt id = 1;
192b724ec41SToby Isaac   PetscPointFunc f0 = user->homogeneous ? f0_trig_homogeneous_u : f0_trig_inhomogeneous_u;
193b724ec41SToby Isaac   PetscErrorCode (*ex)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *) = user->homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u;
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   PetscFunctionBeginUser;
196*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDS(dm, &ds));
197*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetResidual(ds, 0, f0, f1_u));
198*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu));
199*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetExactSolution(ds, 0, ex, user));
200*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLabel(dm, "marker", &label));
201*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) ex, NULL, user, NULL));
202c4762a1bSJed Brown   PetscFunctionReturn(0);
203c4762a1bSJed Brown }
204c4762a1bSJed Brown 
205c4762a1bSJed Brown static PetscErrorCode SetupAdjointProblem(DM dm, AppCtx *user)
206c4762a1bSJed Brown {
20745480ffeSMatthew G. Knepley   PetscDS        ds;
20845480ffeSMatthew G. Knepley   DMLabel        label;
209c4762a1bSJed Brown   const PetscInt id = 1;
210c4762a1bSJed Brown 
211c4762a1bSJed Brown   PetscFunctionBeginUser;
212*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDS(dm, &ds));
213*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetResidual(ds, 0, f0_unity_u, f1_u));
214*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu));
215*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetObjective(ds, 0, obj_error_u));
216*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLabel(dm, "marker", &label));
217*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) zero, NULL, user, NULL));
218c4762a1bSJed Brown   PetscFunctionReturn(0);
219c4762a1bSJed Brown }
220c4762a1bSJed Brown 
221c4762a1bSJed Brown static PetscErrorCode SetupErrorProblem(DM dm, AppCtx *user)
222c4762a1bSJed Brown {
223c4762a1bSJed Brown   PetscDS        prob;
224c4762a1bSJed Brown 
225c4762a1bSJed Brown   PetscFunctionBeginUser;
226*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDS(dm, &prob));
227c4762a1bSJed Brown   PetscFunctionReturn(0);
228c4762a1bSJed Brown }
229c4762a1bSJed Brown 
230c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
231c4762a1bSJed Brown {
232c4762a1bSJed Brown   DM             cdm = dm;
233c4762a1bSJed Brown   PetscFE        fe;
2345be41b18SMatthew G. Knepley   DMPolytopeType ct;
2355be41b18SMatthew G. Knepley   PetscBool      simplex;
2365be41b18SMatthew G. Knepley   PetscInt       dim, cStart;
237c4762a1bSJed Brown   char           prefix[PETSC_MAX_PATH_LEN];
238c4762a1bSJed Brown 
239c4762a1bSJed Brown   PetscFunctionBeginUser;
240*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
241*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
242*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetCellType(dm, cStart, &ct));
2435be41b18SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
244c4762a1bSJed Brown   /* Create finite element */
245*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
246*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe));
247*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) fe, name));
248c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
249*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe));
250*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateDS(dm));
251*5f80ce2aSJacob Faibussowitsch   CHKERRQ((*setup)(dm, user));
252c4762a1bSJed Brown   while (cdm) {
253*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCopyDisc(dm,cdm));
254c4762a1bSJed Brown     /* TODO: Check whether the boundary of coarse meshes is marked */
255*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoarseDM(cdm, &cdm));
256c4762a1bSJed Brown   }
257*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEDestroy(&fe));
258c4762a1bSJed Brown   PetscFunctionReturn(0);
259c4762a1bSJed Brown }
260c4762a1bSJed Brown 
261c4762a1bSJed Brown static PetscErrorCode ComputeSpectral(DM dm, Vec u, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user)
262c4762a1bSJed Brown {
263c4762a1bSJed Brown   MPI_Comm           comm;
264c4762a1bSJed Brown   PetscSection       coordSection, section;
265c4762a1bSJed Brown   Vec                coordinates, uloc;
266c4762a1bSJed Brown   const PetscScalar *coords, *array;
267c4762a1bSJed Brown   PetscInt           p;
268c4762a1bSJed Brown   PetscMPIInt        size, rank;
269c4762a1bSJed Brown 
270c4762a1bSJed Brown   PetscFunctionBeginUser;
271*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm));
272*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm, &size));
273*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm, &rank));
274*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(dm, &uloc));
275*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(dm, u, INSERT_VALUES, uloc));
276*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(dm, u, INSERT_VALUES, uloc));
277*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, uloc, 0.0, NULL, NULL, NULL));
278*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(uloc, NULL, "-sol_view"));
279*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalSection(dm, &section));
280*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(uloc, &array));
281*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(dm, &coordinates));
282*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateSection(dm, &coordSection));
283*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(coordinates, &coords));
284c4762a1bSJed Brown   for (p = 0; p < numPlanes; ++p) {
285c4762a1bSJed Brown     DMLabel         label;
286c4762a1bSJed Brown     char            name[PETSC_MAX_PATH_LEN];
287c4762a1bSJed Brown     Mat             F;
288c4762a1bSJed Brown     Vec             x, y;
289c4762a1bSJed Brown     IS              stratum;
290c4762a1bSJed Brown     PetscReal      *ray, *gray;
291c4762a1bSJed Brown     PetscScalar    *rvals, *svals, *gsvals;
292c4762a1bSJed Brown     PetscInt       *perm, *nperm;
293c4762a1bSJed Brown     PetscInt        n, N, i, j, off, offu;
294c4762a1bSJed Brown     const PetscInt *points;
295c4762a1bSJed Brown 
296*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%D", p));
297*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetLabel(dm, name, &label));
298*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLabelGetStratumIS(label, 1, &stratum));
299*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetLocalSize(stratum, &n));
300*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetIndices(stratum, &points));
301*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc2(n, &ray, n, &svals));
302c4762a1bSJed Brown     for (i = 0; i < n; ++i) {
303*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetOffset(coordSection, points[i], &off));
304*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetOffset(section, points[i], &offu));
305c4762a1bSJed Brown       ray[i]   = PetscRealPart(coords[off+((planeDir[p]+1)%2)]);
306c4762a1bSJed Brown       svals[i] = array[offu];
307c4762a1bSJed Brown     }
308c4762a1bSJed Brown     /* Gather the ray data to proc 0 */
309c4762a1bSJed Brown     if (size > 1) {
310c4762a1bSJed Brown       PetscMPIInt *cnt, *displs, p;
311c4762a1bSJed Brown 
312*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscCalloc2(size, &cnt, size, &displs));
313*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Gather(&n, 1, MPIU_INT, cnt, 1, MPIU_INT, 0, comm));
314c4762a1bSJed Brown       for (p = 1; p < size; ++p) displs[p] = displs[p-1] + cnt[p-1];
315c4762a1bSJed Brown       N = displs[size-1] + cnt[size-1];
316*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc2(N, &gray, N, &gsvals));
317*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Gatherv(ray, n, MPIU_REAL, gray, cnt, displs, MPIU_REAL, 0, comm));
318*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Gatherv(svals, n, MPIU_SCALAR, gsvals, cnt, displs, MPIU_SCALAR, 0, comm));
319*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree2(cnt, displs));
320c4762a1bSJed Brown     } else {
321c4762a1bSJed Brown       N      = n;
322c4762a1bSJed Brown       gray   = ray;
323c4762a1bSJed Brown       gsvals = svals;
324c4762a1bSJed Brown     }
325dd400576SPatrick Sanan     if (rank == 0) {
326c4762a1bSJed Brown       /* Sort point along ray */
327*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc2(N, &perm, N, &nperm));
328c4762a1bSJed Brown       for (i = 0; i < N; ++i) {perm[i] = i;}
329*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSortRealWithPermutation(N, gray, perm));
330c4762a1bSJed Brown       /* Count duplicates and squish mapping */
331c4762a1bSJed Brown       nperm[0] = perm[0];
332c4762a1bSJed Brown       for (i = 1, j = 1; i < N; ++i) {
333c4762a1bSJed Brown         if (PetscAbsReal(gray[perm[i]] - gray[perm[i-1]]) > PETSC_SMALL) nperm[j++] = perm[i];
334c4762a1bSJed Brown       }
335c4762a1bSJed Brown       /* Create FFT structs */
336*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreateFFT(PETSC_COMM_SELF, 1, &j, MATFFTW, &F));
337*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreateVecs(F, &x, &y));
338*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectSetName((PetscObject) y, name));
339*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecGetArray(x, &rvals));
340c4762a1bSJed Brown       for (i = 0, j = 0; i < N; ++i) {
341c4762a1bSJed Brown         if (i > 0 && PetscAbsReal(gray[perm[i]] - gray[perm[i-1]]) < PETSC_SMALL) continue;
342c4762a1bSJed Brown         rvals[j] = gsvals[nperm[j]];
343c4762a1bSJed Brown         ++j;
344c4762a1bSJed Brown       }
345*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree2(perm, nperm));
346*5f80ce2aSJacob Faibussowitsch       if (size > 1) CHKERRQ(PetscFree2(gray, gsvals));
347*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecRestoreArray(x, &rvals));
348c4762a1bSJed Brown       /* Do FFT along the ray */
349*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(F, x, y));
350c4762a1bSJed Brown       /* Chop FFT */
351*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecChop(y, PETSC_SMALL));
352*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecViewFromOptions(x, NULL, "-real_view"));
353*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecViewFromOptions(y, NULL, "-fft_view"));
354*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&x));
355*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&y));
356*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&F));
357c4762a1bSJed Brown     }
358*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISRestoreIndices(stratum, &points));
359*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&stratum));
360*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree2(ray, svals));
361c4762a1bSJed Brown   }
362*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(coordinates, &coords));
363*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(uloc, &array));
364*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(dm, &uloc));
365c4762a1bSJed Brown   PetscFunctionReturn(0);
366c4762a1bSJed Brown }
367c4762a1bSJed Brown 
368c4762a1bSJed Brown int main(int argc, char **argv)
369c4762a1bSJed Brown {
370c4762a1bSJed Brown   DM             dm;   /* Problem specification */
371c4762a1bSJed Brown   SNES           snes; /* Nonlinear solver */
372c4762a1bSJed Brown   Vec            u;    /* Solutions */
373c4762a1bSJed Brown   AppCtx         user; /* User-defined work context */
374c4762a1bSJed Brown   PetscErrorCode ierr;
375c4762a1bSJed Brown 
376c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
377*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user));
378c4762a1bSJed Brown   /* Primal system */
379*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(PETSC_COMM_WORLD, &snes));
380*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
381*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetDM(snes, dm));
382*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user));
383*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(dm, &u));
384*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(u, 0.0));
385*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) u, "potential"));
386*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexSetSNESLocalFEM(dm, &user, &user, &user));
387*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(snes));
388*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes, NULL, u));
389*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetSolution(snes, &u));
390*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(u, NULL, "-potential_view"));
391ad1a7c93SMatthew Knepley   if (user.viewError) {
392ad1a7c93SMatthew Knepley     PetscErrorCode (*sol)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
393ad1a7c93SMatthew Knepley     void            *ctx;
394ad1a7c93SMatthew Knepley     PetscDS          ds;
395ad1a7c93SMatthew Knepley     PetscReal        error;
396ad1a7c93SMatthew Knepley     PetscInt         N;
397ad1a7c93SMatthew Knepley 
398*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetDS(dm, &ds));
399*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSGetExactSolution(ds, 0, &sol, &ctx));
400*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetSize(u, &N));
401*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMComputeL2Diff(dm, 0.0, &sol, &ctx, u, &error));
402*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "N: %D L2 error: %g\n", N, (double)error));
403ad1a7c93SMatthew Knepley   }
404c4762a1bSJed Brown   if (user.spectral) {
405c4762a1bSJed Brown     PetscInt  planeDir[2]   = {0,  1};
406c4762a1bSJed Brown     PetscReal planeCoord[2] = {0., 1.};
407c4762a1bSJed Brown 
408*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ComputeSpectral(dm, u, 2, planeDir, planeCoord, &user));
409c4762a1bSJed Brown   }
410c4762a1bSJed Brown   /* Adjoint system */
411c4762a1bSJed Brown   if (user.adjoint) {
412c4762a1bSJed Brown     DM   dmAdj;
413c4762a1bSJed Brown     SNES snesAdj;
414c4762a1bSJed Brown     Vec  uAdj;
415c4762a1bSJed Brown 
416*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESCreate(PETSC_COMM_WORLD, &snesAdj));
417*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) snesAdj, "adjoint_"));
418*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMClone(dm, &dmAdj));
419*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESSetDM(snesAdj, dmAdj));
420*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SetupDiscretization(dmAdj, "adjoint", SetupAdjointProblem, &user));
421*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCreateGlobalVector(dmAdj, &uAdj));
422*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(uAdj, 0.0));
423*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject) uAdj, "adjoint"));
424*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexSetSNESLocalFEM(dmAdj, &user, &user, &user));
425*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESSetFromOptions(snesAdj));
426*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESSolve(snesAdj, NULL, uAdj));
427*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESGetSolution(snesAdj, &uAdj));
428*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecViewFromOptions(uAdj, NULL, "-adjoint_view"));
429c4762a1bSJed Brown     /* Error representation */
430c4762a1bSJed Brown     {
431c4762a1bSJed Brown       DM        dmErr, dmErrAux, dms[2];
432c4762a1bSJed Brown       Vec       errorEst, errorL2, uErr, uErrLoc, uAdjLoc, uAdjProj;
433c4762a1bSJed Brown       IS       *subis;
434c4762a1bSJed Brown       PetscReal errorEstTot, errorL2Norm, errorL2Tot;
435c4762a1bSJed Brown       PetscInt  N, i;
436b724ec41SToby Isaac       PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *) = {user.homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u};
437c4762a1bSJed Brown       void (*identity[1])(PetscInt, PetscInt, PetscInt,
438c4762a1bSJed Brown                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
439c4762a1bSJed Brown                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
440c4762a1bSJed Brown                           PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]) = {f0_identityaux_u};
441c4762a1bSJed Brown       void            *ctxs[1] = {0};
442c4762a1bSJed Brown 
443c4762a1bSJed Brown       ctxs[0] = &user;
444*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMClone(dm, &dmErr));
445*5f80ce2aSJacob Faibussowitsch       CHKERRQ(SetupDiscretization(dmErr, "error", SetupErrorProblem, &user));
446*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGetGlobalVector(dmErr, &errorEst));
447*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGetGlobalVector(dmErr, &errorL2));
448c4762a1bSJed Brown       /*   Compute auxiliary data (solution and projection of adjoint solution) */
449*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGetLocalVector(dmAdj, &uAdjLoc));
450*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGlobalToLocalBegin(dmAdj, uAdj, INSERT_VALUES, uAdjLoc));
451*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGlobalToLocalEnd(dmAdj, uAdj, INSERT_VALUES, uAdjLoc));
452*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGetGlobalVector(dm, &uAdjProj));
453*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSetAuxiliaryVec(dm, NULL, 0, 0, uAdjLoc));
454*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMProjectField(dm, 0.0, u, identity, INSERT_VALUES, uAdjProj));
455*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSetAuxiliaryVec(dm, NULL, 0, 0, NULL));
456*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMRestoreLocalVector(dmAdj, &uAdjLoc));
457c4762a1bSJed Brown       /*   Attach auxiliary data */
458c4762a1bSJed Brown       dms[0] = dm; dms[1] = dm;
459*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMCreateSuperDM(dms, 2, &subis, &dmErrAux));
460c4762a1bSJed Brown       if (0) {
461c4762a1bSJed Brown         PetscSection sec;
462c4762a1bSJed Brown 
463*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMGetLocalSection(dms[0], &sec));
464*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD));
465*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMGetLocalSection(dms[1], &sec));
466*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD));
467*5f80ce2aSJacob Faibussowitsch         CHKERRQ(DMGetLocalSection(dmErrAux, &sec));
468*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD));
469c4762a1bSJed Brown       }
470*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMViewFromOptions(dmErrAux, NULL, "-dm_err_view"));
471*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISViewFromOptions(subis[0], NULL, "-super_is_view"));
472*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISViewFromOptions(subis[1], NULL, "-super_is_view"));
473*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGetGlobalVector(dmErrAux, &uErr));
474*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecViewFromOptions(u, NULL, "-map_vec_view"));
475*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecViewFromOptions(uAdjProj, NULL, "-map_vec_view"));
476*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecViewFromOptions(uErr, NULL, "-map_vec_view"));
477*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecISCopy(uErr, subis[0], SCATTER_FORWARD, u));
478*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecISCopy(uErr, subis[1], SCATTER_FORWARD, uAdjProj));
479*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMRestoreGlobalVector(dm, &uAdjProj));
480*5f80ce2aSJacob Faibussowitsch       for (i = 0; i < 2; ++i) CHKERRQ(ISDestroy(&subis[i]));
481*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(subis));
482*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGetLocalVector(dmErrAux, &uErrLoc));
483*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGlobalToLocalBegin(dm, uErr, INSERT_VALUES, uErrLoc));
484*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGlobalToLocalEnd(dm, uErr, INSERT_VALUES, uErrLoc));
485*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMRestoreGlobalVector(dmErrAux, &uErr));
486*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSetAuxiliaryVec(dmAdj, NULL, 0, 0, uErrLoc));
487c4762a1bSJed Brown       /*   Compute cellwise error estimate */
488*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSet(errorEst, 0.0));
489*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexComputeCellwiseIntegralFEM(dmAdj, uAdj, errorEst, &user));
490*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSetAuxiliaryVec(dmAdj, NULL, 0, 0, NULL));
491*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMRestoreLocalVector(dmErrAux, &uErrLoc));
492*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMDestroy(&dmErrAux));
493c4762a1bSJed Brown       /*   Plot cellwise error vector */
494*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecViewFromOptions(errorEst, NULL, "-error_view"));
495c4762a1bSJed Brown       /*   Compute ratio of estimate (sum over cells) with actual L_2 error */
496*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMComputeL2Diff(dm, 0.0, funcs, ctxs, u, &errorL2Norm));
497*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexComputeL2DiffVec(dm, 0.0, funcs, ctxs, u, errorL2));
498*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecViewFromOptions(errorL2, NULL, "-l2_error_view"));
499*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(errorL2,  NORM_INFINITY, &errorL2Tot));
500*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(errorEst, NORM_INFINITY, &errorEstTot));
501*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecGetSize(errorEst, &N));
502*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecPointwiseDivide(errorEst, errorEst, errorL2));
503*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectSetName((PetscObject) errorEst, "Error ratio"));
504*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecViewFromOptions(errorEst, NULL, "-error_ratio_view"));
505*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "N: %D L2 error: %g Error Ratio: %g/%g = %g\n", N, (double) errorL2Norm, (double) errorEstTot, (double) PetscSqrtReal(errorL2Tot), (double) errorEstTot/PetscSqrtReal(errorL2Tot)));
506*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMRestoreGlobalVector(dmErr, &errorEst));
507*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMRestoreGlobalVector(dmErr, &errorL2));
508*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMDestroy(&dmErr));
509c4762a1bSJed Brown     }
510*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDestroy(&dmAdj));
511*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&uAdj));
512*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESDestroy(&snesAdj));
513c4762a1bSJed Brown   }
514c4762a1bSJed Brown   /* Cleanup */
515*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
516*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESDestroy(&snes));
517*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
518c4762a1bSJed Brown   ierr = PetscFinalize();
519c4762a1bSJed Brown   return ierr;
520c4762a1bSJed Brown }
521c4762a1bSJed Brown 
522c4762a1bSJed Brown /*TEST
523c4762a1bSJed Brown 
524c4762a1bSJed Brown   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_p1_conv
527c4762a1bSJed Brown     requires: triangle
5285be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2
5295be41b18SMatthew G. Knepley   test:
5305be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9
5315be41b18SMatthew G. Knepley     suffix: 2d_p2_conv
5325be41b18SMatthew G. Knepley     requires: triangle
5335be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2
5345be41b18SMatthew G. Knepley   test:
5355be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9
5365be41b18SMatthew G. Knepley     suffix: 2d_p3_conv
5375be41b18SMatthew G. Knepley     requires: triangle
5385be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2
5395be41b18SMatthew G. Knepley   test:
5405be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9
5415be41b18SMatthew G. Knepley     suffix: 2d_q1_conv
54230602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2
5435be41b18SMatthew G. Knepley   test:
5445be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9
5455be41b18SMatthew G. Knepley     suffix: 2d_q2_conv
54630602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2
5475be41b18SMatthew G. Knepley   test:
5485be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9
5495be41b18SMatthew G. Knepley     suffix: 2d_q3_conv
55030602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2
5515be41b18SMatthew G. Knepley   test:
5525be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9
5535be41b18SMatthew G. Knepley     suffix: 2d_q1_shear_conv
55430602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2
5555be41b18SMatthew G. Knepley   test:
5565be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9
5575be41b18SMatthew G. Knepley     suffix: 2d_q2_shear_conv
55830602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2
5595be41b18SMatthew G. Knepley   test:
5605be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9
5615be41b18SMatthew G. Knepley     suffix: 2d_q3_shear_conv
56230602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2
5635be41b18SMatthew G. Knepley   test:
5640fdc7489SMatthew Knepley     # Using -convest_num_refine 3 we get L_2 convergence rate: 1.7
5655be41b18SMatthew G. Knepley     suffix: 3d_p1_conv
5665be41b18SMatthew G. Knepley     requires: ctetgen
56730602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_refine 1 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1
5685be41b18SMatthew G. Knepley   test:
5695be41b18SMatthew G. Knepley     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 2.8
5705be41b18SMatthew G. Knepley     suffix: 3d_p2_conv
5715be41b18SMatthew G. Knepley     requires: ctetgen
57230602db0SMatthew 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
5735be41b18SMatthew G. Knepley   test:
5745be41b18SMatthew G. Knepley     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 4.0
5755be41b18SMatthew G. Knepley     suffix: 3d_p3_conv
5765be41b18SMatthew G. Knepley     requires: ctetgen
57730602db0SMatthew 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
5785be41b18SMatthew G. Knepley   test:
5795be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.8
5805be41b18SMatthew G. Knepley     suffix: 3d_q1_conv
58130602db0SMatthew 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
5825be41b18SMatthew G. Knepley   test:
5835be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.8
5845be41b18SMatthew G. Knepley     suffix: 3d_q2_conv
58530602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 1
5865be41b18SMatthew G. Knepley   test:
5875be41b18SMatthew G. Knepley     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 3.8
5885be41b18SMatthew G. Knepley     suffix: 3d_q3_conv
58930602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 1
5909139b27eSToby Isaac   test:
5919139b27eSToby Isaac     suffix: 2d_p1_fas_full
5929139b27eSToby Isaac     requires: triangle
593e600fa54SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_refine_hierarchy 5 \
5949139b27eSToby Isaac       -snes_max_it 1 -snes_type fas -snes_fas_levels 5 -snes_fas_type full -snes_fas_full_total \
5959139b27eSToby Isaac         -fas_coarse_snes_monitor -fas_coarse_snes_max_it 1 -fas_coarse_ksp_atol 1.e-13 \
5969139b27eSToby Isaac         -fas_levels_snes_monitor -fas_levels_snes_max_it 1 -fas_levels_snes_type newtonls \
5979139b27eSToby Isaac           -fas_levels_pc_type none -fas_levels_ksp_max_it 2 -fas_levels_ksp_converged_maxits -fas_levels_ksp_type chebyshev \
59873f7197eSJed 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
5999139b27eSToby Isaac   test:
6009139b27eSToby Isaac     suffix: 2d_p1_fas_full_homogeneous
6019139b27eSToby Isaac     requires: triangle
602e600fa54SMatthew G. Knepley     args: -homogeneous -potential_petscspace_degree 1 -dm_refine_hierarchy 5 \
6039139b27eSToby Isaac       -snes_max_it 1 -snes_type fas -snes_fas_levels 5 -snes_fas_type full \
6049139b27eSToby Isaac         -fas_coarse_snes_monitor -fas_coarse_snes_max_it 1 -fas_coarse_ksp_atol 1.e-13 \
6059139b27eSToby Isaac         -fas_levels_snes_monitor -fas_levels_snes_max_it 1 -fas_levels_snes_type newtonls \
6069139b27eSToby Isaac           -fas_levels_pc_type none -fas_levels_ksp_max_it 2 -fas_levels_ksp_converged_maxits -fas_levels_ksp_type chebyshev \
60773f7197eSJed 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
6085be41b18SMatthew G. Knepley 
609c4762a1bSJed Brown   test:
610c4762a1bSJed Brown     suffix: 2d_p1_scalable
6115be41b18SMatthew G. Knepley     requires: triangle
6125be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_refine 3 \
613c4762a1bSJed Brown       -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned \
61473f7197eSJed Brown       -pc_type gamg -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 \
615c4762a1bSJed Brown         -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \
616c4762a1bSJed Brown         -pc_gamg_coarse_eq_limit 1000 \
617c4762a1bSJed Brown         -pc_gamg_square_graph 1 \
618c4762a1bSJed Brown         -pc_gamg_threshold 0.05 \
619c4762a1bSJed Brown         -pc_gamg_threshold_scale .0 \
620c4762a1bSJed Brown         -mg_levels_ksp_type chebyshev \
621c4762a1bSJed Brown         -mg_levels_ksp_max_it 1 \
622c4762a1bSJed Brown         -mg_levels_pc_type jacobi \
623c4762a1bSJed Brown       -matptap_via scalable
624c4762a1bSJed Brown   test:
625c4762a1bSJed Brown     suffix: 2d_p1_gmg_vcycle
626c4762a1bSJed Brown     requires: triangle
6275be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
628c4762a1bSJed Brown           -ksp_rtol 5e-10 -pc_type mg \
629c4762a1bSJed Brown             -mg_levels_ksp_max_it 1 \
630c4762a1bSJed Brown             -mg_levels_esteig_ksp_type cg \
631c4762a1bSJed Brown             -mg_levels_esteig_ksp_max_it 10 \
63273f7197eSJed Brown             -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \
633c4762a1bSJed Brown             -mg_levels_pc_type jacobi
634c4762a1bSJed Brown   test:
635c4762a1bSJed Brown     suffix: 2d_p1_gmg_fcycle
636c4762a1bSJed Brown     requires: triangle
6375be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
638c4762a1bSJed Brown           -ksp_rtol 5e-10 -pc_type mg -pc_mg_type full \
639c4762a1bSJed Brown             -mg_levels_ksp_max_it 2 \
640c4762a1bSJed Brown             -mg_levels_esteig_ksp_type cg \
641c4762a1bSJed Brown             -mg_levels_esteig_ksp_max_it 10 \
64273f7197eSJed Brown             -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \
643c4762a1bSJed Brown             -mg_levels_pc_type jacobi
644c4762a1bSJed Brown   test:
645d6837840SMatthew G. Knepley     suffix: 2d_p1_gmg_vcycle_adapt
646d6837840SMatthew G. Knepley     requires: triangle bamg
6475be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
648d6837840SMatthew G. Knepley           -ksp_rtol 5e-10 -pc_type mg -pc_mg_galerkin -pc_mg_adapt_interp -pc_mg_adapt_interp_coarse_space harmonic -pc_mg_adapt_interp_n 8 \
649d6837840SMatthew G. Knepley             -mg_levels_ksp_max_it 1 \
650d6837840SMatthew G. Knepley             -mg_levels_esteig_ksp_type cg \
651d6837840SMatthew G. Knepley             -mg_levels_esteig_ksp_max_it 10 \
65273f7197eSJed Brown             -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \
653d6837840SMatthew G. Knepley             -mg_levels_pc_type jacobi
654d6837840SMatthew G. Knepley   test:
655c4762a1bSJed Brown     suffix: 2d_p1_spectral_0
656c4762a1bSJed Brown     requires: triangle fftw !complex
6575be41b18SMatthew G. Knepley     args: -dm_plex_box_faces 1,1 -potential_petscspace_degree 1 -dm_refine 6 -spectral -fft_view
658c4762a1bSJed Brown   test:
659c4762a1bSJed Brown     suffix: 2d_p1_spectral_1
660c4762a1bSJed Brown     requires: triangle fftw !complex
661c4762a1bSJed Brown     nsize: 2
662e600fa54SMatthew G. Knepley     args: -dm_plex_box_faces 4,4 -potential_petscspace_degree 1 -spectral -fft_view
663c4762a1bSJed Brown   test:
664c4762a1bSJed Brown     suffix: 2d_p1_adj_0
665c4762a1bSJed Brown     requires: triangle
6665be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_refine 1 -adjoint -adjoint_petscspace_degree 1 -error_petscspace_degree 0
667b1ad28e3SJunchao Zhang   test:
668b1ad28e3SJunchao Zhang     nsize: 2
6693078479eSJunchao Zhang     requires: !sycl kokkos_kernels
670b1ad28e3SJunchao Zhang     suffix: kokkos
671e600fa54SMatthew 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 \
672b1ad28e3SJunchao 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 \
67373f7197eSJed 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 \
67473f7197eSJed Brown          -ksp_monitor -snes_monitor -dm_view -dm_mat_type aijkokkos -dm_vec_type kokkos
675c4762a1bSJed Brown 
676c4762a1bSJed Brown TEST*/
677