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