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