xref: /petsc/src/dm/impls/plex/tests/ex39.c (revision a69119a591a03a9d906b29c0a4e9802e4d7c9795)
1 const char help[] = "A test of H-div conforming discretizations on different cell types.\n";
2 
3 #include <petscdmplex.h>
4 #include <petscds.h>
5 #include <petscsnes.h>
6 #include <petscconvest.h>
7 #include <petscfe.h>
8 #include <petsc/private/petscfeimpl.h>
9 
10 /*
11   We are using the system
12 
13   \vec{u} = \vec{\hat{u}}
14   p = \div{\vec{u}} in low degree approximation space
15   d = \div{\vec{u}} - p == 0 in higher degree approximation space
16 
17   That is, we are using the field d to compute the error between \div{\vec{u}}
18   computed in a space 1 degree higher than p and the value of p which is
19   \div{u} computed in the low degree space. If H-div
20   elements are implemented correctly then this should be identically zero since
21   the divergence of a function in H(div) should be exactly representable in L_2
22   by definition.
23 */
24 static PetscErrorCode zero_func(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
25   PetscInt c;
26   for (c = 0; c < Nc; ++c) u[c] = 0;
27   return 0;
28 }
29 /* Linear Exact Functions
30    \vec{u} = \vec{x};
31    p = dim;
32    */
33 static PetscErrorCode linear_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
34   PetscInt c;
35   for (c = 0; c < Nc; ++c) u[c] = x[c];
36   return 0;
37 }
38 static PetscErrorCode linear_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
39   u[0] = dim;
40   return 0;
41 }
42 
43 /* Sinusoidal Exact Functions
44  * u_i = \sin{2*\pi*x_i}
45  * p = \Sum_{i=1}^{dim} 2*\pi*cos{2*\pi*x_i}
46  * */
47 
48 static PetscErrorCode sinusoid_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
49   PetscInt c;
50   for (c = 0; c < Nc; ++c) u[c] = PetscSinReal(2 * PETSC_PI * x[c]);
51   return 0;
52 }
53 static PetscErrorCode sinusoid_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
54   PetscInt d;
55   u[0] = 0;
56   for (d = 0; d < dim; ++d) u[0] += 2 * PETSC_PI * PetscCosReal(2 * PETSC_PI * x[d]);
57   return 0;
58 }
59 
60 /* Pointwise residual for u = u*. Need one of these for each possible u* */
61 static void f0_v_linear(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
62   PetscInt     i;
63   PetscScalar *u_rhs;
64 
65   PetscCalloc1(dim, &u_rhs);
66   (void)linear_u(dim, t, x, dim, u_rhs, NULL);
67   for (i = 0; i < dim; ++i) f0[i] = u[uOff[0] + i] - u_rhs[i];
68   PetscFree(u_rhs);
69 }
70 
71 static void f0_v_sinusoid(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
72   PetscInt     i;
73   PetscScalar *u_rhs;
74 
75   PetscCalloc1(dim, &u_rhs);
76   (void)sinusoid_u(dim, t, x, dim, u_rhs, NULL);
77   for (i = 0; i < dim; ++i) f0[i] = u[uOff[0] + i] - u_rhs[i];
78   PetscFree(u_rhs);
79 }
80 
81 /* Residual function for enforcing p = \div{u}. */
82 static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
83   PetscInt    i;
84   PetscScalar divu;
85 
86   divu = 0.;
87   for (i = 0; i < dim; ++i) divu += u_x[uOff_x[0] + i * dim + i];
88   f0[0] = u[uOff[1]] - divu;
89 }
90 
91 /* Residual function for p_err = \div{u} - p. */
92 static void f0_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
93   PetscInt    i;
94   PetscScalar divu;
95 
96   divu = 0.;
97   for (i = 0; i < dim; ++i) divu += u_x[uOff_x[0] + i * dim + i];
98   f0[0] = u[uOff[2]] - u[uOff[1]] + divu;
99 }
100 
101 /* Boundary residual for the embedding system. Need one for each form of
102  * solution. These enforce u = \hat{u} at the boundary. */
103 static void f0_bd_u_sinusoid(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
104   PetscInt     d;
105   PetscScalar *u_rhs;
106   PetscCalloc1(dim, &u_rhs);
107   (void)sinusoid_u(dim, t, x, dim, u_rhs, NULL);
108 
109   for (d = 0; d < dim; ++d) f0[d] = u_rhs[d];
110   PetscFree(u_rhs);
111 }
112 
113 static void f0_bd_u_linear(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
114   PetscInt     d;
115   PetscScalar *u_rhs;
116   PetscCalloc1(dim, &u_rhs);
117   (void)linear_u(dim, t, x, dim, u_rhs, NULL);
118 
119   for (d = 0; d < dim; ++d) f0[d] = u_rhs[d];
120   PetscFree(u_rhs);
121 }
122 /* Jacobian functions. For the following, v is the test function associated with
123  * u, q the test function associated with p, and w the test function associated
124  * with d. */
125 /* <v, u> */
126 static void g0_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) {
127   PetscInt c;
128   for (c = 0; c < dim; ++c) g0[c * dim + c] = 1.0;
129 }
130 
131 /* <q, p> */
132 static void g0_qp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) {
133   PetscInt d;
134   for (d = 0; d < dim; ++d) g0[d * dim + d] = 1.0;
135 }
136 
137 /* -<q, \div{u}> For the embedded system. This is different from the method of
138  * manufactured solution because instead of computing <q,\div{u}> - <q,f> we
139  * need <q,p> - <q,\div{u}.*/
140 static void g1_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) {
141   PetscInt d;
142   for (d = 0; d < dim; ++d) g1[d * dim + d] = -1.0;
143 }
144 
145 /* <w, p> This is only used by the embedded system. Where we need to compute
146  * <w,d> - <w,p> + <w, \div{u}>*/
147 static void g0_wp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) {
148   PetscInt d;
149   for (d = 0; d < dim; ++d) g0[d * dim + d] = -1.0;
150 }
151 
152 /* <w, d> */
153 static void g0_wd(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) {
154   PetscInt c;
155   for (c = 0; c < dim; ++c) g0[c * dim + c] = 1.0;
156 }
157 
158 /* <w, \div{u}> for the embedded system. */
159 static void g1_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) {
160   PetscInt d;
161   for (d = 0; d < dim; ++d) g1[d * dim + d] = 1.0;
162 }
163 
164 /* Enum and string array for selecting mesh perturbation options */
165 typedef enum {
166   NONE         = 0,
167   PERTURB      = 1,
168   SKEW         = 2,
169   SKEW_PERTURB = 3
170 } Transform;
171 const char *const TransformTypes[] = {"none", "perturb", "skew", "skew_perturb", "Perturbation", "", NULL};
172 
173 /* Enum and string array for selecting the form of the exact solution*/
174 typedef enum {
175   LINEAR     = 0,
176   SINUSOIDAL = 1
177 } Solution;
178 const char *const SolutionTypes[] = {"linear", "sinusoidal", "Solution", "", NULL};
179 
180 typedef struct {
181   Transform mesh_transform;
182   Solution  sol_form;
183 } UserCtx;
184 
185 /* Process command line options and initialize the UserCtx struct */
186 static PetscErrorCode ProcessOptions(MPI_Comm comm, UserCtx *user) {
187   PetscFunctionBegin;
188   /* Default to  2D, unperturbed triangle mesh and Linear solution.*/
189   user->mesh_transform = NONE;
190   user->sol_form       = LINEAR;
191 
192   PetscOptionsBegin(comm, "", "H-div Test Options", "DMPLEX");
193   PetscCall(PetscOptionsEnum("-mesh_transform", "Method used to perturb the mesh vertices. Options are skew, perturb, skew_perturb,or none", "ex39.c", TransformTypes, (PetscEnum)user->mesh_transform, (PetscEnum *)&user->mesh_transform, NULL));
194   PetscCall(PetscOptionsEnum("-sol_form", "Form of the exact solution. Options are Linear or Sinusoidal", "ex39.c", SolutionTypes, (PetscEnum)user->sol_form, (PetscEnum *)&user->sol_form, NULL));
195   PetscOptionsEnd();
196   PetscFunctionReturn(0);
197 }
198 
199 /* Perturb the position of each mesh vertex by a small amount.*/
200 static PetscErrorCode PerturbMesh(DM *mesh, PetscScalar *coordVals, PetscInt npoints, PetscInt dim) {
201   PetscInt    i, j, k;
202   PetscReal   minCoords[3], maxCoords[3], maxPert[3], randVal, amp;
203   PetscRandom ran;
204 
205   PetscFunctionBegin;
206   PetscCall(DMGetCoordinateDim(*mesh, &dim));
207   PetscCall(DMGetLocalBoundingBox(*mesh, minCoords, maxCoords));
208   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &ran));
209 
210   /* Compute something approximately equal to half an edge length. This is the
211    * most we can perturb points and guarantee that there won't be any topology
212    * issues. */
213   for (k = 0; k < dim; ++k) maxPert[k] = 0.025 * (maxCoords[k] - minCoords[k]) / (PetscPowReal(npoints, 1. / dim) - 1);
214   /* For each mesh vertex */
215   for (i = 0; i < npoints; ++i) {
216     /* For each coordinate of the vertex */
217     for (j = 0; j < dim; ++j) {
218       /* Generate a random amplitude in [-0.5*maxPert, 0.5*maxPert] */
219       PetscCall(PetscRandomGetValueReal(ran, &randVal));
220       amp = maxPert[j] * (randVal - 0.5);
221       /* Add the perturbation to the vertex*/
222       coordVals[dim * i + j] += amp;
223     }
224   }
225 
226   PetscRandomDestroy(&ran);
227   PetscFunctionReturn(0);
228 }
229 
230 /* Apply a global skew transformation to the mesh. */
231 static PetscErrorCode SkewMesh(DM *mesh, PetscScalar *coordVals, PetscInt npoints, PetscInt dim) {
232   PetscInt     i, j, k, l;
233   PetscScalar *transMat;
234   PetscScalar  tmpcoord[3];
235   PetscRandom  ran;
236   PetscReal    randVal;
237 
238   PetscFunctionBegin;
239   PetscCall(PetscCalloc1(dim * dim, &transMat));
240   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &ran));
241 
242   /* Make a matrix representing a skew transformation */
243   for (i = 0; i < dim; ++i) {
244     for (j = 0; j < dim; ++j) {
245       PetscCall(PetscRandomGetValueReal(ran, &randVal));
246       if (i == j) transMat[i * dim + j] = 1.;
247       else if (j < i) transMat[i * dim + j] = 2 * (j + i) * randVal;
248       else transMat[i * dim + j] = 0;
249     }
250   }
251 
252   /* Multiply each coordinate vector by our tranformation.*/
253   for (i = 0; i < npoints; ++i) {
254     for (j = 0; j < dim; ++j) {
255       tmpcoord[j] = 0;
256       for (k = 0; k < dim; ++k) tmpcoord[j] += coordVals[dim * i + k] * transMat[dim * k + j];
257     }
258     for (l = 0; l < dim; ++l) coordVals[dim * i + l] = tmpcoord[l];
259   }
260   PetscCall(PetscFree(transMat));
261   PetscCall(PetscRandomDestroy(&ran));
262   PetscFunctionReturn(0);
263 }
264 
265 /* Accesses the mesh coordinate array and performs the transformation operations
266  * specified by the user options */
267 static PetscErrorCode TransformMesh(UserCtx *user, DM *mesh) {
268   PetscInt     dim, npoints;
269   PetscScalar *coordVals;
270   Vec          coords;
271 
272   PetscFunctionBegin;
273   PetscCall(DMGetCoordinates(*mesh, &coords));
274   PetscCall(VecGetArray(coords, &coordVals));
275   PetscCall(VecGetLocalSize(coords, &npoints));
276   PetscCall(DMGetCoordinateDim(*mesh, &dim));
277   npoints = npoints / dim;
278 
279   switch (user->mesh_transform) {
280   case PERTURB: PetscCall(PerturbMesh(mesh, coordVals, npoints, dim)); break;
281   case SKEW: PetscCall(SkewMesh(mesh, coordVals, npoints, dim)); break;
282   case SKEW_PERTURB:
283     PetscCall(SkewMesh(mesh, coordVals, npoints, dim));
284     PetscCall(PerturbMesh(mesh, coordVals, npoints, dim));
285     break;
286   default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "invalid mesh transformation");
287   }
288   PetscCall(VecRestoreArray(coords, &coordVals));
289   PetscCall(DMSetCoordinates(*mesh, coords));
290   PetscFunctionReturn(0);
291 }
292 
293 static PetscErrorCode CreateMesh(MPI_Comm comm, UserCtx *user, DM *mesh) {
294   PetscFunctionBegin;
295   PetscCall(DMCreate(comm, mesh));
296   PetscCall(DMSetType(*mesh, DMPLEX));
297   PetscCall(DMSetFromOptions(*mesh));
298 
299   /* Perform any mesh transformations if specified by user */
300   if (user->mesh_transform != NONE) PetscCall(TransformMesh(user, mesh));
301 
302   /* Get any other mesh options from the command line */
303   PetscCall(DMSetApplicationContext(*mesh, user));
304   PetscCall(DMViewFromOptions(*mesh, NULL, "-dm_view"));
305   PetscFunctionReturn(0);
306 }
307 
308 /* Setup the system of equations that we wish to solve */
309 static PetscErrorCode SetupProblem(DM dm, UserCtx *user) {
310   PetscDS        prob;
311   DMLabel        label;
312   const PetscInt id = 1;
313 
314   PetscFunctionBegin;
315   PetscCall(DMGetDS(dm, &prob));
316   /* All of these are independent of the user's choice of solution */
317   PetscCall(PetscDSSetResidual(prob, 1, f0_q, NULL));
318   PetscCall(PetscDSSetResidual(prob, 2, f0_w, NULL));
319   PetscCall(PetscDSSetJacobian(prob, 0, 0, g0_vu, NULL, NULL, NULL));
320   PetscCall(PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL, NULL));
321   PetscCall(PetscDSSetJacobian(prob, 1, 1, g0_qp, NULL, NULL, NULL));
322   PetscCall(PetscDSSetJacobian(prob, 2, 0, NULL, g1_wu, NULL, NULL));
323   PetscCall(PetscDSSetJacobian(prob, 2, 1, g0_wp, NULL, NULL, NULL));
324   PetscCall(PetscDSSetJacobian(prob, 2, 2, g0_wd, NULL, NULL, NULL));
325 
326   /* Field 2 is the error between \div{u} and pressure in a higher dimensional
327    * space. If all is right this should be machine zero. */
328   PetscCall(PetscDSSetExactSolution(prob, 2, zero_func, NULL));
329 
330   switch (user->sol_form) {
331   case LINEAR:
332     PetscCall(PetscDSSetResidual(prob, 0, f0_v_linear, NULL));
333     PetscCall(PetscDSSetBdResidual(prob, 0, f0_bd_u_linear, NULL));
334     PetscCall(PetscDSSetExactSolution(prob, 0, linear_u, NULL));
335     PetscCall(PetscDSSetExactSolution(prob, 1, linear_p, NULL));
336     break;
337   case SINUSOIDAL:
338     PetscCall(PetscDSSetResidual(prob, 0, f0_v_sinusoid, NULL));
339     PetscCall(PetscDSSetBdResidual(prob, 0, f0_bd_u_sinusoid, NULL));
340     PetscCall(PetscDSSetExactSolution(prob, 0, sinusoid_u, NULL));
341     PetscCall(PetscDSSetExactSolution(prob, 1, sinusoid_p, NULL));
342     break;
343   default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "invalid solution form");
344   }
345 
346   PetscCall(DMGetLabel(dm, "marker", &label));
347   PetscCall(PetscDSAddBoundary(prob, DM_BC_NATURAL, "Boundary Integral", label, 1, &id, 0, 0, NULL, (void (*)(void))NULL, NULL, user, NULL));
348   PetscFunctionReturn(0);
349 }
350 
351 /* Create the finite element spaces we will use for this system */
352 static PetscErrorCode SetupDiscretization(DM mesh, PetscErrorCode (*setup)(DM, UserCtx *), UserCtx *user) {
353   DM        cdm = mesh;
354   PetscFE   fevel, fepres, fedivErr;
355   PetscInt  dim;
356   PetscBool simplex;
357 
358   PetscFunctionBegin;
359   PetscCall(DMGetDimension(mesh, &dim));
360   PetscCall(DMPlexIsSimplex(mesh, &simplex));
361   /* Create FE objects and give them names so that options can be set from
362    * command line */
363   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh), dim, dim, simplex, "velocity_", -1, &fevel));
364   PetscCall(PetscObjectSetName((PetscObject)fevel, "velocity"));
365 
366   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh), dim, 1, simplex, "pressure_", -1, &fepres));
367   PetscCall(PetscObjectSetName((PetscObject)fepres, "pressure"));
368 
369   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh), dim, 1, simplex, "divErr_", -1, &fedivErr));
370   PetscCall(PetscObjectSetName((PetscObject)fedivErr, "divErr"));
371 
372   PetscCall(PetscFECopyQuadrature(fevel, fepres));
373   PetscCall(PetscFECopyQuadrature(fevel, fedivErr));
374 
375   /* Associate the FE objects with the mesh and setup the system */
376   PetscCall(DMSetField(mesh, 0, NULL, (PetscObject)fevel));
377   PetscCall(DMSetField(mesh, 1, NULL, (PetscObject)fepres));
378   PetscCall(DMSetField(mesh, 2, NULL, (PetscObject)fedivErr));
379   PetscCall(DMCreateDS(mesh));
380   PetscCall((*setup)(mesh, user));
381 
382   while (cdm) {
383     PetscCall(DMCopyDisc(mesh, cdm));
384     PetscCall(DMGetCoarseDM(cdm, &cdm));
385   }
386 
387   /* The Mesh now owns the fields, so we can destroy the FEs created in this
388    * function */
389   PetscCall(PetscFEDestroy(&fevel));
390   PetscCall(PetscFEDestroy(&fepres));
391   PetscCall(PetscFEDestroy(&fedivErr));
392   PetscCall(DMDestroy(&cdm));
393   PetscFunctionReturn(0);
394 }
395 
396 int main(int argc, char **argv) {
397   PetscInt        i;
398   UserCtx         user;
399   DM              mesh;
400   SNES            snes;
401   Vec             computed, divErr;
402   PetscReal       divErrNorm;
403   IS             *fieldIS;
404   PetscBool       exampleSuccess = PETSC_FALSE;
405   const PetscReal errTol         = 10. * PETSC_SMALL;
406 
407   char stdFormat[] = "L2 Norm of the Divergence Error is: %g\n H(div) elements working correctly: %s\n";
408 
409   /* Initialize PETSc */
410   PetscFunctionBeginUser;
411   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
412   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
413 
414   /* Set up the system, we need to create a solver and a mesh and then assign
415    * the correct spaces into the mesh*/
416   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
417   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &mesh));
418   PetscCall(SNESSetDM(snes, mesh));
419   PetscCall(SetupDiscretization(mesh, SetupProblem, &user));
420   PetscCall(DMPlexSetSNESLocalFEM(mesh, &user, &user, &user));
421   PetscCall(SNESSetFromOptions(snes));
422 
423   /* Grab field IS so that we can view the solution by field */
424   PetscCall(DMCreateFieldIS(mesh, NULL, NULL, &fieldIS));
425 
426   /* Create a vector to store the SNES solution, solve the system and grab the
427    * solution from SNES */
428   PetscCall(DMCreateGlobalVector(mesh, &computed));
429   PetscCall(PetscObjectSetName((PetscObject)computed, "computedSolution"));
430   PetscCall(VecSet(computed, 0.0));
431   PetscCall(SNESSolve(snes, NULL, computed));
432   PetscCall(SNESGetSolution(snes, &computed));
433   PetscCall(VecViewFromOptions(computed, NULL, "-computedSolution_view"));
434 
435   /* Now we pull out the portion of the vector corresponding to the 3rd field
436    * which is the error between \div{u} computed in a higher dimensional space
437    * and p = \div{u} computed in a low dimension space. We report the L2 norm of
438    * this vector which should be zero if the H(div) spaces are implemented
439    * correctly. */
440   PetscCall(VecGetSubVector(computed, fieldIS[2], &divErr));
441   PetscCall(VecNorm(divErr, NORM_2, &divErrNorm));
442   PetscCall(VecRestoreSubVector(computed, fieldIS[2], &divErr));
443   exampleSuccess = (PetscBool)(divErrNorm <= errTol);
444 
445   PetscCall(PetscPrintf(PETSC_COMM_WORLD, stdFormat, divErrNorm, exampleSuccess ? "true" : "false"));
446 
447   /* Tear down */
448   PetscCall(VecDestroy(&divErr));
449   PetscCall(VecDestroy(&computed));
450   for (i = 0; i < 3; ++i) PetscCall(ISDestroy(&fieldIS[i]));
451   PetscCall(PetscFree(fieldIS));
452   PetscCall(SNESDestroy(&snes));
453   PetscCall(DMDestroy(&mesh));
454   PetscCall(PetscFinalize());
455   return 0;
456 }
457 
458 /*TEST
459   testset:
460     suffix: 2d_bdm
461     requires: triangle
462     args: -velocity_petscfe_default_quadrature_order 1 \
463       -velocity_petscspace_degree 1 \
464       -velocity_petscdualspace_type bdm \
465       -divErr_petscspace_degree 1 \
466       -divErr_petscdualspace_lagrange_continuity false \
467       -snes_error_if_not_converged \
468       -ksp_rtol 1e-10 \
469       -ksp_error_if_not_converged \
470       -pc_type fieldsplit\
471       -pc_fieldsplit_detect_saddle_point\
472       -pc_fieldsplit_type schur\
473       -pc_fieldsplit_schur_precondition full
474     test:
475       suffix: linear
476       args: -sol_form linear -mesh_transform none
477     test:
478       suffix: sinusoidal
479       args: -sol_form sinusoidal -mesh_transform none
480     test:
481       suffix: sinusoidal_skew
482       args: -sol_form sinusoidal -mesh_transform skew
483     test:
484       suffix: sinusoidal_perturb
485       args: -sol_form sinusoidal -mesh_transform perturb
486     test:
487       suffix: sinusoidal_skew_perturb
488       args: -sol_form sinusoidal -mesh_transform skew_perturb
489 
490   testset:
491     TODO: broken
492     suffix: 2d_bdmq
493     args: -dm_plex_simplex false \
494       -velocity_petscspace_degree 1 \
495       -velocity_petscdualspace_type bdm \
496       -velocity_petscdualspace_lagrange_tensor 1 \
497       -divErr_petscspace_degree 1 \
498       -divErr_petscdualspace_lagrange_continuity false \
499       -snes_error_if_not_converged \
500       -ksp_rtol 1e-10 \
501       -ksp_error_if_not_converged \
502       -pc_type fieldsplit\
503       -pc_fieldsplit_detect_saddle_point\
504       -pc_fieldsplit_type schur\
505       -pc_fieldsplit_schur_precondition full
506     test:
507       suffix: linear
508       args: -sol_form linear -mesh_transform none
509     test:
510       suffix: sinusoidal
511       args: -sol_form sinusoidal -mesh_transform none
512     test:
513       suffix: sinusoidal_skew
514       args: -sol_form sinusoidal -mesh_transform skew
515     test:
516       suffix: sinusoidal_perturb
517       args: -sol_form sinusoidal -mesh_transform perturb
518     test:
519       suffix: sinusoidal_skew_perturb
520       args: -sol_form sinusoidal -mesh_transform skew_perturb
521 
522   testset:
523     suffix: 3d_bdm
524     requires: ctetgen
525     args: -dm_plex_dim 3 \
526       -velocity_petscspace_degree 1 \
527       -velocity_petscdualspace_type bdm \
528       -divErr_petscspace_degree 1 \
529       -divErr_petscdualspace_lagrange_continuity false \
530       -snes_error_if_not_converged \
531       -ksp_rtol 1e-10 \
532       -ksp_error_if_not_converged \
533       -pc_type fieldsplit \
534       -pc_fieldsplit_detect_saddle_point \
535       -pc_fieldsplit_type schur \
536       -pc_fieldsplit_schur_precondition full
537     test:
538       suffix: linear
539       args: -sol_form linear -mesh_transform none
540     test:
541       suffix: sinusoidal
542       args: -sol_form sinusoidal -mesh_transform none
543     test:
544       suffix: sinusoidal_skew
545       args: -sol_form sinusoidal -mesh_transform skew
546     test:
547       suffix: sinusoidal_perturb
548       args: -sol_form sinusoidal -mesh_transform perturb
549     test:
550       suffix: sinusoidal_skew_perturb
551       args: -sol_form sinusoidal -mesh_transform skew_perturb
552 
553   testset:
554     TODO: broken
555     suffix: 3d_bdmq
556     requires: ctetgen
557     args: -dm_plex_dim 3 \
558       -dm_plex_simplex false \
559       -velocity_petscspace_degree 1 \
560       -velocity_petscdualspace_type bdm \
561       -velocity_petscdualspace_lagrange_tensor 1 \
562       -divErr_petscspace_degree 1 \
563       -divErr_petscdualspace_lagrange_continuity false \
564       -snes_error_if_not_converged \
565       -ksp_rtol 1e-10 \
566       -ksp_error_if_not_converged \
567       -pc_type fieldsplit \
568       -pc_fieldsplit_detect_saddle_point \
569       -pc_fieldsplit_type schur \
570       -pc_fieldsplit_schur_precondition full
571     test:
572       suffix: linear
573       args: -sol_form linear -mesh_transform none
574     test:
575       suffix: sinusoidal
576       args: -sol_form sinusoidal -mesh_transform none
577     test:
578       suffix: sinusoidal_skew
579       args: -sol_form sinusoidal -mesh_transform skew
580     test:
581       suffix: sinusoidal_perturb
582       args: -sol_form sinusoidal -mesh_transform perturb
583     test:
584       suffix: sinusoidal_skew_perturb
585       args: -sol_form sinusoidal -mesh_transform skew_perturb
586 
587   test:
588     suffix: quad_rt_0
589     args: -dm_plex_simplex false -mesh_transform skew \
590           -divErr_petscspace_degree 1 \
591           -divErr_petscdualspace_lagrange_continuity false \
592           -snes_error_if_not_converged \
593           -ksp_rtol 1e-10 \
594           -ksp_error_if_not_converged \
595           -pc_type fieldsplit\
596           -pc_fieldsplit_detect_saddle_point\
597           -pc_fieldsplit_type schur\
598           -pc_fieldsplit_schur_precondition full \
599           -velocity_petscfe_default_quadrature_order 1 \
600           -velocity_petscspace_type sum \
601           -velocity_petscspace_variables 2 \
602           -velocity_petscspace_components 2 \
603           -velocity_petscspace_sum_spaces 2 \
604           -velocity_petscspace_sum_concatenate true \
605           -velocity_sumcomp_0_petscspace_variables 2 \
606           -velocity_sumcomp_0_petscspace_type tensor \
607           -velocity_sumcomp_0_petscspace_tensor_spaces 2 \
608           -velocity_sumcomp_0_petscspace_tensor_uniform false \
609           -velocity_sumcomp_0_tensorcomp_0_petscspace_degree 1 \
610           -velocity_sumcomp_0_tensorcomp_1_petscspace_degree 0 \
611           -velocity_sumcomp_1_petscspace_variables 2 \
612           -velocity_sumcomp_1_petscspace_type tensor \
613           -velocity_sumcomp_1_petscspace_tensor_spaces 2 \
614           -velocity_sumcomp_1_petscspace_tensor_uniform false \
615           -velocity_sumcomp_1_tensorcomp_0_petscspace_degree 0 \
616           -velocity_sumcomp_1_tensorcomp_1_petscspace_degree 1 \
617           -velocity_petscdualspace_form_degree -1 \
618           -velocity_petscdualspace_order 1 \
619           -velocity_petscdualspace_lagrange_trimmed true
620 TEST*/
621