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