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