xref: /petsc/src/snes/tutorials/ex13.c (revision b724ec41a692ea93fcd25473f3bd97088e04f9b6)
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 */
17*b724ec41SToby Isaac   PetscBool homogeneous;
18c4762a1bSJed Brown } AppCtx;
19c4762a1bSJed Brown 
20c4762a1bSJed Brown static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
21c4762a1bSJed Brown {
22c4762a1bSJed Brown   *u = 0.0;
23c4762a1bSJed Brown   return 0;
24c4762a1bSJed Brown }
25c4762a1bSJed Brown 
26*b724ec41SToby Isaac static PetscErrorCode trig_inhomogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
27c4762a1bSJed Brown {
28c4762a1bSJed Brown   PetscInt d;
29c4762a1bSJed Brown   *u = 0.0;
30c4762a1bSJed Brown   for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0*PETSC_PI*x[d]);
31c4762a1bSJed Brown   return 0;
32c4762a1bSJed Brown }
33c4762a1bSJed Brown 
34*b724ec41SToby Isaac static PetscErrorCode trig_homogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
35*b724ec41SToby Isaac {
36*b724ec41SToby Isaac   PetscInt d;
37*b724ec41SToby Isaac   *u = 1.0;
38*b724ec41SToby Isaac   for (d = 0; d < dim; ++d) *u *= PetscSinReal(2.0*PETSC_PI*x[d]);
39*b724ec41SToby Isaac   return 0;
40*b724ec41SToby Isaac }
41*b724ec41SToby Isaac 
42c4762a1bSJed Brown /* Compute integral of (residual of solution)*(adjoint solution - projection of adjoint solution) */
43c4762a1bSJed Brown static void obj_error_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
44c4762a1bSJed Brown                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
45c4762a1bSJed Brown                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
46c4762a1bSJed Brown                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar obj[])
47c4762a1bSJed Brown {
48c4762a1bSJed Brown   obj[0] = a[aOff[0]]*(u[0] - a[aOff[1]]);
49c4762a1bSJed Brown }
50c4762a1bSJed Brown 
51*b724ec41SToby Isaac static void f0_trig_inhomogeneous_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
52c4762a1bSJed Brown                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
53c4762a1bSJed Brown                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
54c4762a1bSJed Brown                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
55c4762a1bSJed Brown {
56c4762a1bSJed Brown   PetscInt d;
57c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f0[0] += -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]);
58c4762a1bSJed Brown }
59c4762a1bSJed Brown 
60*b724ec41SToby Isaac static void f0_trig_homogeneous_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
61*b724ec41SToby Isaac                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
62*b724ec41SToby Isaac                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
63*b724ec41SToby Isaac                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
64*b724ec41SToby Isaac {
65*b724ec41SToby Isaac   PetscInt d;
66*b724ec41SToby Isaac   for (d = 0; d < dim; ++d) {
67*b724ec41SToby Isaac     PetscScalar v = 1.;
68*b724ec41SToby Isaac     for (PetscInt e = 0; e < dim; e++) {
69*b724ec41SToby Isaac       if (e == d) {
70*b724ec41SToby Isaac         v *= -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]);
71*b724ec41SToby Isaac       } else {
72*b724ec41SToby Isaac         v *= PetscSinReal(2.0*PETSC_PI*x[d]);
73*b724ec41SToby Isaac       }
74*b724ec41SToby Isaac     }
75*b724ec41SToby Isaac     f0[0] += v;
76*b724ec41SToby Isaac   }
77*b724ec41SToby Isaac }
78*b724ec41SToby Isaac 
79c4762a1bSJed Brown static void f0_unity_u(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, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
83c4762a1bSJed Brown {
84c4762a1bSJed Brown   f0[0] = 1.0;
85c4762a1bSJed Brown }
86c4762a1bSJed Brown 
87c4762a1bSJed Brown static void f0_identityaux_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
88c4762a1bSJed Brown                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
89c4762a1bSJed Brown                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
90c4762a1bSJed Brown                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
91c4762a1bSJed Brown {
92c4762a1bSJed Brown   f0[0] = a[0];
93c4762a1bSJed Brown }
94c4762a1bSJed Brown 
95c4762a1bSJed Brown static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
96c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
97c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
98c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
99c4762a1bSJed Brown {
100c4762a1bSJed Brown   PetscInt d;
101c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
102c4762a1bSJed Brown }
103c4762a1bSJed Brown 
104c4762a1bSJed Brown static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
105c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
106c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
107c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
108c4762a1bSJed Brown {
109c4762a1bSJed Brown   PetscInt d;
110c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
111c4762a1bSJed Brown }
112c4762a1bSJed Brown 
113c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
114c4762a1bSJed Brown {
115c4762a1bSJed Brown   PetscErrorCode ierr;
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   PetscFunctionBeginUser;
118c4762a1bSJed Brown   options->shear    = PETSC_FALSE;
119c4762a1bSJed Brown   options->spectral = PETSC_FALSE;
120c4762a1bSJed Brown   options->adjoint  = PETSC_FALSE;
121*b724ec41SToby Isaac   options->homogeneous = PETSC_FALSE;
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");CHKERRQ(ierr);
124c4762a1bSJed Brown   ierr = PetscOptionsBool("-shear", "Shear the domain", "ex13.c", options->shear, &options->shear, NULL);CHKERRQ(ierr);
125c4762a1bSJed Brown   ierr = PetscOptionsBool("-spectral", "Look at the spectrum along planes of the solution", "ex13.c", options->spectral, &options->spectral, NULL);CHKERRQ(ierr);
126c4762a1bSJed Brown   ierr = PetscOptionsBool("-adjoint", "Solve the adjoint problem", "ex13.c", options->adjoint, &options->adjoint, NULL);CHKERRQ(ierr);
127*b724ec41SToby Isaac   ierr = PetscOptionsBool("-homogeneous", "Use homogeneous boundary conditions", "ex13.c", options->homogeneous, &options->homogeneous, NULL);CHKERRQ(ierr);
128c4762a1bSJed Brown   ierr = PetscOptionsEnd();
129c4762a1bSJed Brown   PetscFunctionReturn(0);
130c4762a1bSJed Brown }
131c4762a1bSJed Brown 
132c4762a1bSJed Brown static PetscErrorCode CreateSpectralPlanes(DM dm, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user)
133c4762a1bSJed Brown {
134c4762a1bSJed Brown   PetscSection       coordSection;
135c4762a1bSJed Brown   Vec                coordinates;
136c4762a1bSJed Brown   const PetscScalar *coords;
137c4762a1bSJed Brown   PetscInt           dim, p, vStart, vEnd, v;
138c4762a1bSJed Brown   PetscErrorCode     ierr;
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   PetscFunctionBeginUser;
141c4762a1bSJed Brown   ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr);
142c4762a1bSJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
143c4762a1bSJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
144c4762a1bSJed Brown   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
145c4762a1bSJed Brown   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
146c4762a1bSJed Brown   for (p = 0; p < numPlanes; ++p) {
147c4762a1bSJed Brown     DMLabel label;
148c4762a1bSJed Brown     char    name[PETSC_MAX_PATH_LEN];
149c4762a1bSJed Brown 
150c4762a1bSJed Brown     ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%D", p);CHKERRQ(ierr);
151c4762a1bSJed Brown     ierr = DMCreateLabel(dm, name);CHKERRQ(ierr);
152c4762a1bSJed Brown     ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr);
153c4762a1bSJed Brown     ierr = DMLabelAddStratum(label, 1);CHKERRQ(ierr);
154c4762a1bSJed Brown     for (v = vStart; v < vEnd; ++v) {
155c4762a1bSJed Brown       PetscInt off;
156c4762a1bSJed Brown 
157c4762a1bSJed Brown       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
158c4762a1bSJed Brown       if (PetscAbsReal(planeCoord[p] - PetscRealPart(coords[off+planeDir[p]])) < PETSC_SMALL) {
159c4762a1bSJed Brown         ierr = DMLabelSetValue(label, v, 1);CHKERRQ(ierr);
160c4762a1bSJed Brown       }
161c4762a1bSJed Brown     }
162c4762a1bSJed Brown   }
163c4762a1bSJed Brown   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
164c4762a1bSJed Brown   PetscFunctionReturn(0);
165c4762a1bSJed Brown }
166c4762a1bSJed Brown 
167c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
168c4762a1bSJed Brown {
169c4762a1bSJed Brown   PetscErrorCode ierr;
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   PetscFunctionBeginUser;
172c4762a1bSJed Brown   /* Create box mesh */
1735be41b18SMatthew G. Knepley   ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
174c4762a1bSJed Brown   /* TODO: This should be pulled into the library */
175c4762a1bSJed Brown   {
176c4762a1bSJed Brown     char      convType[256];
177c4762a1bSJed Brown     PetscBool flg;
178c4762a1bSJed Brown 
179c4762a1bSJed Brown     ierr = PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");CHKERRQ(ierr);
180c4762a1bSJed Brown     ierr = PetscOptionsFList("-dm_plex_convert_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);CHKERRQ(ierr);
181c4762a1bSJed Brown     ierr = PetscOptionsEnd();
182c4762a1bSJed Brown     if (flg) {
183c4762a1bSJed Brown       DM dmConv;
184c4762a1bSJed Brown 
185c4762a1bSJed Brown       ierr = DMConvert(*dm,convType,&dmConv);CHKERRQ(ierr);
186c4762a1bSJed Brown       if (dmConv) {
187c4762a1bSJed Brown         ierr = DMDestroy(dm);CHKERRQ(ierr);
188c4762a1bSJed Brown         *dm  = dmConv;
189c4762a1bSJed Brown       }
190c4762a1bSJed Brown     }
191c4762a1bSJed Brown   }
192c4762a1bSJed Brown   if (user->shear) {ierr = DMPlexShearGeometry(*dm, DM_X, NULL);CHKERRQ(ierr);}
193c4762a1bSJed Brown   ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr);
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
196c4762a1bSJed Brown   ierr = DMSetApplicationContext(*dm, user);CHKERRQ(ierr);
197c4762a1bSJed Brown   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
198c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
199c4762a1bSJed Brown   if (user->spectral) {
200c4762a1bSJed Brown     PetscInt  planeDir[2]   = {0,  1};
201c4762a1bSJed Brown     PetscReal planeCoord[2] = {0., 1.};
202c4762a1bSJed Brown 
203c4762a1bSJed Brown     ierr = CreateSpectralPlanes(*dm, 2, planeDir, planeCoord, user);CHKERRQ(ierr);
204c4762a1bSJed Brown   }
205c4762a1bSJed Brown   PetscFunctionReturn(0);
206c4762a1bSJed Brown }
207c4762a1bSJed Brown 
208c4762a1bSJed Brown static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
209c4762a1bSJed Brown {
210c4762a1bSJed Brown   PetscDS        prob;
211c4762a1bSJed Brown   const PetscInt id = 1;
212*b724ec41SToby Isaac   PetscPointFunc f0 = user->homogeneous ? f0_trig_homogeneous_u : f0_trig_inhomogeneous_u;
213*b724ec41SToby Isaac   PetscErrorCode (*ex)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *) = user->homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u;
214c4762a1bSJed Brown   PetscErrorCode ierr;
215c4762a1bSJed Brown 
216c4762a1bSJed Brown   PetscFunctionBeginUser;
217c4762a1bSJed Brown   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
218*b724ec41SToby Isaac   ierr = PetscDSSetResidual(prob, 0, f0, f1_u);CHKERRQ(ierr);
219*b724ec41SToby Isaac   ierr = PetscDSSetExactSolution(prob, 0, ex, user);CHKERRQ(ierr);
220*b724ec41SToby Isaac   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) ex, NULL, 1, &id, user);CHKERRQ(ierr);
221c4762a1bSJed Brown   ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr);
222c4762a1bSJed Brown   PetscFunctionReturn(0);
223c4762a1bSJed Brown }
224c4762a1bSJed Brown 
225c4762a1bSJed Brown static PetscErrorCode SetupAdjointProblem(DM dm, AppCtx *user)
226c4762a1bSJed Brown {
227c4762a1bSJed Brown   PetscDS        prob;
228c4762a1bSJed Brown   const PetscInt id = 1;
229c4762a1bSJed Brown   PetscErrorCode ierr;
230c4762a1bSJed Brown 
231c4762a1bSJed Brown   PetscFunctionBeginUser;
232c4762a1bSJed Brown   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
233c4762a1bSJed Brown   ierr = PetscDSSetResidual(prob, 0, f0_unity_u, f1_u);CHKERRQ(ierr);
234c4762a1bSJed Brown   ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr);
235c4762a1bSJed Brown   ierr = PetscDSSetObjective(prob, 0, obj_error_u);CHKERRQ(ierr);
23656cf3b9cSMatthew G. Knepley   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) zero, NULL, 1, &id, user);CHKERRQ(ierr);
237c4762a1bSJed Brown   PetscFunctionReturn(0);
238c4762a1bSJed Brown }
239c4762a1bSJed Brown 
240c4762a1bSJed Brown static PetscErrorCode SetupErrorProblem(DM dm, AppCtx *user)
241c4762a1bSJed Brown {
242c4762a1bSJed Brown   PetscDS        prob;
243c4762a1bSJed Brown   PetscErrorCode ierr;
244c4762a1bSJed Brown 
245c4762a1bSJed Brown   PetscFunctionBeginUser;
246c4762a1bSJed Brown   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
247c4762a1bSJed Brown   PetscFunctionReturn(0);
248c4762a1bSJed Brown }
249c4762a1bSJed Brown 
250c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
251c4762a1bSJed Brown {
252c4762a1bSJed Brown   DM             cdm = dm;
253c4762a1bSJed Brown   PetscFE        fe;
2545be41b18SMatthew G. Knepley   DMPolytopeType ct;
2555be41b18SMatthew G. Knepley   PetscBool      simplex;
2565be41b18SMatthew G. Knepley   PetscInt       dim, cStart;
257c4762a1bSJed Brown   char           prefix[PETSC_MAX_PATH_LEN];
258c4762a1bSJed Brown   PetscErrorCode ierr;
259c4762a1bSJed Brown 
260c4762a1bSJed Brown   PetscFunctionBeginUser;
2615be41b18SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2625be41b18SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr);
2635be41b18SMatthew G. Knepley   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
2645be41b18SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
265c4762a1bSJed Brown   /* Create finite element */
266c4762a1bSJed Brown   ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name);CHKERRQ(ierr);
267260ffc31SMatthew G. Knepley   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe);CHKERRQ(ierr);
268c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr);
269c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
270c4762a1bSJed Brown   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr);
271c4762a1bSJed Brown   ierr = DMCreateDS(dm);CHKERRQ(ierr);
272c4762a1bSJed Brown   ierr = (*setup)(dm, user);CHKERRQ(ierr);
273c4762a1bSJed Brown   while (cdm) {
274c4762a1bSJed Brown     ierr = DMCopyDisc(dm,cdm);CHKERRQ(ierr);
275c4762a1bSJed Brown     /* TODO: Check whether the boundary of coarse meshes is marked */
276c4762a1bSJed Brown     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
277c4762a1bSJed Brown   }
278c4762a1bSJed Brown   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
279c4762a1bSJed Brown   PetscFunctionReturn(0);
280c4762a1bSJed Brown }
281c4762a1bSJed Brown 
282c4762a1bSJed Brown static PetscErrorCode ComputeSpectral(DM dm, Vec u, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user)
283c4762a1bSJed Brown {
284c4762a1bSJed Brown   MPI_Comm           comm;
285c4762a1bSJed Brown   PetscSection       coordSection, section;
286c4762a1bSJed Brown   Vec                coordinates, uloc;
287c4762a1bSJed Brown   const PetscScalar *coords, *array;
288c4762a1bSJed Brown   PetscInt           p;
289c4762a1bSJed Brown   PetscMPIInt        size, rank;
290c4762a1bSJed Brown   PetscErrorCode     ierr;
291c4762a1bSJed Brown 
292c4762a1bSJed Brown   PetscFunctionBeginUser;
293c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
294ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
295ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
296c4762a1bSJed Brown   ierr = DMGetLocalVector(dm, &uloc);CHKERRQ(ierr);
297c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(dm, u, INSERT_VALUES, uloc);CHKERRQ(ierr);
298c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(dm, u, INSERT_VALUES, uloc);CHKERRQ(ierr);
299c4762a1bSJed Brown   ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, uloc, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
300c4762a1bSJed Brown   ierr = VecViewFromOptions(uloc, NULL, "-sol_view");CHKERRQ(ierr);
301c4762a1bSJed Brown   ierr = DMGetLocalSection(dm, &section);CHKERRQ(ierr);
302c4762a1bSJed Brown   ierr = VecGetArrayRead(uloc, &array);CHKERRQ(ierr);
303c4762a1bSJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
304c4762a1bSJed Brown   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
305c4762a1bSJed Brown   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
306c4762a1bSJed Brown   for (p = 0; p < numPlanes; ++p) {
307c4762a1bSJed Brown     DMLabel         label;
308c4762a1bSJed Brown     char            name[PETSC_MAX_PATH_LEN];
309c4762a1bSJed Brown     Mat             F;
310c4762a1bSJed Brown     Vec             x, y;
311c4762a1bSJed Brown     IS              stratum;
312c4762a1bSJed Brown     PetscReal      *ray, *gray;
313c4762a1bSJed Brown     PetscScalar    *rvals, *svals, *gsvals;
314c4762a1bSJed Brown     PetscInt       *perm, *nperm;
315c4762a1bSJed Brown     PetscInt        n, N, i, j, off, offu;
316c4762a1bSJed Brown     const PetscInt *points;
317c4762a1bSJed Brown 
318c4762a1bSJed Brown     ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%D", p);CHKERRQ(ierr);
319c4762a1bSJed Brown     ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr);
320c4762a1bSJed Brown     ierr = DMLabelGetStratumIS(label, 1, &stratum);CHKERRQ(ierr);
321c4762a1bSJed Brown     ierr = ISGetLocalSize(stratum, &n);CHKERRQ(ierr);
322c4762a1bSJed Brown     ierr = ISGetIndices(stratum, &points);CHKERRQ(ierr);
323c4762a1bSJed Brown     ierr = PetscMalloc2(n, &ray, n, &svals);CHKERRQ(ierr);
324c4762a1bSJed Brown     for (i = 0; i < n; ++i) {
325c4762a1bSJed Brown       ierr = PetscSectionGetOffset(coordSection, points[i], &off);CHKERRQ(ierr);
326c4762a1bSJed Brown       ierr = PetscSectionGetOffset(section, points[i], &offu);CHKERRQ(ierr);
327c4762a1bSJed Brown       ray[i]   = PetscRealPart(coords[off+((planeDir[p]+1)%2)]);
328c4762a1bSJed Brown       svals[i] = array[offu];
329c4762a1bSJed Brown     }
330c4762a1bSJed Brown     /* Gather the ray data to proc 0 */
331c4762a1bSJed Brown     if (size > 1) {
332c4762a1bSJed Brown       PetscMPIInt *cnt, *displs, p;
333c4762a1bSJed Brown 
334c4762a1bSJed Brown       ierr = PetscCalloc2(size, &cnt, size, &displs);CHKERRQ(ierr);
335ffc4695bSBarry Smith       ierr = MPI_Gather(&n, 1, MPIU_INT, cnt, 1, MPIU_INT, 0, comm);CHKERRMPI(ierr);
336c4762a1bSJed Brown       for (p = 1; p < size; ++p) displs[p] = displs[p-1] + cnt[p-1];
337c4762a1bSJed Brown       N = displs[size-1] + cnt[size-1];
338c4762a1bSJed Brown       ierr = PetscMalloc2(N, &gray, N, &gsvals);CHKERRQ(ierr);
339ffc4695bSBarry Smith       ierr = MPI_Gatherv(ray, n, MPIU_REAL, gray, cnt, displs, MPIU_REAL, 0, comm);CHKERRMPI(ierr);
340ffc4695bSBarry Smith       ierr = MPI_Gatherv(svals, n, MPIU_SCALAR, gsvals, cnt, displs, MPIU_SCALAR, 0, comm);CHKERRMPI(ierr);
341c4762a1bSJed Brown       ierr = PetscFree2(cnt, displs);CHKERRQ(ierr);
342c4762a1bSJed Brown     } else {
343c4762a1bSJed Brown       N      = n;
344c4762a1bSJed Brown       gray   = ray;
345c4762a1bSJed Brown       gsvals = svals;
346c4762a1bSJed Brown     }
347c4762a1bSJed Brown     if (!rank) {
348c4762a1bSJed Brown       /* Sort point along ray */
349c4762a1bSJed Brown       ierr = PetscMalloc2(N, &perm, N, &nperm);CHKERRQ(ierr);
350c4762a1bSJed Brown       for (i = 0; i < N; ++i) {perm[i] = i;}
351c4762a1bSJed Brown       ierr = PetscSortRealWithPermutation(N, gray, perm);CHKERRQ(ierr);
352c4762a1bSJed Brown       /* Count duplicates and squish mapping */
353c4762a1bSJed Brown       nperm[0] = perm[0];
354c4762a1bSJed Brown       for (i = 1, j = 1; i < N; ++i) {
355c4762a1bSJed Brown         if (PetscAbsReal(gray[perm[i]] - gray[perm[i-1]]) > PETSC_SMALL) nperm[j++] = perm[i];
356c4762a1bSJed Brown       }
357c4762a1bSJed Brown       /* Create FFT structs */
358c4762a1bSJed Brown       ierr = MatCreateFFT(PETSC_COMM_SELF, 1, &j, MATFFTW, &F);CHKERRQ(ierr);
359c4762a1bSJed Brown       ierr = MatCreateVecs(F, &x, &y);CHKERRQ(ierr);
360c4762a1bSJed Brown       ierr = PetscObjectSetName((PetscObject) y, name);CHKERRQ(ierr);
361c4762a1bSJed Brown       ierr = VecGetArray(x, &rvals);CHKERRQ(ierr);
362c4762a1bSJed Brown       for (i = 0, j = 0; i < N; ++i) {
363c4762a1bSJed Brown         if (i > 0 && PetscAbsReal(gray[perm[i]] - gray[perm[i-1]]) < PETSC_SMALL) continue;
364c4762a1bSJed Brown         rvals[j] = gsvals[nperm[j]];
365c4762a1bSJed Brown         ++j;
366c4762a1bSJed Brown       }
367c4762a1bSJed Brown       ierr = PetscFree2(perm, nperm);CHKERRQ(ierr);
368c4762a1bSJed Brown       if (size > 1) {ierr = PetscFree2(gray, gsvals);CHKERRQ(ierr);}
369c4762a1bSJed Brown       ierr = VecRestoreArray(x, &rvals);CHKERRQ(ierr);
370c4762a1bSJed Brown       /* Do FFT along the ray */
371c4762a1bSJed Brown       ierr = MatMult(F, x, y);CHKERRQ(ierr);
372c4762a1bSJed Brown       /* Chop FFT */
373c4762a1bSJed Brown       ierr = VecChop(y, PETSC_SMALL);CHKERRQ(ierr);
374c4762a1bSJed Brown       ierr = VecViewFromOptions(x, NULL, "-real_view");CHKERRQ(ierr);
375c4762a1bSJed Brown       ierr = VecViewFromOptions(y, NULL, "-fft_view");CHKERRQ(ierr);
376c4762a1bSJed Brown       ierr = VecDestroy(&x);CHKERRQ(ierr);
377c4762a1bSJed Brown       ierr = VecDestroy(&y);CHKERRQ(ierr);
378c4762a1bSJed Brown       ierr = MatDestroy(&F);CHKERRQ(ierr);
379c4762a1bSJed Brown     }
380c4762a1bSJed Brown     ierr = ISRestoreIndices(stratum, &points);CHKERRQ(ierr);
381c4762a1bSJed Brown     ierr = ISDestroy(&stratum);CHKERRQ(ierr);
382c4762a1bSJed Brown     ierr = PetscFree2(ray, svals);CHKERRQ(ierr);
383c4762a1bSJed Brown   }
384c4762a1bSJed Brown   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
385c4762a1bSJed Brown   ierr = VecRestoreArrayRead(uloc, &array);CHKERRQ(ierr);
386c4762a1bSJed Brown   ierr = DMRestoreLocalVector(dm, &uloc);CHKERRQ(ierr);
387c4762a1bSJed Brown   PetscFunctionReturn(0);
388c4762a1bSJed Brown }
389c4762a1bSJed Brown 
390c4762a1bSJed Brown int main(int argc, char **argv)
391c4762a1bSJed Brown {
392c4762a1bSJed Brown   DM             dm;   /* Problem specification */
393c4762a1bSJed Brown   SNES           snes; /* Nonlinear solver */
394c4762a1bSJed Brown   Vec            u;    /* Solutions */
395c4762a1bSJed Brown   AppCtx         user; /* User-defined work context */
396c4762a1bSJed Brown   PetscErrorCode ierr;
397c4762a1bSJed Brown 
398c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
399c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
400c4762a1bSJed Brown   /* Primal system */
401c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
402c4762a1bSJed Brown   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
403c4762a1bSJed Brown   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
404c4762a1bSJed Brown   ierr = SetupDiscretization(dm, "potential", SetupPrimalProblem, &user);CHKERRQ(ierr);
405c4762a1bSJed Brown   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
406c4762a1bSJed Brown   ierr = VecSet(u, 0.0);CHKERRQ(ierr);
407c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) u, "potential");CHKERRQ(ierr);
408c4762a1bSJed Brown   ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr);
409c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
410c4762a1bSJed Brown   ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
411c4762a1bSJed Brown   ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
412c4762a1bSJed Brown   ierr = VecViewFromOptions(u, NULL, "-potential_view");CHKERRQ(ierr);
413c4762a1bSJed Brown   if (user.spectral) {
414c4762a1bSJed Brown     PetscInt  planeDir[2]   = {0,  1};
415c4762a1bSJed Brown     PetscReal planeCoord[2] = {0., 1.};
416c4762a1bSJed Brown 
417c4762a1bSJed Brown     ierr = ComputeSpectral(dm, u, 2, planeDir, planeCoord, &user);CHKERRQ(ierr);
418c4762a1bSJed Brown   }
419c4762a1bSJed Brown   /* Adjoint system */
420c4762a1bSJed Brown   if (user.adjoint) {
421c4762a1bSJed Brown     DM   dmAdj;
422c4762a1bSJed Brown     SNES snesAdj;
423c4762a1bSJed Brown     Vec  uAdj;
424c4762a1bSJed Brown 
425c4762a1bSJed Brown     ierr = SNESCreate(PETSC_COMM_WORLD, &snesAdj);CHKERRQ(ierr);
426c4762a1bSJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject) snesAdj, "adjoint_");CHKERRQ(ierr);
427c4762a1bSJed Brown     ierr = DMClone(dm, &dmAdj);CHKERRQ(ierr);
428c4762a1bSJed Brown     ierr = SNESSetDM(snesAdj, dmAdj);CHKERRQ(ierr);
429c4762a1bSJed Brown     ierr = SetupDiscretization(dmAdj, "adjoint", SetupAdjointProblem, &user);CHKERRQ(ierr);
430c4762a1bSJed Brown     ierr = DMCreateGlobalVector(dmAdj, &uAdj);CHKERRQ(ierr);
431c4762a1bSJed Brown     ierr = VecSet(uAdj, 0.0);CHKERRQ(ierr);
432c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject) uAdj, "adjoint");CHKERRQ(ierr);
433c4762a1bSJed Brown     ierr = DMPlexSetSNESLocalFEM(dmAdj, &user, &user, &user);CHKERRQ(ierr);
434c4762a1bSJed Brown     ierr = SNESSetFromOptions(snesAdj);CHKERRQ(ierr);
435c4762a1bSJed Brown     ierr = SNESSolve(snesAdj, NULL, uAdj);CHKERRQ(ierr);
436c4762a1bSJed Brown     ierr = SNESGetSolution(snesAdj, &uAdj);CHKERRQ(ierr);
437c4762a1bSJed Brown     ierr = VecViewFromOptions(uAdj, NULL, "-adjoint_view");CHKERRQ(ierr);
438c4762a1bSJed Brown     /* Error representation */
439c4762a1bSJed Brown     {
440c4762a1bSJed Brown       DM        dmErr, dmErrAux, dms[2];
441c4762a1bSJed Brown       Vec       errorEst, errorL2, uErr, uErrLoc, uAdjLoc, uAdjProj;
442c4762a1bSJed Brown       IS       *subis;
443c4762a1bSJed Brown       PetscReal errorEstTot, errorL2Norm, errorL2Tot;
444c4762a1bSJed Brown       PetscInt  N, i;
445*b724ec41SToby Isaac       PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *) = {user.homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u};
446c4762a1bSJed Brown       void (*identity[1])(PetscInt, PetscInt, PetscInt,
447c4762a1bSJed Brown                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
448c4762a1bSJed Brown                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
449c4762a1bSJed Brown                           PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]) = {f0_identityaux_u};
450c4762a1bSJed Brown       void            *ctxs[1] = {0};
451c4762a1bSJed Brown 
452c4762a1bSJed Brown       ctxs[0] = &user;
453c4762a1bSJed Brown       ierr = DMClone(dm, &dmErr);CHKERRQ(ierr);
454c4762a1bSJed Brown       ierr = SetupDiscretization(dmErr, "error", SetupErrorProblem, &user);CHKERRQ(ierr);
455c4762a1bSJed Brown       ierr = DMGetGlobalVector(dmErr, &errorEst);CHKERRQ(ierr);
456c4762a1bSJed Brown       ierr = DMGetGlobalVector(dmErr, &errorL2);CHKERRQ(ierr);
457c4762a1bSJed Brown       /*   Compute auxiliary data (solution and projection of adjoint solution) */
458c4762a1bSJed Brown       ierr = DMGetLocalVector(dmAdj, &uAdjLoc);CHKERRQ(ierr);
459c4762a1bSJed Brown       ierr = DMGlobalToLocalBegin(dmAdj, uAdj, INSERT_VALUES, uAdjLoc);CHKERRQ(ierr);
460c4762a1bSJed Brown       ierr = DMGlobalToLocalEnd(dmAdj, uAdj, INSERT_VALUES, uAdjLoc);CHKERRQ(ierr);
461c4762a1bSJed Brown       ierr = DMGetGlobalVector(dm, &uAdjProj);CHKERRQ(ierr);
462c4762a1bSJed Brown       ierr = PetscObjectCompose((PetscObject) dm, "dmAux", (PetscObject) dmAdj);CHKERRQ(ierr);
463c4762a1bSJed Brown       ierr = PetscObjectCompose((PetscObject) dm, "A", (PetscObject) uAdjLoc);CHKERRQ(ierr);
464c4762a1bSJed Brown       ierr = DMProjectField(dm, 0.0, u, identity, INSERT_VALUES, uAdjProj);CHKERRQ(ierr);
465c4762a1bSJed Brown       ierr = PetscObjectCompose((PetscObject) dm, "dmAux", NULL);CHKERRQ(ierr);
466c4762a1bSJed Brown       ierr = PetscObjectCompose((PetscObject) dm, "A", NULL);CHKERRQ(ierr);
467c4762a1bSJed Brown       ierr = DMRestoreLocalVector(dmAdj, &uAdjLoc);CHKERRQ(ierr);
468c4762a1bSJed Brown       /*   Attach auxiliary data */
469c4762a1bSJed Brown       dms[0] = dm; dms[1] = dm;
470c4762a1bSJed Brown       ierr = DMCreateSuperDM(dms, 2, &subis, &dmErrAux);CHKERRQ(ierr);
471c4762a1bSJed Brown       if (0) {
472c4762a1bSJed Brown         PetscSection sec;
473c4762a1bSJed Brown 
474c4762a1bSJed Brown         ierr = DMGetLocalSection(dms[0], &sec);CHKERRQ(ierr);
475c4762a1bSJed Brown         ierr = PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
476c4762a1bSJed Brown         ierr = DMGetLocalSection(dms[1], &sec);CHKERRQ(ierr);
477c4762a1bSJed Brown         ierr = PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
478c4762a1bSJed Brown         ierr = DMGetLocalSection(dmErrAux, &sec);CHKERRQ(ierr);
479c4762a1bSJed Brown         ierr = PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
480c4762a1bSJed Brown       }
481c4762a1bSJed Brown       ierr = DMViewFromOptions(dmErrAux, NULL, "-dm_err_view");CHKERRQ(ierr);
482c4762a1bSJed Brown       ierr = ISViewFromOptions(subis[0], NULL, "-super_is_view");CHKERRQ(ierr);
483c4762a1bSJed Brown       ierr = ISViewFromOptions(subis[1], NULL, "-super_is_view");CHKERRQ(ierr);
484c4762a1bSJed Brown       ierr = DMGetGlobalVector(dmErrAux, &uErr);CHKERRQ(ierr);
485c4762a1bSJed Brown       ierr = VecViewFromOptions(u, NULL, "-map_vec_view");CHKERRQ(ierr);
486c4762a1bSJed Brown       ierr = VecViewFromOptions(uAdjProj, NULL, "-map_vec_view");CHKERRQ(ierr);
487c4762a1bSJed Brown       ierr = VecViewFromOptions(uErr, NULL, "-map_vec_view");CHKERRQ(ierr);
488c4762a1bSJed Brown       ierr = VecISCopy(uErr, subis[0], SCATTER_FORWARD, u);CHKERRQ(ierr);
489c4762a1bSJed Brown       ierr = VecISCopy(uErr, subis[1], SCATTER_FORWARD, uAdjProj);CHKERRQ(ierr);
490c4762a1bSJed Brown       ierr = DMRestoreGlobalVector(dm, &uAdjProj);CHKERRQ(ierr);
491c4762a1bSJed Brown       for (i = 0; i < 2; ++i) {ierr = ISDestroy(&subis[i]);CHKERRQ(ierr);}
492c4762a1bSJed Brown       ierr = PetscFree(subis);CHKERRQ(ierr);
493c4762a1bSJed Brown       ierr = DMGetLocalVector(dmErrAux, &uErrLoc);CHKERRQ(ierr);
494c4762a1bSJed Brown       ierr = DMGlobalToLocalBegin(dm, uErr, INSERT_VALUES, uErrLoc);CHKERRQ(ierr);
495c4762a1bSJed Brown       ierr = DMGlobalToLocalEnd(dm, uErr, INSERT_VALUES, uErrLoc);CHKERRQ(ierr);
496c4762a1bSJed Brown       ierr = DMRestoreGlobalVector(dmErrAux, &uErr);CHKERRQ(ierr);
497c4762a1bSJed Brown       ierr = PetscObjectCompose((PetscObject) dmAdj, "dmAux", (PetscObject) dmErrAux);CHKERRQ(ierr);
498c4762a1bSJed Brown       ierr = PetscObjectCompose((PetscObject) dmAdj, "A", (PetscObject) uErrLoc);CHKERRQ(ierr);
499c4762a1bSJed Brown       /*   Compute cellwise error estimate */
500c4762a1bSJed Brown       ierr = VecSet(errorEst, 0.0);CHKERRQ(ierr);
501c4762a1bSJed Brown       ierr = DMPlexComputeCellwiseIntegralFEM(dmAdj, uAdj, errorEst, &user);CHKERRQ(ierr);
502c4762a1bSJed Brown       ierr = PetscObjectCompose((PetscObject) dmAdj, "dmAux", NULL);CHKERRQ(ierr);
503c4762a1bSJed Brown       ierr = PetscObjectCompose((PetscObject) dmAdj, "A", NULL);CHKERRQ(ierr);
504c4762a1bSJed Brown       ierr = DMRestoreLocalVector(dmErrAux, &uErrLoc);CHKERRQ(ierr);
505c4762a1bSJed Brown       ierr = DMDestroy(&dmErrAux);CHKERRQ(ierr);
506c4762a1bSJed Brown       /*   Plot cellwise error vector */
507c4762a1bSJed Brown       ierr = VecViewFromOptions(errorEst, NULL, "-error_view");CHKERRQ(ierr);
508c4762a1bSJed Brown       /*   Compute ratio of estimate (sum over cells) with actual L_2 error */
509c4762a1bSJed Brown       ierr = DMComputeL2Diff(dm, 0.0, funcs, ctxs, u, &errorL2Norm);CHKERRQ(ierr);
510c4762a1bSJed Brown       ierr = DMPlexComputeL2DiffVec(dm, 0.0, funcs, ctxs, u, errorL2);CHKERRQ(ierr);
511c4762a1bSJed Brown       ierr = VecViewFromOptions(errorL2, NULL, "-l2_error_view");CHKERRQ(ierr);
512c4762a1bSJed Brown       ierr = VecNorm(errorL2,  NORM_INFINITY, &errorL2Tot);CHKERRQ(ierr);
513c4762a1bSJed Brown       ierr = VecNorm(errorEst, NORM_INFINITY, &errorEstTot);CHKERRQ(ierr);
514c4762a1bSJed Brown       ierr = VecGetSize(errorEst, &N);CHKERRQ(ierr);
515c4762a1bSJed Brown       ierr = VecPointwiseDivide(errorEst, errorEst, errorL2);CHKERRQ(ierr);
516c4762a1bSJed Brown       ierr = PetscObjectSetName((PetscObject) errorEst, "Error ratio");CHKERRQ(ierr);
517c4762a1bSJed Brown       ierr = VecViewFromOptions(errorEst, NULL, "-error_ratio_view");CHKERRQ(ierr);
518c4762a1bSJed 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);
519c4762a1bSJed Brown       ierr = DMRestoreGlobalVector(dmErr, &errorEst);CHKERRQ(ierr);
520c4762a1bSJed Brown       ierr = DMRestoreGlobalVector(dmErr, &errorL2);CHKERRQ(ierr);
521c4762a1bSJed Brown       ierr = DMDestroy(&dmErr);CHKERRQ(ierr);
522c4762a1bSJed Brown     }
523c4762a1bSJed Brown     ierr = DMDestroy(&dmAdj);CHKERRQ(ierr);
524c4762a1bSJed Brown     ierr = VecDestroy(&uAdj);CHKERRQ(ierr);
525c4762a1bSJed Brown     ierr = SNESDestroy(&snesAdj);CHKERRQ(ierr);
526c4762a1bSJed Brown   }
527c4762a1bSJed Brown   /* Cleanup */
528c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
529c4762a1bSJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
530c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
531c4762a1bSJed Brown   ierr = PetscFinalize();
532c4762a1bSJed Brown   return ierr;
533c4762a1bSJed Brown }
534c4762a1bSJed Brown 
535c4762a1bSJed Brown /*TEST
536c4762a1bSJed Brown 
537c4762a1bSJed Brown   test:
5385be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9
5395be41b18SMatthew G. Knepley     suffix: 2d_p1_conv
540c4762a1bSJed Brown     requires: triangle
5415be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2
5425be41b18SMatthew G. Knepley   test:
5435be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9
5445be41b18SMatthew G. Knepley     suffix: 2d_p2_conv
5455be41b18SMatthew G. Knepley     requires: triangle
5465be41b18SMatthew G. Knepley     args: -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_p3_conv
5505be41b18SMatthew G. Knepley     requires: triangle
5515be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2
5525be41b18SMatthew G. Knepley   test:
5535be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9
5545be41b18SMatthew G. Knepley     suffix: 2d_q1_conv
5555be41b18SMatthew G. Knepley     args: -dm_plex_box_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2
5565be41b18SMatthew G. Knepley   test:
5575be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9
5585be41b18SMatthew G. Knepley     suffix: 2d_q2_conv
5595be41b18SMatthew G. Knepley     args: -dm_plex_box_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2
5605be41b18SMatthew G. Knepley   test:
5615be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9
5625be41b18SMatthew G. Knepley     suffix: 2d_q3_conv
5635be41b18SMatthew G. Knepley     args: -dm_plex_box_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2
5645be41b18SMatthew G. Knepley   test:
5655be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9
5665be41b18SMatthew G. Knepley     suffix: 2d_q1_shear_conv
5675be41b18SMatthew G. Knepley     args: -dm_plex_box_simplex 0 -shear -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2
5685be41b18SMatthew G. Knepley   test:
5695be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9
5705be41b18SMatthew G. Knepley     suffix: 2d_q2_shear_conv
5715be41b18SMatthew G. Knepley     args: -dm_plex_box_simplex 0 -shear -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2
5725be41b18SMatthew G. Knepley   test:
5735be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9
5745be41b18SMatthew G. Knepley     suffix: 2d_q3_shear_conv
5755be41b18SMatthew G. Knepley     args: -dm_plex_box_simplex 0 -shear -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2
5765be41b18SMatthew G. Knepley   test:
5770fdc7489SMatthew Knepley     # Using -convest_num_refine 3 we get L_2 convergence rate: 1.7
5785be41b18SMatthew G. Knepley     suffix: 3d_p1_conv
5795be41b18SMatthew G. Knepley     requires: ctetgen
5800fdc7489SMatthew Knepley     args: -dm_plex_box_dim 3 -dm_refine 1 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1
5815be41b18SMatthew G. Knepley   test:
5825be41b18SMatthew G. Knepley     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 2.8
5835be41b18SMatthew G. Knepley     suffix: 3d_p2_conv
5845be41b18SMatthew G. Knepley     requires: ctetgen
5855be41b18SMatthew G. Knepley     args: -dm_plex_box_dim 3 -dm_plex_box_faces 2,2,2 -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: 4.0
5885be41b18SMatthew G. Knepley     suffix: 3d_p3_conv
5895be41b18SMatthew G. Knepley     requires: ctetgen
5905be41b18SMatthew G. Knepley     args: -dm_plex_box_dim 3 -dm_plex_box_faces 2,2,2 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 1
5915be41b18SMatthew G. Knepley   test:
5925be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.8
5935be41b18SMatthew G. Knepley     suffix: 3d_q1_conv
5945be41b18SMatthew G. Knepley     args: -dm_plex_box_dim 3 -dm_plex_box_simplex 0 -dm_refine 1 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1
5955be41b18SMatthew G. Knepley   test:
5965be41b18SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.8
5975be41b18SMatthew G. Knepley     suffix: 3d_q2_conv
5985be41b18SMatthew G. Knepley     args: -dm_plex_box_dim 3 -dm_plex_box_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 1
5995be41b18SMatthew G. Knepley   test:
6005be41b18SMatthew G. Knepley     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 3.8
6015be41b18SMatthew G. Knepley     suffix: 3d_q3_conv
6025be41b18SMatthew G. Knepley     args: -dm_plex_box_dim 3 -dm_plex_box_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 1
6035be41b18SMatthew G. Knepley 
604c4762a1bSJed Brown   test:
605c4762a1bSJed Brown     suffix: 2d_p1_scalable
6065be41b18SMatthew G. Knepley     requires: triangle
6075be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_refine 3 \
608c4762a1bSJed Brown       -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned \
609c4762a1bSJed Brown       -pc_type gamg \
610c4762a1bSJed Brown         -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \
611c4762a1bSJed Brown         -pc_gamg_coarse_eq_limit 1000 \
612c4762a1bSJed Brown         -pc_gamg_square_graph 1 \
613c4762a1bSJed Brown         -pc_gamg_threshold 0.05 \
614c4762a1bSJed Brown         -pc_gamg_threshold_scale .0 \
615c4762a1bSJed Brown       -mg_levels_ksp_type chebyshev \
616c4762a1bSJed Brown         -mg_levels_ksp_max_it 1 \
617c4762a1bSJed Brown         -mg_levels_esteig_ksp_type cg \
618c4762a1bSJed Brown         -mg_levels_esteig_ksp_max_it 10 \
619c4762a1bSJed Brown         -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \
620c4762a1bSJed Brown         -mg_levels_pc_type jacobi \
621c4762a1bSJed Brown       -matptap_via scalable
622c4762a1bSJed Brown   test:
623c4762a1bSJed Brown     suffix: 2d_p1_gmg_vcycle
624c4762a1bSJed Brown     requires: triangle
6255be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
626c4762a1bSJed Brown           -ksp_rtol 5e-10 -pc_type mg \
627c4762a1bSJed Brown             -mg_levels_ksp_max_it 1 \
628c4762a1bSJed Brown             -mg_levels_esteig_ksp_type cg \
629c4762a1bSJed Brown             -mg_levels_esteig_ksp_max_it 10 \
630c4762a1bSJed Brown             -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \
631c4762a1bSJed Brown             -mg_levels_pc_type jacobi
632c4762a1bSJed Brown   test:
633c4762a1bSJed Brown     suffix: 2d_p1_gmg_fcycle
634c4762a1bSJed Brown     requires: triangle
6355be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
636c4762a1bSJed Brown           -ksp_rtol 5e-10 -pc_type mg -pc_mg_type full \
637c4762a1bSJed Brown             -mg_levels_ksp_max_it 2 \
638c4762a1bSJed Brown             -mg_levels_esteig_ksp_type cg \
639c4762a1bSJed Brown             -mg_levels_esteig_ksp_max_it 10 \
640c4762a1bSJed Brown             -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \
641c4762a1bSJed Brown             -mg_levels_pc_type jacobi
642c4762a1bSJed Brown   test:
643d6837840SMatthew G. Knepley     suffix: 2d_p1_gmg_vcycle_adapt
644d6837840SMatthew G. Knepley     requires: triangle bamg
6455be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \
646d6837840SMatthew 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 \
647d6837840SMatthew G. Knepley             -mg_levels_ksp_max_it 1 \
648d6837840SMatthew G. Knepley             -mg_levels_esteig_ksp_type cg \
649d6837840SMatthew G. Knepley             -mg_levels_esteig_ksp_max_it 10 \
650d6837840SMatthew G. Knepley             -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \
651d6837840SMatthew G. Knepley             -mg_levels_pc_type jacobi
652d6837840SMatthew G. Knepley   test:
653c4762a1bSJed Brown     suffix: 2d_p1_spectral_0
654c4762a1bSJed Brown     requires: triangle fftw !complex
6555be41b18SMatthew G. Knepley     args: -dm_plex_box_faces 1,1 -potential_petscspace_degree 1 -dm_refine 6 -spectral -fft_view
656c4762a1bSJed Brown   test:
657c4762a1bSJed Brown     suffix: 2d_p1_spectral_1
658c4762a1bSJed Brown     requires: triangle fftw !complex
659c4762a1bSJed Brown     nsize: 2
6605be41b18SMatthew G. Knepley     args: -dm_plex_box_faces 4,4 -dm_distribute -potential_petscspace_degree 1 -spectral -fft_view
661c4762a1bSJed Brown   test:
662c4762a1bSJed Brown     suffix: 2d_p1_adj_0
663c4762a1bSJed Brown     requires: triangle
6645be41b18SMatthew G. Knepley     args: -potential_petscspace_degree 1 -dm_refine 1 -adjoint -adjoint_petscspace_degree 1 -error_petscspace_degree 0
665b1ad28e3SJunchao Zhang   test:
666b1ad28e3SJunchao Zhang     nsize: 2
667b1ad28e3SJunchao Zhang     requires: kokkos_kernels
668b1ad28e3SJunchao Zhang     suffix: kokkos
669b1ad28e3SJunchao Zhang     args: -dm_plex_box_dim 3 -dm_plex_box_faces 2,3,6 -dm_distribute -petscpartitioner_type simple -dm_plex_box_simplex 0 -potential_petscspace_degree 1 \
670b1ad28e3SJunchao 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 \
671b1ad28e3SJunchao Zhang          -pc_gamg_threshold_scale .5 -mg_levels_ksp_type chebyshev -mg_levels_ksp_max_it 2 -mg_levels_esteig_ksp_type cg -mg_levels_esteig_ksp_max_it 10 \
672b1ad28e3SJunchao Zhang          -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 -mg_levels_pc_type jacobi -ksp_monitor -snes_monitor -dm_view -dm_mat_type aijkokkos -dm_vec_type kokkos
673c4762a1bSJed Brown 
674c4762a1bSJed Brown TEST*/
675