xref: /petsc/src/dm/impls/plex/tests/ex39.c (revision 800f99ff9e85495c69e9e5819c0be0dbd8cbc57c)
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   PetscFunctionBegin;
203   /* Default to  2D, unperturbed triangle mesh and Linear solution.*/
204   user->mesh_transform = NONE;
205   user->sol_form       = LINEAR;
206 
207   PetscOptionsBegin(comm,"","H-div Test Options","DMPLEX");
208   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));
209   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));
210   PetscOptionsEnd();
211   PetscFunctionReturn(0);
212 }
213 
214 /* Perturb the position of each mesh vertex by a small amount.*/
215 static PetscErrorCode PerturbMesh(DM *mesh,PetscScalar *coordVals,PetscInt npoints,PetscInt dim)
216 {
217   PetscInt       i,j,k;
218   PetscReal      minCoords[3],maxCoords[3],maxPert[3],randVal,amp;
219   PetscRandom    ran;
220 
221   PetscFunctionBegin;
222   PetscCall(DMGetCoordinateDim(*mesh,&dim));
223   PetscCall(DMGetLocalBoundingBox(*mesh,minCoords,maxCoords));
224   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&ran));
225 
226   /* Compute something approximately equal to half an edge length. This is the
227    * most we can perturb points and guarantee that there won't be any topology
228    * issues. */
229   for (k = 0; k < dim; ++k) maxPert[k] = 0.025 * (maxCoords[k] - minCoords[k]) / (PetscPowReal(npoints,1. / dim) - 1);
230   /* For each mesh vertex */
231   for (i = 0; i < npoints; ++i) {
232     /* For each coordinate of the vertex */
233     for (j = 0; j < dim; ++j) {
234       /* Generate a random amplitude in [-0.5*maxPert, 0.5*maxPert] */
235       PetscCall(PetscRandomGetValueReal(ran,&randVal));
236       amp  = maxPert[j] * (randVal - 0.5);
237       /* Add the perturbation to the vertex*/
238       coordVals[dim * i + j] += amp;
239     }
240   }
241 
242   PetscRandomDestroy(&ran);
243   PetscFunctionReturn(0);
244 }
245 
246 /* Apply a global skew transformation to the mesh. */
247 static PetscErrorCode SkewMesh(DM * mesh,PetscScalar * coordVals,PetscInt npoints,PetscInt dim)
248 {
249   PetscInt       i,j,k,l;
250   PetscScalar    * transMat;
251   PetscScalar    tmpcoord[3];
252   PetscRandom    ran;
253   PetscReal      randVal;
254 
255   PetscFunctionBegin;
256   PetscCall(PetscCalloc1(dim * dim,&transMat));
257   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&ran));
258 
259   /* Make a matrix representing a skew transformation */
260   for (i = 0; i < dim; ++i) {
261     for (j = 0; j < dim; ++j) {
262       PetscCall(PetscRandomGetValueReal(ran,&randVal));
263       if (i == j) transMat[i * dim + j] = 1.;
264       else if (j < i) transMat[i * dim + j] = 2 * (j + i)*randVal;
265       else transMat[i * dim + j] = 0;
266     }
267   }
268 
269   /* Multiply each coordinate vector by our tranformation.*/
270   for (i = 0; i < npoints; ++i) {
271     for (j = 0; j < dim; ++j) {
272       tmpcoord[j] = 0;
273       for (k = 0; k < dim; ++k) tmpcoord[j] += coordVals[dim * i + k] * transMat[dim * k + j];
274     }
275     for (l = 0; l < dim; ++l) coordVals[dim * i + l] = tmpcoord[l];
276   }
277   PetscCall(PetscFree(transMat));
278   PetscCall(PetscRandomDestroy(&ran));
279   PetscFunctionReturn(0);
280 }
281 
282 /* Accesses the mesh coordinate array and performs the transformation operations
283  * specified by the user options */
284 static PetscErrorCode TransformMesh(UserCtx * user,DM * mesh)
285 {
286   PetscInt       dim,npoints;
287   PetscScalar    * coordVals;
288   Vec            coords;
289 
290   PetscFunctionBegin;
291   PetscCall(DMGetCoordinates(*mesh,&coords));
292   PetscCall(VecGetArray(coords,&coordVals));
293   PetscCall(VecGetLocalSize(coords,&npoints));
294   PetscCall(DMGetCoordinateDim(*mesh,&dim));
295   npoints = npoints / dim;
296 
297   switch (user->mesh_transform) {
298   case PERTURB:
299     PetscCall(PerturbMesh(mesh,coordVals,npoints,dim));
300     break;
301   case SKEW:
302     PetscCall(SkewMesh(mesh,coordVals,npoints,dim));
303     break;
304   case SKEW_PERTURB:
305     PetscCall(SkewMesh(mesh,coordVals,npoints,dim));
306     PetscCall(PerturbMesh(mesh,coordVals,npoints,dim));
307     break;
308   default:
309     PetscFunctionReturn(-1);
310   }
311   PetscCall(VecRestoreArray(coords,&coordVals));
312   PetscCall(DMSetCoordinates(*mesh,coords));
313   PetscFunctionReturn(0);
314 }
315 
316 static PetscErrorCode CreateMesh(MPI_Comm comm,UserCtx * user,DM * mesh)
317 {
318   PetscFunctionBegin;
319   PetscCall(DMCreate(comm, mesh));
320   PetscCall(DMSetType(*mesh, DMPLEX));
321   PetscCall(DMSetFromOptions(*mesh));
322 
323   /* Perform any mesh transformations if specified by user */
324   if (user->mesh_transform != NONE) {
325     PetscCall(TransformMesh(user,mesh));
326   }
327 
328   /* Get any other mesh options from the command line */
329   PetscCall(DMSetApplicationContext(*mesh,user));
330   PetscCall(DMViewFromOptions(*mesh,NULL,"-dm_view"));
331   PetscFunctionReturn(0);
332 }
333 
334 /* Setup the system of equations that we wish to solve */
335 static PetscErrorCode SetupProblem(DM dm,UserCtx * user)
336 {
337   PetscDS        prob;
338   DMLabel        label;
339   const PetscInt id=1;
340 
341   PetscFunctionBegin;
342   PetscCall(DMGetDS(dm,&prob));
343   /* All of these are independent of the user's choice of solution */
344   PetscCall(PetscDSSetResidual(prob,1,f0_q,NULL));
345   PetscCall(PetscDSSetResidual(prob,2,f0_w,NULL));
346   PetscCall(PetscDSSetJacobian(prob,0,0,g0_vu,NULL,NULL,NULL));
347   PetscCall(PetscDSSetJacobian(prob,1,0,NULL,g1_qu,NULL,NULL));
348   PetscCall(PetscDSSetJacobian(prob,1,1,g0_qp,NULL,NULL,NULL));
349   PetscCall(PetscDSSetJacobian(prob,2,0,NULL,g1_wu,NULL,NULL));
350   PetscCall(PetscDSSetJacobian(prob,2,1,g0_wp,NULL,NULL,NULL));
351   PetscCall(PetscDSSetJacobian(prob,2,2,g0_wd,NULL,NULL,NULL));
352 
353   /* Field 2 is the error between \div{u} and pressure in a higher dimensional
354    * space. If all is right this should be machine zero. */
355   PetscCall(PetscDSSetExactSolution(prob,2,zero_func,NULL));
356 
357   switch (user->sol_form) {
358   case LINEAR:
359     PetscCall(PetscDSSetResidual(prob,0,f0_v_linear,NULL));
360     PetscCall(PetscDSSetBdResidual(prob,0,f0_bd_u_linear,NULL));
361     PetscCall(PetscDSSetExactSolution(prob,0,linear_u,NULL));
362     PetscCall(PetscDSSetExactSolution(prob,1,linear_p,NULL));
363     break;
364   case SINUSOIDAL:
365     PetscCall(PetscDSSetResidual(prob,0,f0_v_sinusoid,NULL));
366     PetscCall(PetscDSSetBdResidual(prob,0,f0_bd_u_sinusoid,NULL));
367     PetscCall(PetscDSSetExactSolution(prob,0,sinusoid_u,NULL));
368     PetscCall(PetscDSSetExactSolution(prob,1,sinusoid_p,NULL));
369     break;
370   default:
371     PetscFunctionReturn(-1);
372   }
373 
374   PetscCall(DMGetLabel(dm, "marker", &label));
375   PetscCall(PetscDSAddBoundary(prob,DM_BC_NATURAL,"Boundary Integral",label,1,&id,0,0,NULL,(void (*)(void))NULL,NULL,user,NULL));
376   PetscFunctionReturn(0);
377 }
378 
379 /* Create the finite element spaces we will use for this system */
380 static PetscErrorCode SetupDiscretization(DM mesh,PetscErrorCode (*setup)(DM,UserCtx*),UserCtx *user)
381 {
382   DM             cdm = mesh;
383   PetscFE        fevel,fepres,fedivErr;
384   PetscInt       dim;
385   PetscBool      simplex;
386 
387   PetscFunctionBegin;
388   PetscCall(DMGetDimension(mesh, &dim));
389   PetscCall(DMPlexIsSimplex(mesh, &simplex));
390   /* Create FE objects and give them names so that options can be set from
391    * command line */
392   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject) mesh),dim,dim,simplex,"velocity_",-1,&fevel));
393   PetscCall(PetscObjectSetName((PetscObject) fevel,"velocity"));
394 
395   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject) mesh),dim,1,simplex,"pressure_",-1,&fepres));
396   PetscCall(PetscObjectSetName((PetscObject) fepres,"pressure"));
397 
398   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh),dim,1,simplex,"divErr_",-1,&fedivErr));
399   PetscCall(PetscObjectSetName((PetscObject) fedivErr,"divErr"));
400 
401   PetscCall(PetscFECopyQuadrature(fevel,fepres));
402   PetscCall(PetscFECopyQuadrature(fevel,fedivErr));
403 
404   /* Associate the FE objects with the mesh and setup the system */
405   PetscCall(DMSetField(mesh,0,NULL,(PetscObject) fevel));
406   PetscCall(DMSetField(mesh,1,NULL,(PetscObject) fepres));
407   PetscCall(DMSetField(mesh,2,NULL,(PetscObject) fedivErr));
408   PetscCall(DMCreateDS(mesh));
409   PetscCall((*setup)(mesh,user));
410 
411   while (cdm) {
412     PetscCall(DMCopyDisc(mesh,cdm));
413     PetscCall(DMGetCoarseDM(cdm,&cdm));
414   }
415 
416   /* The Mesh now owns the fields, so we can destroy the FEs created in this
417    * function */
418   PetscCall(PetscFEDestroy(&fevel));
419   PetscCall(PetscFEDestroy(&fepres));
420   PetscCall(PetscFEDestroy(&fedivErr));
421   PetscCall(DMDestroy(&cdm));
422   PetscFunctionReturn(0);
423 }
424 
425 int main(int argc,char **argv)
426 {
427   PetscInt        i;
428   UserCtx         user;
429   DM              mesh;
430   SNES            snes;
431   Vec             computed,divErr;
432   PetscReal       divErrNorm;
433   IS              * fieldIS;
434   PetscBool       exampleSuccess = PETSC_FALSE;
435   const PetscReal errTol         = 10. * PETSC_SMALL;
436 
437   char stdFormat[] = "L2 Norm of the Divergence Error is: %g\n H(div) elements working correctly: %s\n";
438 
439   /* Initialize PETSc */
440   PetscFunctionBeginUser;
441   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
442   PetscCall(ProcessOptions(PETSC_COMM_WORLD,&user));
443 
444   /* Set up the system, we need to create a solver and a mesh and then assign
445    * the correct spaces into the mesh*/
446   PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes));
447   PetscCall(CreateMesh(PETSC_COMM_WORLD,&user,&mesh));
448   PetscCall(SNESSetDM(snes,mesh));
449   PetscCall(SetupDiscretization(mesh,SetupProblem,&user));
450   PetscCall(DMPlexSetSNESLocalFEM(mesh,&user,&user,&user));
451   PetscCall(SNESSetFromOptions(snes));
452 
453   /* Grab field IS so that we can view the solution by field */
454   PetscCall(DMCreateFieldIS(mesh,NULL,NULL,&fieldIS));
455 
456   /* Create a vector to store the SNES solution, solve the system and grab the
457    * solution from SNES */
458   PetscCall(DMCreateGlobalVector(mesh,&computed));
459   PetscCall(PetscObjectSetName((PetscObject) computed,"computedSolution"));
460   PetscCall(VecSet(computed,0.0));
461   PetscCall(SNESSolve(snes,NULL,computed));
462   PetscCall(SNESGetSolution(snes,&computed));
463   PetscCall(VecViewFromOptions(computed,NULL,"-computedSolution_view"));
464 
465   /* Now we pull out the portion of the vector corresponding to the 3rd field
466    * which is the error between \div{u} computed in a higher dimensional space
467    * and p = \div{u} computed in a low dimension space. We report the L2 norm of
468    * this vector which should be zero if the H(div) spaces are implemented
469    * correctly. */
470   PetscCall(VecGetSubVector(computed,fieldIS[2],&divErr));
471   PetscCall(VecNorm(divErr,NORM_2,&divErrNorm));
472   PetscCall(VecRestoreSubVector(computed,fieldIS[2],&divErr));
473   exampleSuccess = (PetscBool)(divErrNorm <= errTol);
474 
475   PetscCall(PetscPrintf(PETSC_COMM_WORLD,stdFormat,divErrNorm,exampleSuccess ? "true" : "false"));
476 
477   /* Tear down */
478   PetscCall(VecDestroy(&divErr));
479   PetscCall(VecDestroy(&computed));
480   for (i = 0; i < 3; ++i) {
481     PetscCall(ISDestroy(&fieldIS[i]));
482   }
483   PetscCall(PetscFree(fieldIS));
484   PetscCall(SNESDestroy(&snes));
485   PetscCall(DMDestroy(&mesh));
486   PetscCall(PetscFinalize());
487   return 0;
488 }
489 
490 /*TEST
491   testset:
492     suffix: 2d_bdm
493     requires: triangle
494     args: -velocity_petscfe_default_quadrature_order 1 \
495       -velocity_petscspace_degree 1 \
496       -velocity_petscdualspace_type bdm \
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     TODO: broken
524     suffix: 2d_bdmq
525     args: -dm_plex_simplex false \
526       -velocity_petscspace_degree 1 \
527       -velocity_petscdualspace_type bdm \
528       -velocity_petscdualspace_lagrange_tensor 1 \
529       -divErr_petscspace_degree 1 \
530       -divErr_petscdualspace_lagrange_continuity false \
531       -snes_error_if_not_converged \
532       -ksp_rtol 1e-10 \
533       -ksp_error_if_not_converged \
534       -pc_type fieldsplit\
535       -pc_fieldsplit_detect_saddle_point\
536       -pc_fieldsplit_type schur\
537       -pc_fieldsplit_schur_precondition full
538     test:
539       suffix: linear
540       args: -sol_form linear -mesh_transform none
541     test:
542       suffix: sinusoidal
543       args: -sol_form sinusoidal -mesh_transform none
544     test:
545       suffix: sinusoidal_skew
546       args: -sol_form sinusoidal -mesh_transform skew
547     test:
548       suffix: sinusoidal_perturb
549       args: -sol_form sinusoidal -mesh_transform perturb
550     test:
551       suffix: sinusoidal_skew_perturb
552       args: -sol_form sinusoidal -mesh_transform skew_perturb
553 
554   testset:
555     suffix: 3d_bdm
556     requires: ctetgen
557     args: -dm_plex_dim 3 \
558       -velocity_petscspace_degree 1 \
559       -velocity_petscdualspace_type bdm \
560       -divErr_petscspace_degree 1 \
561       -divErr_petscdualspace_lagrange_continuity false \
562       -snes_error_if_not_converged \
563       -ksp_rtol 1e-10 \
564       -ksp_error_if_not_converged \
565       -pc_type fieldsplit \
566       -pc_fieldsplit_detect_saddle_point \
567       -pc_fieldsplit_type schur \
568       -pc_fieldsplit_schur_precondition full
569     test:
570       suffix: linear
571       args: -sol_form linear -mesh_transform none
572     test:
573       suffix: sinusoidal
574       args: -sol_form sinusoidal -mesh_transform none
575     test:
576       suffix: sinusoidal_skew
577       args: -sol_form sinusoidal -mesh_transform skew
578     test:
579       suffix: sinusoidal_perturb
580       args: -sol_form sinusoidal -mesh_transform perturb
581     test:
582       suffix: sinusoidal_skew_perturb
583       args: -sol_form sinusoidal -mesh_transform skew_perturb
584 
585   testset:
586     TODO: broken
587     suffix: 3d_bdmq
588     requires: ctetgen
589     args: -dm_plex_dim 3 \
590       -dm_plex_simplex false \
591       -velocity_petscspace_degree 1 \
592       -velocity_petscdualspace_type bdm \
593       -velocity_petscdualspace_lagrange_tensor 1 \
594       -divErr_petscspace_degree 1 \
595       -divErr_petscdualspace_lagrange_continuity false \
596       -snes_error_if_not_converged \
597       -ksp_rtol 1e-10 \
598       -ksp_error_if_not_converged \
599       -pc_type fieldsplit \
600       -pc_fieldsplit_detect_saddle_point \
601       -pc_fieldsplit_type schur \
602       -pc_fieldsplit_schur_precondition full
603     test:
604       suffix: linear
605       args: -sol_form linear -mesh_transform none
606     test:
607       suffix: sinusoidal
608       args: -sol_form sinusoidal -mesh_transform none
609     test:
610       suffix: sinusoidal_skew
611       args: -sol_form sinusoidal -mesh_transform skew
612     test:
613       suffix: sinusoidal_perturb
614       args: -sol_form sinusoidal -mesh_transform perturb
615     test:
616       suffix: sinusoidal_skew_perturb
617       args: -sol_form sinusoidal -mesh_transform skew_perturb
618 
619   test:
620     suffix: quad_rt_0
621     args: -dm_plex_simplex false -mesh_transform skew \
622           -divErr_petscspace_degree 1 \
623           -divErr_petscdualspace_lagrange_continuity false \
624           -snes_error_if_not_converged \
625           -ksp_rtol 1e-10 \
626           -ksp_error_if_not_converged \
627           -pc_type fieldsplit\
628           -pc_fieldsplit_detect_saddle_point\
629           -pc_fieldsplit_type schur\
630           -pc_fieldsplit_schur_precondition full \
631           -velocity_petscfe_default_quadrature_order 1 \
632           -velocity_petscspace_type sum \
633           -velocity_petscspace_variables 2 \
634           -velocity_petscspace_components 2 \
635           -velocity_petscspace_sum_spaces 2 \
636           -velocity_petscspace_sum_concatenate true \
637           -velocity_sumcomp_0_petscspace_variables 2 \
638           -velocity_sumcomp_0_petscspace_type tensor \
639           -velocity_sumcomp_0_petscspace_tensor_spaces 2 \
640           -velocity_sumcomp_0_petscspace_tensor_uniform false \
641           -velocity_sumcomp_0_tensorcomp_0_petscspace_degree 1 \
642           -velocity_sumcomp_0_tensorcomp_1_petscspace_degree 0 \
643           -velocity_sumcomp_1_petscspace_variables 2 \
644           -velocity_sumcomp_1_petscspace_type tensor \
645           -velocity_sumcomp_1_petscspace_tensor_spaces 2 \
646           -velocity_sumcomp_1_petscspace_tensor_uniform false \
647           -velocity_sumcomp_1_tensorcomp_0_petscspace_degree 0 \
648           -velocity_sumcomp_1_tensorcomp_1_petscspace_degree 1 \
649           -velocity_petscdualspace_form_degree -1 \
650           -velocity_petscdualspace_order 1 \
651           -velocity_petscdualspace_lagrange_trimmed true
652 TEST*/
653