xref: /petsc/src/dm/impls/plex/tests/ex39.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
130602db0SMatthew G. Knepley const char help[] = "A test of H-div conforming discretizations on different cell types.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown #include <petscds.h>
5c4762a1bSJed Brown #include <petscsnes.h>
6c4762a1bSJed Brown #include <petscconvest.h>
7c4762a1bSJed Brown #include <petscfe.h>
8c4762a1bSJed Brown #include <petsc/private/petscfeimpl.h>
9c4762a1bSJed Brown 
10c4762a1bSJed Brown /*
110e3d61c9SBarry Smith   We are using the system
120e3d61c9SBarry Smith 
130e3d61c9SBarry Smith   \vec{u} = \vec{\hat{u}}
140e3d61c9SBarry Smith   p = \div{\vec{u}} in low degree approximation space
150e3d61c9SBarry Smith   d = \div{\vec{u}} - p == 0 in higher degree approximation space
160e3d61c9SBarry Smith 
170e3d61c9SBarry Smith   That is, we are using the field d to compute the error between \div{\vec{u}}
180e3d61c9SBarry Smith   computed in a space 1 degree higher than p and the value of p which is
190e3d61c9SBarry Smith   \div{u} computed in the low degree space. If H-div
200e3d61c9SBarry Smith   elements are implemented correctly then this should be identically zero since
210e3d61c9SBarry Smith   the divergence of a function in H(div) should be exactly representable in L_2
220e3d61c9SBarry Smith   by definition.
23c4762a1bSJed Brown */
24c4762a1bSJed Brown static PetscErrorCode zero_func(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar *u,void *ctx)
25c4762a1bSJed Brown {
26c4762a1bSJed Brown   PetscInt c;
27c4762a1bSJed Brown   for (c = 0; c < Nc; ++c) u[c] = 0;
28c4762a1bSJed Brown   return 0;
29c4762a1bSJed Brown }
30c4762a1bSJed Brown /* Linear Exact Functions
31c4762a1bSJed Brown    \vec{u} = \vec{x};
32c4762a1bSJed Brown    p = dim;
33c4762a1bSJed Brown    */
34c4762a1bSJed Brown static PetscErrorCode linear_u(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar *u,void *ctx)
35c4762a1bSJed Brown {
36c4762a1bSJed Brown   PetscInt c;
37c4762a1bSJed Brown   for (c = 0; c < Nc; ++c) u[c] = x[c];
38c4762a1bSJed Brown   return 0;
39c4762a1bSJed Brown }
40c4762a1bSJed Brown static PetscErrorCode linear_p(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar *u,void *ctx)
41c4762a1bSJed Brown {
42c4762a1bSJed Brown   u[0] = dim;
43c4762a1bSJed Brown   return 0;
44c4762a1bSJed Brown }
45c4762a1bSJed Brown 
46c4762a1bSJed Brown /* Sinusoidal Exact Functions
47c4762a1bSJed Brown  * u_i = \sin{2*\pi*x_i}
48c4762a1bSJed Brown  * p = \Sum_{i=1}^{dim} 2*\pi*cos{2*\pi*x_i}
49c4762a1bSJed Brown  * */
50c4762a1bSJed Brown 
51c4762a1bSJed Brown static PetscErrorCode sinusoid_u(PetscInt dim,PetscReal time,const PetscReal
52c4762a1bSJed Brown                                  x[],PetscInt Nc,PetscScalar *u,void *ctx)
53c4762a1bSJed Brown {
54c4762a1bSJed Brown   PetscInt c;
55c4762a1bSJed Brown   for (c = 0; c< Nc; ++c) u[c] = PetscSinReal(2*PETSC_PI*x[c]);
56c4762a1bSJed Brown   return 0;
57c4762a1bSJed Brown }
58c4762a1bSJed Brown static PetscErrorCode sinusoid_p(PetscInt dim,PetscReal time,const PetscReal
59c4762a1bSJed Brown                                  x[],PetscInt Nc,PetscScalar *u,void *ctx)
60c4762a1bSJed Brown {
61c4762a1bSJed Brown   PetscInt d;
62c4762a1bSJed Brown   u[0]=0;
63c4762a1bSJed Brown   for (d=0; d<dim; ++d) u[0] += 2*PETSC_PI*PetscCosReal(2*PETSC_PI*x[d]);
64c4762a1bSJed Brown   return 0;
65c4762a1bSJed Brown }
66c4762a1bSJed Brown 
67c4762a1bSJed Brown /* Pointwise residual for u = u*. Need one of these for each possible u* */
68c4762a1bSJed 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[])
69c4762a1bSJed Brown {
70c4762a1bSJed Brown   PetscInt    i;
71c4762a1bSJed Brown   PetscScalar *u_rhs;
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   PetscCalloc1(dim,&u_rhs);
74c4762a1bSJed Brown   (void) linear_u(dim,t,x,dim,u_rhs,NULL);
75c4762a1bSJed Brown   for (i = 0; i < dim; ++i) f0[i] = u[uOff[0]+i]-u_rhs[i];
76c4762a1bSJed Brown   PetscFree(u_rhs);
77c4762a1bSJed Brown }
78c4762a1bSJed Brown 
79c4762a1bSJed 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[])
80c4762a1bSJed Brown {
81c4762a1bSJed Brown   PetscInt    i;
82c4762a1bSJed Brown   PetscScalar *u_rhs;
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   PetscCalloc1(dim,&u_rhs);
85c4762a1bSJed Brown   (void) sinusoid_u(dim,t,x,dim,u_rhs,NULL);
86c4762a1bSJed Brown   for (i = 0; i < dim; ++i) f0[i] = u[uOff[0]+i]-u_rhs[i];
87c4762a1bSJed Brown   PetscFree(u_rhs);
88c4762a1bSJed Brown }
89c4762a1bSJed Brown 
90c4762a1bSJed Brown /* Residual function for enforcing p = \div{u}. */
91c4762a1bSJed 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[])
92c4762a1bSJed Brown {
93c4762a1bSJed Brown   PetscInt    i;
94c4762a1bSJed Brown   PetscScalar divu;
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   divu = 0.;
97c4762a1bSJed Brown   for (i = 0; i< dim; ++i) divu += u_x[uOff_x[0]+i*dim+i];
98c4762a1bSJed Brown   f0[0] = u[uOff[1]] - divu;
99c4762a1bSJed Brown }
100c4762a1bSJed Brown 
101c4762a1bSJed Brown /* Residual function for p_err = \div{u} - p. */
102c4762a1bSJed 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[])
103c4762a1bSJed Brown {
104c4762a1bSJed Brown   PetscInt    i;
105c4762a1bSJed Brown   PetscScalar divu;
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   divu = 0.;
108c4762a1bSJed Brown   for (i = 0; i < dim; ++i) divu += u_x[uOff_x[0] + i*dim +i];
109c4762a1bSJed Brown   f0[0] = u[uOff[2]] - u[uOff[1]] + divu;
110c4762a1bSJed Brown }
111c4762a1bSJed Brown 
112c4762a1bSJed Brown /* Boundary residual for the embedding system. Need one for each form of
113c4762a1bSJed Brown  * solution. These enforce u = \hat{u} at the boundary. */
114c4762a1bSJed 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[])
115c4762a1bSJed Brown {
116c4762a1bSJed Brown   PetscInt    d;
117c4762a1bSJed Brown   PetscScalar *u_rhs;
118c4762a1bSJed Brown   PetscCalloc1(dim,&u_rhs);
119c4762a1bSJed Brown   (void) sinusoid_u(dim,t,x,dim,u_rhs,NULL);
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   for (d=0; d<dim; ++d) f0[d] = u_rhs[d];
122c4762a1bSJed Brown   PetscFree(u_rhs);
123c4762a1bSJed Brown 
124c4762a1bSJed Brown }
125c4762a1bSJed Brown 
126c4762a1bSJed 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[])
127c4762a1bSJed Brown {
128c4762a1bSJed Brown   PetscInt    d;
129c4762a1bSJed Brown   PetscScalar *u_rhs;
130c4762a1bSJed Brown   PetscCalloc1(dim,&u_rhs);
131c4762a1bSJed Brown   (void) linear_u(dim,t,x,dim,u_rhs,NULL);
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   for (d=0; d<dim; ++d) f0[d] = u_rhs[d];
134c4762a1bSJed Brown   PetscFree(u_rhs);
135c4762a1bSJed Brown }
136c4762a1bSJed Brown /* Jacobian functions. For the following, v is the test function associated with
137c4762a1bSJed Brown  * u, q the test function associated with p, and w the test function associated
138c4762a1bSJed Brown  * with d. */
139c4762a1bSJed Brown /* <v, u> */
140c4762a1bSJed 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[])
141c4762a1bSJed Brown {
142c4762a1bSJed Brown   PetscInt c;
143c4762a1bSJed Brown   for (c = 0; c < dim; ++c) g0[c * dim + c] = 1.0;
144c4762a1bSJed Brown }
145c4762a1bSJed Brown 
146c4762a1bSJed Brown /* <q, p> */
147c4762a1bSJed 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[])
148c4762a1bSJed Brown {
149c4762a1bSJed Brown   PetscInt d;
150c4762a1bSJed Brown   for (d=0; d< dim; ++d) g0[d * dim + d] = 1.0;
151c4762a1bSJed Brown }
152c4762a1bSJed Brown 
153c4762a1bSJed Brown /* -<q, \div{u}> For the embedded system. This is different from the method of
154c4762a1bSJed Brown  * manufactured solution because instead of computing <q,\div{u}> - <q,f> we
155c4762a1bSJed Brown  * need <q,p> - <q,\div{u}.*/
156c4762a1bSJed 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[])
157c4762a1bSJed Brown {
158c4762a1bSJed Brown   PetscInt d;
159c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g1[d * dim + d] = -1.0;
160c4762a1bSJed Brown }
161c4762a1bSJed Brown 
162c4762a1bSJed Brown /* <w, p> This is only used by the embedded system. Where we need to compute
163c4762a1bSJed Brown  * <w,d> - <w,p> + <w, \div{u}>*/
164c4762a1bSJed 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[])
165c4762a1bSJed Brown {
166c4762a1bSJed Brown   PetscInt d;
167c4762a1bSJed Brown   for (d=0; d< dim; ++d) g0[d * dim + d] = -1.0;
168c4762a1bSJed Brown }
169c4762a1bSJed Brown 
170c4762a1bSJed Brown /* <w, d> */
171c4762a1bSJed 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[])
172c4762a1bSJed Brown {
173c4762a1bSJed Brown   PetscInt c;
174c4762a1bSJed Brown   for (c = 0; c < dim; ++c) g0[c*dim+c] = 1.0;
175c4762a1bSJed Brown }
176c4762a1bSJed Brown 
177c4762a1bSJed Brown /* <w, \div{u}> for the embedded system. */
178c4762a1bSJed 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[])
179c4762a1bSJed Brown {
180c4762a1bSJed Brown   PetscInt d;
181c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g1[d * dim + d] = 1.0;
182c4762a1bSJed Brown }
183c4762a1bSJed Brown 
184c4762a1bSJed Brown /* Enum and string array for selecting mesh perturbation options */
185c4762a1bSJed Brown typedef enum { NONE = 0,PERTURB = 1,SKEW = 2,SKEW_PERTURB = 3 } Transform;
186c4762a1bSJed Brown const char* const TransformTypes[] = {"none","perturb","skew","skew_perturb","Perturbation","",NULL};
187c4762a1bSJed Brown 
188c4762a1bSJed Brown /* Enum and string array for selecting the form of the exact solution*/
189c4762a1bSJed Brown typedef enum
190c4762a1bSJed Brown { LINEAR = 0,SINUSOIDAL = 1 } Solution;
191c4762a1bSJed Brown const char* const SolutionTypes[] = {"linear","sinusoidal","Solution","",NULL};
192c4762a1bSJed Brown 
193c4762a1bSJed Brown typedef struct
194c4762a1bSJed Brown {
195c4762a1bSJed Brown   Transform mesh_transform;
196c4762a1bSJed Brown   Solution  sol_form;
197c4762a1bSJed Brown } UserCtx;
198c4762a1bSJed Brown 
199c4762a1bSJed Brown /* Process command line options and initialize the UserCtx struct */
200c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm,UserCtx * user)
201c4762a1bSJed Brown {
202c4762a1bSJed Brown   PetscFunctionBegin;
203c4762a1bSJed Brown   /* Default to  2D, unperturbed triangle mesh and Linear solution.*/
204c4762a1bSJed Brown   user->mesh_transform = NONE;
205c4762a1bSJed Brown   user->sol_form       = LINEAR;
206c4762a1bSJed Brown 
207d0609cedSBarry Smith   PetscOptionsBegin(comm,"","H-div Test Options","DMPLEX");
2089566063dSJacob Faibussowitsch   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));
2099566063dSJacob Faibussowitsch   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));
210d0609cedSBarry Smith   PetscOptionsEnd();
211c4762a1bSJed Brown   PetscFunctionReturn(0);
212c4762a1bSJed Brown }
213c4762a1bSJed Brown 
214c4762a1bSJed Brown /* Perturb the position of each mesh vertex by a small amount.*/
215c4762a1bSJed Brown static PetscErrorCode PerturbMesh(DM *mesh,PetscScalar *coordVals,PetscInt npoints,PetscInt dim)
216c4762a1bSJed Brown {
217c4762a1bSJed Brown   PetscInt       i,j,k;
218d092c84bSBrandon Whitchurch   PetscReal      minCoords[3],maxCoords[3],maxPert[3],randVal,amp;
219c4762a1bSJed Brown   PetscRandom    ran;
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   PetscFunctionBegin;
2229566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(*mesh,&dim));
2239566063dSJacob Faibussowitsch   PetscCall(DMGetLocalBoundingBox(*mesh,minCoords,maxCoords));
2249566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&ran));
225c4762a1bSJed Brown 
226c4762a1bSJed Brown   /* Compute something approximately equal to half an edge length. This is the
227a5b23f4aSJose E. Roman    * most we can perturb points and guarantee that there won't be any topology
228c4762a1bSJed Brown    * issues. */
229d092c84bSBrandon Whitchurch   for (k = 0; k < dim; ++k) maxPert[k] = 0.025 * (maxCoords[k] - minCoords[k]) / (PetscPowReal(npoints,1. / dim) - 1);
230c4762a1bSJed Brown   /* For each mesh vertex */
231c4762a1bSJed Brown   for (i = 0; i < npoints; ++i) {
232c4762a1bSJed Brown     /* For each coordinate of the vertex */
233c4762a1bSJed Brown     for (j = 0; j < dim; ++j) {
234c4762a1bSJed Brown       /* Generate a random amplitude in [-0.5*maxPert, 0.5*maxPert] */
2359566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValueReal(ran,&randVal));
236c4762a1bSJed Brown       amp  = maxPert[j] * (randVal - 0.5);
237c4762a1bSJed Brown       /* Add the perturbation to the vertex*/
238d092c84bSBrandon Whitchurch       coordVals[dim * i + j] += amp;
239c4762a1bSJed Brown     }
240c4762a1bSJed Brown   }
241c4762a1bSJed Brown 
242c4762a1bSJed Brown   PetscRandomDestroy(&ran);
243c4762a1bSJed Brown   PetscFunctionReturn(0);
244c4762a1bSJed Brown }
245c4762a1bSJed Brown 
246c4762a1bSJed Brown /* Apply a global skew transformation to the mesh. */
247c4762a1bSJed Brown static PetscErrorCode SkewMesh(DM * mesh,PetscScalar * coordVals,PetscInt npoints,PetscInt dim)
248c4762a1bSJed Brown {
249c4762a1bSJed Brown   PetscInt       i,j,k,l;
250c4762a1bSJed Brown   PetscScalar    * transMat;
251c4762a1bSJed Brown   PetscScalar    tmpcoord[3];
252c4762a1bSJed Brown   PetscRandom    ran;
253c4762a1bSJed Brown   PetscReal      randVal;
254c4762a1bSJed Brown 
255c4762a1bSJed Brown   PetscFunctionBegin;
2569566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(dim * dim,&transMat));
2579566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&ran));
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   /* Make a matrix representing a skew transformation */
260c4762a1bSJed Brown   for (i = 0; i < dim; ++i) {
261c4762a1bSJed Brown     for (j = 0; j < dim; ++j) {
2629566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValueReal(ran,&randVal));
263d092c84bSBrandon Whitchurch       if (i == j) transMat[i * dim + j] = 1.;
264c4762a1bSJed Brown       else if (j < i) transMat[i * dim + j] = 2 * (j + i)*randVal;
265c4762a1bSJed Brown       else transMat[i * dim + j] = 0;
266c4762a1bSJed Brown     }
267c4762a1bSJed Brown   }
268c4762a1bSJed Brown 
269c4762a1bSJed Brown   /* Multiply each coordinate vector by our tranformation.*/
270c4762a1bSJed Brown   for (i = 0; i < npoints; ++i) {
271c4762a1bSJed Brown     for (j = 0; j < dim; ++j) {
272c4762a1bSJed Brown       tmpcoord[j] = 0;
273c4762a1bSJed Brown       for (k = 0; k < dim; ++k) tmpcoord[j] += coordVals[dim * i + k] * transMat[dim * k + j];
274c4762a1bSJed Brown     }
275c4762a1bSJed Brown     for (l = 0; l < dim; ++l) coordVals[dim * i + l] = tmpcoord[l];
276c4762a1bSJed Brown   }
2779566063dSJacob Faibussowitsch   PetscCall(PetscFree(transMat));
2789566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&ran));
279c4762a1bSJed Brown   PetscFunctionReturn(0);
280c4762a1bSJed Brown }
281c4762a1bSJed Brown 
282c4762a1bSJed Brown /* Accesses the mesh coordinate array and performs the transformation operations
283c4762a1bSJed Brown  * specified by the user options */
284c4762a1bSJed Brown static PetscErrorCode TransformMesh(UserCtx * user,DM * mesh)
285c4762a1bSJed Brown {
286c4762a1bSJed Brown   PetscInt       dim,npoints;
287c4762a1bSJed Brown   PetscScalar    * coordVals;
288c4762a1bSJed Brown   Vec            coords;
289c4762a1bSJed Brown 
290c4762a1bSJed Brown   PetscFunctionBegin;
2919566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(*mesh,&coords));
2929566063dSJacob Faibussowitsch   PetscCall(VecGetArray(coords,&coordVals));
2939566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(coords,&npoints));
2949566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(*mesh,&dim));
295c4762a1bSJed Brown   npoints = npoints / dim;
296c4762a1bSJed Brown 
297c4762a1bSJed Brown   switch (user->mesh_transform) {
298c4762a1bSJed Brown   case PERTURB:
2999566063dSJacob Faibussowitsch     PetscCall(PerturbMesh(mesh,coordVals,npoints,dim));
300c4762a1bSJed Brown     break;
301c4762a1bSJed Brown   case SKEW:
3029566063dSJacob Faibussowitsch     PetscCall(SkewMesh(mesh,coordVals,npoints,dim));
303c4762a1bSJed Brown     break;
304c4762a1bSJed Brown   case SKEW_PERTURB:
3059566063dSJacob Faibussowitsch     PetscCall(SkewMesh(mesh,coordVals,npoints,dim));
3069566063dSJacob Faibussowitsch     PetscCall(PerturbMesh(mesh,coordVals,npoints,dim));
307c4762a1bSJed Brown     break;
308c4762a1bSJed Brown   default:
309c4762a1bSJed Brown     PetscFunctionReturn(-1);
310c4762a1bSJed Brown   }
3119566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(coords,&coordVals));
3129566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinates(*mesh,coords));
313c4762a1bSJed Brown   PetscFunctionReturn(0);
314c4762a1bSJed Brown }
315c4762a1bSJed Brown 
316c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm,UserCtx * user,DM * mesh)
317c4762a1bSJed Brown {
318c4762a1bSJed Brown   PetscFunctionBegin;
3199566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, mesh));
3209566063dSJacob Faibussowitsch   PetscCall(DMSetType(*mesh, DMPLEX));
3219566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*mesh));
322c4762a1bSJed Brown 
323c4762a1bSJed Brown   /* Perform any mesh transformations if specified by user */
324c4762a1bSJed Brown   if (user->mesh_transform != NONE) {
3259566063dSJacob Faibussowitsch     PetscCall(TransformMesh(user,mesh));
326c4762a1bSJed Brown   }
327c4762a1bSJed Brown 
328c4762a1bSJed Brown   /* Get any other mesh options from the command line */
3299566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(*mesh,user));
3309566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*mesh,NULL,"-dm_view"));
331c4762a1bSJed Brown   PetscFunctionReturn(0);
332c4762a1bSJed Brown }
333c4762a1bSJed Brown 
334c4762a1bSJed Brown /* Setup the system of equations that we wish to solve */
335c4762a1bSJed Brown static PetscErrorCode SetupProblem(DM dm,UserCtx * user)
336c4762a1bSJed Brown {
337c4762a1bSJed Brown   PetscDS        prob;
33845480ffeSMatthew G. Knepley   DMLabel        label;
339c4762a1bSJed Brown   const PetscInt id=1;
340c4762a1bSJed Brown 
341c4762a1bSJed Brown   PetscFunctionBegin;
3429566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm,&prob));
343c4762a1bSJed Brown   /* All of these are independent of the user's choice of solution */
3449566063dSJacob Faibussowitsch   PetscCall(PetscDSSetResidual(prob,1,f0_q,NULL));
3459566063dSJacob Faibussowitsch   PetscCall(PetscDSSetResidual(prob,2,f0_w,NULL));
3469566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob,0,0,g0_vu,NULL,NULL,NULL));
3479566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob,1,0,NULL,g1_qu,NULL,NULL));
3489566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob,1,1,g0_qp,NULL,NULL,NULL));
3499566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob,2,0,NULL,g1_wu,NULL,NULL));
3509566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob,2,1,g0_wp,NULL,NULL,NULL));
3519566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob,2,2,g0_wd,NULL,NULL,NULL));
352c4762a1bSJed Brown 
353c4762a1bSJed Brown   /* Field 2 is the error between \div{u} and pressure in a higher dimensional
354c4762a1bSJed Brown    * space. If all is right this should be machine zero. */
3559566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(prob,2,zero_func,NULL));
356c4762a1bSJed Brown 
357c4762a1bSJed Brown   switch (user->sol_form) {
358c4762a1bSJed Brown   case LINEAR:
3599566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(prob,0,f0_v_linear,NULL));
3609566063dSJacob Faibussowitsch     PetscCall(PetscDSSetBdResidual(prob,0,f0_bd_u_linear,NULL));
3619566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(prob,0,linear_u,NULL));
3629566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(prob,1,linear_p,NULL));
363c4762a1bSJed Brown     break;
364c4762a1bSJed Brown   case SINUSOIDAL:
3659566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(prob,0,f0_v_sinusoid,NULL));
3669566063dSJacob Faibussowitsch     PetscCall(PetscDSSetBdResidual(prob,0,f0_bd_u_sinusoid,NULL));
3679566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(prob,0,sinusoid_u,NULL));
3689566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(prob,1,sinusoid_p,NULL));
369c4762a1bSJed Brown     break;
370c4762a1bSJed Brown   default:
371c4762a1bSJed Brown     PetscFunctionReturn(-1);
372c4762a1bSJed Brown   }
373c4762a1bSJed Brown 
3749566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
3759566063dSJacob Faibussowitsch   PetscCall(PetscDSAddBoundary(prob,DM_BC_NATURAL,"Boundary Integral",label,1,&id,0,0,NULL,(void (*)(void))NULL,NULL,user,NULL));
376c4762a1bSJed Brown   PetscFunctionReturn(0);
377c4762a1bSJed Brown }
378c4762a1bSJed Brown 
379c4762a1bSJed Brown /* Create the finite element spaces we will use for this system */
380c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM mesh,PetscErrorCode (*setup)(DM,UserCtx*),UserCtx *user)
381c4762a1bSJed Brown {
382c4762a1bSJed Brown   DM             cdm = mesh;
383c4762a1bSJed Brown   PetscFE        fevel,fepres,fedivErr;
38430602db0SMatthew G. Knepley   PetscInt       dim;
38530602db0SMatthew G. Knepley   PetscBool      simplex;
386c4762a1bSJed Brown 
387c4762a1bSJed Brown   PetscFunctionBegin;
3889566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(mesh, &dim));
3899566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(mesh, &simplex));
390c4762a1bSJed Brown   /* Create FE objects and give them names so that options can be set from
391c4762a1bSJed Brown    * command line */
3929566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject) mesh),dim,dim,simplex,"velocity_",-1,&fevel));
3939566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) fevel,"velocity"));
394c4762a1bSJed Brown 
3959566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject) mesh),dim,1,simplex,"pressure_",-1,&fepres));
3969566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) fepres,"pressure"));
397c4762a1bSJed Brown 
398d0609cedSBarry Smith   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh),dim,1,simplex,"divErr_",-1,&fedivErr));
3999566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) fedivErr,"divErr"));
400c4762a1bSJed Brown 
4019566063dSJacob Faibussowitsch   PetscCall(PetscFECopyQuadrature(fevel,fepres));
4029566063dSJacob Faibussowitsch   PetscCall(PetscFECopyQuadrature(fevel,fedivErr));
403c4762a1bSJed Brown 
404c4762a1bSJed Brown   /* Associate the FE objects with the mesh and setup the system */
4059566063dSJacob Faibussowitsch   PetscCall(DMSetField(mesh,0,NULL,(PetscObject) fevel));
4069566063dSJacob Faibussowitsch   PetscCall(DMSetField(mesh,1,NULL,(PetscObject) fepres));
4079566063dSJacob Faibussowitsch   PetscCall(DMSetField(mesh,2,NULL,(PetscObject) fedivErr));
4089566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(mesh));
4099566063dSJacob Faibussowitsch   PetscCall((*setup)(mesh,user));
410c4762a1bSJed Brown 
411c4762a1bSJed Brown   while (cdm) {
4129566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(mesh,cdm));
4139566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm,&cdm));
414c4762a1bSJed Brown   }
415c4762a1bSJed Brown 
416c4762a1bSJed Brown   /* The Mesh now owns the fields, so we can destroy the FEs created in this
417c4762a1bSJed Brown    * function */
4189566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fevel));
4199566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fepres));
4209566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fedivErr));
4219566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&cdm));
422c4762a1bSJed Brown   PetscFunctionReturn(0);
423c4762a1bSJed Brown }
424c4762a1bSJed Brown 
425c4762a1bSJed Brown int main(int argc,char **argv)
426c4762a1bSJed Brown {
427c4762a1bSJed Brown   PetscInt        i;
428c4762a1bSJed Brown   UserCtx         user;
429c4762a1bSJed Brown   DM              mesh;
430c4762a1bSJed Brown   SNES            snes;
431c4762a1bSJed Brown   Vec             computed,divErr;
432c4762a1bSJed Brown   PetscReal       divErrNorm;
433c4762a1bSJed Brown   IS              * fieldIS;
434c4762a1bSJed Brown   PetscBool       exampleSuccess = PETSC_FALSE;
435d092c84bSBrandon Whitchurch   const PetscReal errTol         = 10. * PETSC_SMALL;
436c4762a1bSJed Brown 
437c4762a1bSJed Brown   char stdFormat[] = "L2 Norm of the Divergence Error is: %g\n H(div) elements working correctly: %s\n";
438c4762a1bSJed Brown 
439c4762a1bSJed Brown   /* Initialize PETSc */
440*327415f7SBarry Smith   PetscFunctionBeginUser;
4419566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
4429566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD,&user));
443c4762a1bSJed Brown 
444c4762a1bSJed Brown   /* Set up the system, we need to create a solver and a mesh and then assign
445c4762a1bSJed Brown    * the correct spaces into the mesh*/
4469566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes));
4479566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD,&user,&mesh));
4489566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes,mesh));
4499566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(mesh,SetupProblem,&user));
4509566063dSJacob Faibussowitsch   PetscCall(DMPlexSetSNESLocalFEM(mesh,&user,&user,&user));
4519566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
452c4762a1bSJed Brown 
453c4762a1bSJed Brown   /* Grab field IS so that we can view the solution by field */
4549566063dSJacob Faibussowitsch   PetscCall(DMCreateFieldIS(mesh,NULL,NULL,&fieldIS));
455c4762a1bSJed Brown 
456c4762a1bSJed Brown   /* Create a vector to store the SNES solution, solve the system and grab the
457c4762a1bSJed Brown    * solution from SNES */
4589566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(mesh,&computed));
4599566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) computed,"computedSolution"));
4609566063dSJacob Faibussowitsch   PetscCall(VecSet(computed,0.0));
4619566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes,NULL,computed));
4629566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes,&computed));
4639566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(computed,NULL,"-computedSolution_view"));
464c4762a1bSJed Brown 
465c4762a1bSJed Brown   /* Now we pull out the portion of the vector corresponding to the 3rd field
466c4762a1bSJed Brown    * which is the error between \div{u} computed in a higher dimensional space
467c4762a1bSJed Brown    * and p = \div{u} computed in a low dimension space. We report the L2 norm of
468c4762a1bSJed Brown    * this vector which should be zero if the H(div) spaces are implemented
469c4762a1bSJed Brown    * correctly. */
4709566063dSJacob Faibussowitsch   PetscCall(VecGetSubVector(computed,fieldIS[2],&divErr));
4719566063dSJacob Faibussowitsch   PetscCall(VecNorm(divErr,NORM_2,&divErrNorm));
4729566063dSJacob Faibussowitsch   PetscCall(VecRestoreSubVector(computed,fieldIS[2],&divErr));
473c4762a1bSJed Brown   exampleSuccess = (PetscBool)(divErrNorm <= errTol);
474c4762a1bSJed Brown 
4759566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,stdFormat,divErrNorm,exampleSuccess ? "true" : "false"));
476c4762a1bSJed Brown 
477c4762a1bSJed Brown   /* Tear down */
4789566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&divErr));
4799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&computed));
480c4762a1bSJed Brown   for (i = 0; i < 3; ++i) {
4819566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&fieldIS[i]));
482c4762a1bSJed Brown   }
4839566063dSJacob Faibussowitsch   PetscCall(PetscFree(fieldIS));
4849566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
4859566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&mesh));
4869566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
487b122ec5aSJacob Faibussowitsch   return 0;
488c4762a1bSJed Brown }
489c4762a1bSJed Brown 
490c4762a1bSJed Brown /*TEST
491c4762a1bSJed Brown   testset:
492c4762a1bSJed Brown     suffix: 2d_bdm
493c4762a1bSJed Brown     requires: triangle
49430602db0SMatthew G. Knepley     args: -velocity_petscfe_default_quadrature_order 1 \
495c4762a1bSJed Brown       -velocity_petscspace_degree 1 \
496c4762a1bSJed Brown       -velocity_petscdualspace_type bdm \
497c4762a1bSJed Brown       -divErr_petscspace_degree 1 \
498c4762a1bSJed Brown       -divErr_petscdualspace_lagrange_continuity false \
499c4762a1bSJed Brown       -snes_error_if_not_converged \
500c4762a1bSJed Brown       -ksp_rtol 1e-10 \
501c4762a1bSJed Brown       -ksp_error_if_not_converged \
502c4762a1bSJed Brown       -pc_type fieldsplit\
503c4762a1bSJed Brown       -pc_fieldsplit_detect_saddle_point\
504c4762a1bSJed Brown       -pc_fieldsplit_type schur\
505c4762a1bSJed Brown       -pc_fieldsplit_schur_precondition full
506c4762a1bSJed Brown     test:
507c4762a1bSJed Brown       suffix: linear
508c4762a1bSJed Brown       args: -sol_form linear -mesh_transform none
509c4762a1bSJed Brown     test:
510c4762a1bSJed Brown       suffix: sinusoidal
511c4762a1bSJed Brown       args: -sol_form sinusoidal -mesh_transform none
512c4762a1bSJed Brown     test:
513c4762a1bSJed Brown       suffix: sinusoidal_skew
514c4762a1bSJed Brown       args: -sol_form sinusoidal -mesh_transform skew
515c4762a1bSJed Brown     test:
516c4762a1bSJed Brown       suffix: sinusoidal_perturb
517c4762a1bSJed Brown       args: -sol_form sinusoidal -mesh_transform perturb
518c4762a1bSJed Brown     test:
519c4762a1bSJed Brown       suffix: sinusoidal_skew_perturb
520c4762a1bSJed Brown       args: -sol_form sinusoidal -mesh_transform skew_perturb
521c4762a1bSJed Brown 
522c4762a1bSJed Brown   testset:
523c4762a1bSJed Brown     TODO: broken
524c4762a1bSJed Brown     suffix: 2d_bdmq
52530602db0SMatthew G. Knepley     args: -dm_plex_simplex false \
526c4762a1bSJed Brown       -velocity_petscspace_degree 1 \
527c4762a1bSJed Brown       -velocity_petscdualspace_type bdm \
528d092c84bSBrandon Whitchurch       -velocity_petscdualspace_lagrange_tensor 1 \
529c4762a1bSJed Brown       -divErr_petscspace_degree 1 \
530c4762a1bSJed Brown       -divErr_petscdualspace_lagrange_continuity false \
531c4762a1bSJed Brown       -snes_error_if_not_converged \
532c4762a1bSJed Brown       -ksp_rtol 1e-10 \
533c4762a1bSJed Brown       -ksp_error_if_not_converged \
534c4762a1bSJed Brown       -pc_type fieldsplit\
535c4762a1bSJed Brown       -pc_fieldsplit_detect_saddle_point\
536c4762a1bSJed Brown       -pc_fieldsplit_type schur\
537c4762a1bSJed Brown       -pc_fieldsplit_schur_precondition full
538c4762a1bSJed Brown     test:
539c4762a1bSJed Brown       suffix: linear
540c4762a1bSJed Brown       args: -sol_form linear -mesh_transform none
541c4762a1bSJed Brown     test:
542c4762a1bSJed Brown       suffix: sinusoidal
543c4762a1bSJed Brown       args: -sol_form sinusoidal -mesh_transform none
544c4762a1bSJed Brown     test:
545c4762a1bSJed Brown       suffix: sinusoidal_skew
546c4762a1bSJed Brown       args: -sol_form sinusoidal -mesh_transform skew
547c4762a1bSJed Brown     test:
548c4762a1bSJed Brown       suffix: sinusoidal_perturb
549c4762a1bSJed Brown       args: -sol_form sinusoidal -mesh_transform perturb
550c4762a1bSJed Brown     test:
551c4762a1bSJed Brown       suffix: sinusoidal_skew_perturb
552c4762a1bSJed Brown       args: -sol_form sinusoidal -mesh_transform skew_perturb
553c4762a1bSJed Brown 
554c4762a1bSJed Brown   testset:
555c4762a1bSJed Brown     suffix: 3d_bdm
556c4762a1bSJed Brown     requires: ctetgen
55730602db0SMatthew G. Knepley     args: -dm_plex_dim 3 \
558c4762a1bSJed Brown       -velocity_petscspace_degree 1 \
559c4762a1bSJed Brown       -velocity_petscdualspace_type bdm \
560c4762a1bSJed Brown       -divErr_petscspace_degree 1 \
561c4762a1bSJed Brown       -divErr_petscdualspace_lagrange_continuity false \
562c4762a1bSJed Brown       -snes_error_if_not_converged \
563c4762a1bSJed Brown       -ksp_rtol 1e-10 \
564c4762a1bSJed Brown       -ksp_error_if_not_converged \
565c4762a1bSJed Brown       -pc_type fieldsplit \
566c4762a1bSJed Brown       -pc_fieldsplit_detect_saddle_point \
567c4762a1bSJed Brown       -pc_fieldsplit_type schur \
568c4762a1bSJed Brown       -pc_fieldsplit_schur_precondition full
569c4762a1bSJed Brown     test:
570c4762a1bSJed Brown       suffix: linear
571c4762a1bSJed Brown       args: -sol_form linear -mesh_transform none
572c4762a1bSJed Brown     test:
573c4762a1bSJed Brown       suffix: sinusoidal
574c4762a1bSJed Brown       args: -sol_form sinusoidal -mesh_transform none
575c4762a1bSJed Brown     test:
576c4762a1bSJed Brown       suffix: sinusoidal_skew
577c4762a1bSJed Brown       args: -sol_form sinusoidal -mesh_transform skew
578c4762a1bSJed Brown     test:
579c4762a1bSJed Brown       suffix: sinusoidal_perturb
580c4762a1bSJed Brown       args: -sol_form sinusoidal -mesh_transform perturb
581c4762a1bSJed Brown     test:
582c4762a1bSJed Brown       suffix: sinusoidal_skew_perturb
583c4762a1bSJed Brown       args: -sol_form sinusoidal -mesh_transform skew_perturb
584c4762a1bSJed Brown 
585c4762a1bSJed Brown   testset:
586c4762a1bSJed Brown     TODO: broken
587c4762a1bSJed Brown     suffix: 3d_bdmq
588c4762a1bSJed Brown     requires: ctetgen
58930602db0SMatthew G. Knepley     args: -dm_plex_dim 3 \
59030602db0SMatthew G. Knepley       -dm_plex_simplex false \
591c4762a1bSJed Brown       -velocity_petscspace_degree 1 \
592c4762a1bSJed Brown       -velocity_petscdualspace_type bdm \
593d092c84bSBrandon Whitchurch       -velocity_petscdualspace_lagrange_tensor 1 \
594c4762a1bSJed Brown       -divErr_petscspace_degree 1 \
595c4762a1bSJed Brown       -divErr_petscdualspace_lagrange_continuity false \
596c4762a1bSJed Brown       -snes_error_if_not_converged \
597c4762a1bSJed Brown       -ksp_rtol 1e-10 \
598c4762a1bSJed Brown       -ksp_error_if_not_converged \
599c4762a1bSJed Brown       -pc_type fieldsplit \
600c4762a1bSJed Brown       -pc_fieldsplit_detect_saddle_point \
601c4762a1bSJed Brown       -pc_fieldsplit_type schur \
602c4762a1bSJed Brown       -pc_fieldsplit_schur_precondition full
603c4762a1bSJed Brown     test:
604c4762a1bSJed Brown       suffix: linear
605c4762a1bSJed Brown       args: -sol_form linear -mesh_transform none
606c4762a1bSJed Brown     test:
607c4762a1bSJed Brown       suffix: sinusoidal
608c4762a1bSJed Brown       args: -sol_form sinusoidal -mesh_transform none
609c4762a1bSJed Brown     test:
610c4762a1bSJed Brown       suffix: sinusoidal_skew
611c4762a1bSJed Brown       args: -sol_form sinusoidal -mesh_transform skew
612c4762a1bSJed Brown     test:
613c4762a1bSJed Brown       suffix: sinusoidal_perturb
614c4762a1bSJed Brown       args: -sol_form sinusoidal -mesh_transform perturb
615c4762a1bSJed Brown     test:
616c4762a1bSJed Brown       suffix: sinusoidal_skew_perturb
617c4762a1bSJed Brown       args: -sol_form sinusoidal -mesh_transform skew_perturb
618d092c84bSBrandon Whitchurch 
619d092c84bSBrandon Whitchurch   test:
620d092c84bSBrandon Whitchurch     suffix: quad_rt_0
62130602db0SMatthew G. Knepley     args: -dm_plex_simplex false -mesh_transform skew \
622d092c84bSBrandon Whitchurch           -divErr_petscspace_degree 1 \
623d092c84bSBrandon Whitchurch           -divErr_petscdualspace_lagrange_continuity false \
624d092c84bSBrandon Whitchurch           -snes_error_if_not_converged \
625d092c84bSBrandon Whitchurch           -ksp_rtol 1e-10 \
626d092c84bSBrandon Whitchurch           -ksp_error_if_not_converged \
627d092c84bSBrandon Whitchurch           -pc_type fieldsplit\
628d092c84bSBrandon Whitchurch           -pc_fieldsplit_detect_saddle_point\
629d092c84bSBrandon Whitchurch           -pc_fieldsplit_type schur\
630d092c84bSBrandon Whitchurch           -pc_fieldsplit_schur_precondition full \
631d092c84bSBrandon Whitchurch           -velocity_petscfe_default_quadrature_order 1 \
632d092c84bSBrandon Whitchurch           -velocity_petscspace_type sum \
633d092c84bSBrandon Whitchurch           -velocity_petscspace_variables 2 \
634d092c84bSBrandon Whitchurch           -velocity_petscspace_components 2 \
635d092c84bSBrandon Whitchurch           -velocity_petscspace_sum_spaces 2 \
636d092c84bSBrandon Whitchurch           -velocity_petscspace_sum_concatenate true \
637417c287bSToby Isaac           -velocity_sumcomp_0_petscspace_variables 2 \
638417c287bSToby Isaac           -velocity_sumcomp_0_petscspace_type tensor \
639417c287bSToby Isaac           -velocity_sumcomp_0_petscspace_tensor_spaces 2 \
640417c287bSToby Isaac           -velocity_sumcomp_0_petscspace_tensor_uniform false \
641417c287bSToby Isaac           -velocity_sumcomp_0_tensorcomp_0_petscspace_degree 1 \
642417c287bSToby Isaac           -velocity_sumcomp_0_tensorcomp_1_petscspace_degree 0 \
643417c287bSToby Isaac           -velocity_sumcomp_1_petscspace_variables 2 \
644417c287bSToby Isaac           -velocity_sumcomp_1_petscspace_type tensor \
645417c287bSToby Isaac           -velocity_sumcomp_1_petscspace_tensor_spaces 2 \
646417c287bSToby Isaac           -velocity_sumcomp_1_petscspace_tensor_uniform false \
647417c287bSToby Isaac           -velocity_sumcomp_1_tensorcomp_0_petscspace_degree 0 \
648417c287bSToby Isaac           -velocity_sumcomp_1_tensorcomp_1_petscspace_degree 1 \
649d092c84bSBrandon Whitchurch           -velocity_petscdualspace_form_degree -1 \
650d092c84bSBrandon Whitchurch           -velocity_petscdualspace_order 1 \
651d092c84bSBrandon Whitchurch           -velocity_petscdualspace_lagrange_trimmed true
652c4762a1bSJed Brown TEST*/
653