xref: /petsc/src/snes/tests/ex8.c (revision 557cf1953a2187a2a87a3adbc9420a6cc217e91d)
1*557cf195SMatthew G. Knepley static char help[] = "Test adaptive interpolation of functions of a given polynomial order\n\n";
2*557cf195SMatthew G. Knepley 
3*557cf195SMatthew G. Knepley #include <petscdmplex.h>
4*557cf195SMatthew G. Knepley #include <petscsnes.h>
5*557cf195SMatthew G. Knepley 
6*557cf195SMatthew G. Knepley /*
7*557cf195SMatthew G. Knepley   What properties does the adapted interpolator have?
8*557cf195SMatthew G. Knepley 
9*557cf195SMatthew G. Knepley 1) If we adapt to quadratics, we can get lower interpolation error for quadratics (than local interpolation) when using a linear basis
10*557cf195SMatthew G. Knepley 
11*557cf195SMatthew G. Knepley $ ./ex8 -dm_refine 2 -petscspace_degree 1 -qorder 1 -dim 2 -porder 2 -K 2 -num_comp 1 -use_poly 1
12*557cf195SMatthew G. Knepley Function tests FAIL for order 2 at tolerance 1e-10 error 0.00273757
13*557cf195SMatthew G. Knepley Function tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.0721688
14*557cf195SMatthew G. Knepley Interpolation tests FAIL for order 2 at tolerance 1e-10 error 0.00284555
15*557cf195SMatthew G. Knepley Interpolation tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.0721688
16*557cf195SMatthew G. Knepley  Adapting interpolator using polynomials
17*557cf195SMatthew G. Knepley The number of input vectors 4 < 7 the maximum number of column entries
18*557cf195SMatthew G. Knepley   Interpolation poly tests FAIL for order 2 at tolerance 1e-10 error 0.00659864
19*557cf195SMatthew G. Knepley   Interpolation poly tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.0836582
20*557cf195SMatthew G. Knepley   Interpolation trig (0, 0) tests FAIL for order 1 at tolerance 1e-10 error 0.476194
21*557cf195SMatthew G. Knepley   Interpolation trig (0, 0) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22144
22*557cf195SMatthew G. Knepley   Interpolation trig (0, 1) tests FAIL for order 1 at tolerance 1e-10 error 1.39768
23*557cf195SMatthew G. Knepley   Interpolation trig (0, 1) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22144
24*557cf195SMatthew G. Knepley   Interpolation trig (1, 0) tests FAIL for order 2 at tolerance 1e-10 error 1.07315
25*557cf195SMatthew G. Knepley   Interpolation trig (1, 0) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55403
26*557cf195SMatthew G. Knepley   Interpolation trig (1, 1) tests FAIL for order 2 at tolerance 1e-10 error 1.07315
27*557cf195SMatthew G. Knepley   Interpolation trig (1, 1) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55403
28*557cf195SMatthew G. Knepley 
29*557cf195SMatthew G. Knepley $ ./ex8 -dm_refine 2 -petscspace_degree 1 -qorder 1 -dim 2 -porder 2 -K 3 -num_comp 1 -use_poly 1
30*557cf195SMatthew G. Knepley Function tests FAIL for order 2 at tolerance 1e-10 error 0.00273757
31*557cf195SMatthew G. Knepley Function tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.0721688
32*557cf195SMatthew G. Knepley Interpolation tests FAIL for order 2 at tolerance 1e-10 error 0.00284555
33*557cf195SMatthew G. Knepley Interpolation tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.0721688
34*557cf195SMatthew G. Knepley  Adapting interpolator using polynomials
35*557cf195SMatthew G. Knepley The number of input vectors 6 < 7 the maximum number of column entries
36*557cf195SMatthew G. Knepley   Interpolation poly tests FAIL for order 2 at tolerance 1e-10 error 0.00194055
37*557cf195SMatthew G. Knepley   Interpolation poly tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.0525591
38*557cf195SMatthew G. Knepley   Interpolation trig (0, 0) tests FAIL for order 1 at tolerance 1e-10 error 0.476255
39*557cf195SMatthew G. Knepley   Interpolation trig (0, 0) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22132
40*557cf195SMatthew G. Knepley   Interpolation trig (0, 1) tests FAIL for order 1 at tolerance 1e-10 error 1.39785
41*557cf195SMatthew G. Knepley   Interpolation trig (0, 1) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22119
42*557cf195SMatthew G. Knepley   Interpolation trig (1, 0) tests FAIL for order 2 at tolerance 1e-10 error 1.0727
43*557cf195SMatthew G. Knepley   Interpolation trig (1, 0) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55364
44*557cf195SMatthew G. Knepley   Interpolation trig (1, 1) tests FAIL for order 2 at tolerance 1e-10 error 1.0727
45*557cf195SMatthew G. Knepley   Interpolation trig (1, 1) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55364
46*557cf195SMatthew G. Knepley   Interpolation trig (2, 0) tests FAIL for order 3 at tolerance 1e-10 error 0.705258
47*557cf195SMatthew G. Knepley   Interpolation trig (2, 0) tests FAIL for order 3 derivatives at tolerance 1e-10 error 6.82037
48*557cf195SMatthew G. Knepley   Interpolation trig (2, 1) tests FAIL for order 3 at tolerance 1e-10 error 0.705258
49*557cf195SMatthew G. Knepley   Interpolation trig (2, 1) tests FAIL for order 3 derivatives at tolerance 1e-10 error 6.82037
50*557cf195SMatthew G. Knepley 
51*557cf195SMatthew G. Knepley 2) We can more accurately capture low harmonics
52*557cf195SMatthew G. Knepley 
53*557cf195SMatthew G. Knepley If we adapt polynomials, we can be exact
54*557cf195SMatthew G. Knepley 
55*557cf195SMatthew G. Knepley $ ./ex8 -dm_refine 2 -petscspace_degree 1 -qorder 1 -dim 2 -porder 1 -K 2 -num_comp 1 -use_poly 1
56*557cf195SMatthew G. Knepley Function tests pass for order 1 at tolerance 1e-10
57*557cf195SMatthew G. Knepley Function tests pass for order 1 derivatives at tolerance 1e-10
58*557cf195SMatthew G. Knepley Interpolation tests pass for order 1 at tolerance 1e-10
59*557cf195SMatthew G. Knepley Interpolation tests pass for order 1 derivatives at tolerance 1e-10
60*557cf195SMatthew G. Knepley  Adapting interpolator using polynomials
61*557cf195SMatthew G. Knepley The number of input vectors 4 < 7 the maximum number of column entries
62*557cf195SMatthew G. Knepley   Interpolation poly tests pass for order 1 at tolerance 1e-10
63*557cf195SMatthew G. Knepley   Interpolation poly tests pass for order 1 derivatives at tolerance 1e-10
64*557cf195SMatthew G. Knepley   Interpolation trig (0, 0) tests FAIL for order 1 at tolerance 1e-10 error 0.476194
65*557cf195SMatthew G. Knepley   Interpolation trig (0, 0) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22144
66*557cf195SMatthew G. Knepley   Interpolation trig (0, 1) tests FAIL for order 1 at tolerance 1e-10 error 1.39768
67*557cf195SMatthew G. Knepley   Interpolation trig (0, 1) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22144
68*557cf195SMatthew G. Knepley   Interpolation trig (1, 0) tests FAIL for order 2 at tolerance 1e-10 error 1.07315
69*557cf195SMatthew G. Knepley   Interpolation trig (1, 0) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55403
70*557cf195SMatthew G. Knepley   Interpolation trig (1, 1) tests FAIL for order 2 at tolerance 1e-10 error 1.07315
71*557cf195SMatthew G. Knepley   Interpolation trig (1, 1) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55403
72*557cf195SMatthew G. Knepley 
73*557cf195SMatthew G. Knepley and least for small K,
74*557cf195SMatthew G. Knepley 
75*557cf195SMatthew G. Knepley $ ./ex8 -dm_refine 2 -petscspace_degree 1 -qorder 1 -dim 2 -porder 1 -K 4 -num_comp 1 -use_poly 1
76*557cf195SMatthew G. Knepley Function tests pass for order 1 at tolerance 1e-10
77*557cf195SMatthew G. Knepley Function tests pass for order 1 derivatives at tolerance 1e-10
78*557cf195SMatthew G. Knepley Interpolation tests pass for order 1 at tolerance 1e-10
79*557cf195SMatthew G. Knepley Interpolation tests pass for order 1 derivatives at tolerance 1e-10
80*557cf195SMatthew G. Knepley  Adapting interpolator using polynomials
81*557cf195SMatthew G. Knepley   Interpolation poly tests FAIL for order 1 at tolerance 1e-10 error 0.0015351
82*557cf195SMatthew G. Knepley   Interpolation poly tests FAIL for order 1 derivatives at tolerance 1e-10 error 0.0427369
83*557cf195SMatthew G. Knepley   Interpolation trig (0, 0) tests FAIL for order 1 at tolerance 1e-10 error 0.476359
84*557cf195SMatthew G. Knepley   Interpolation trig (0, 0) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22115
85*557cf195SMatthew G. Knepley   Interpolation trig (0, 1) tests FAIL for order 1 at tolerance 1e-10 error 1.3981
86*557cf195SMatthew G. Knepley   Interpolation trig (0, 1) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22087
87*557cf195SMatthew G. Knepley   Interpolation trig (1, 0) tests FAIL for order 2 at tolerance 1e-10 error 1.07228
88*557cf195SMatthew G. Knepley   Interpolation trig (1, 0) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55238
89*557cf195SMatthew G. Knepley   Interpolation trig (1, 1) tests FAIL for order 2 at tolerance 1e-10 error 1.07228
90*557cf195SMatthew G. Knepley   Interpolation trig (1, 1) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55238
91*557cf195SMatthew G. Knepley   Interpolation trig (2, 0) tests FAIL for order 3 at tolerance 1e-10 error 0.704947
92*557cf195SMatthew G. Knepley   Interpolation trig (2, 0) tests FAIL for order 3 derivatives at tolerance 1e-10 error 6.82254
93*557cf195SMatthew G. Knepley   Interpolation trig (2, 1) tests FAIL for order 3 at tolerance 1e-10 error 0.704948
94*557cf195SMatthew G. Knepley   Interpolation trig (2, 1) tests FAIL for order 3 derivatives at tolerance 1e-10 error 6.82254
95*557cf195SMatthew G. Knepley   Interpolation trig (3, 0) tests FAIL for order 4 at tolerance 1e-10 error 0.893279
96*557cf195SMatthew G. Knepley   Interpolation trig (3, 0) tests FAIL for order 4 derivatives at tolerance 1e-10 error 8.93718
97*557cf195SMatthew G. Knepley   Interpolation trig (3, 1) tests FAIL for order 4 at tolerance 1e-10 error 0.89328
98*557cf195SMatthew G. Knepley   Interpolation trig (3, 1) tests FAIL for order 4 derivatives at tolerance 1e-10 error 8.93717
99*557cf195SMatthew G. Knepley 
100*557cf195SMatthew G. Knepley but adapting to harmonics gives alright polynomials errors and much better harmonics errors.
101*557cf195SMatthew G. Knepley 
102*557cf195SMatthew G. Knepley $ ./ex8 -dm_refine 2 -petscspace_degree 1 -qorder 1 -dim 2 -porder 1 -K 4 -num_comp 1 -use_poly 0
103*557cf195SMatthew G. Knepley Function tests pass for order 1 at tolerance 1e-10
104*557cf195SMatthew G. Knepley Function tests pass for order 1 derivatives at tolerance 1e-10
105*557cf195SMatthew G. Knepley Interpolation tests pass for order 1 at tolerance 1e-10
106*557cf195SMatthew G. Knepley Interpolation tests pass for order 1 derivatives at tolerance 1e-10
107*557cf195SMatthew G. Knepley  Adapting interpolator using harmonics
108*557cf195SMatthew G. Knepley   Interpolation poly tests FAIL for order 1 at tolerance 1e-10 error 0.0720606
109*557cf195SMatthew G. Knepley   Interpolation poly tests FAIL for order 1 derivatives at tolerance 1e-10 error 1.97779
110*557cf195SMatthew G. Knepley   Interpolation trig (0, 0) tests FAIL for order 1 at tolerance 1e-10 error 0.0398055
111*557cf195SMatthew G. Knepley   Interpolation trig (0, 0) tests FAIL for order 1 derivatives at tolerance 1e-10 error 0.995963
112*557cf195SMatthew G. Knepley   Interpolation trig (0, 1) tests FAIL for order 1 at tolerance 1e-10 error 0.0398051
113*557cf195SMatthew G. Knepley   Interpolation trig (0, 1) tests FAIL for order 1 derivatives at tolerance 1e-10 error 0.995964
114*557cf195SMatthew G. Knepley   Interpolation trig (1, 0) tests FAIL for order 2 at tolerance 1e-10 error 0.0238441
115*557cf195SMatthew G. Knepley   Interpolation trig (1, 0) tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.888611
116*557cf195SMatthew G. Knepley   Interpolation trig (1, 1) tests FAIL for order 2 at tolerance 1e-10 error 0.0238346
117*557cf195SMatthew G. Knepley   Interpolation trig (1, 1) tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.888612
118*557cf195SMatthew G. Knepley   Interpolation trig (2, 0) tests FAIL for order 3 at tolerance 1e-10 error 0.0537968
119*557cf195SMatthew G. Knepley   Interpolation trig (2, 0) tests FAIL for order 3 derivatives at tolerance 1e-10 error 1.57665
120*557cf195SMatthew G. Knepley   Interpolation trig (2, 1) tests FAIL for order 3 at tolerance 1e-10 error 0.0537779
121*557cf195SMatthew G. Knepley   Interpolation trig (2, 1) tests FAIL for order 3 derivatives at tolerance 1e-10 error 1.57666
122*557cf195SMatthew G. Knepley   Interpolation trig (3, 0) tests FAIL for order 4 at tolerance 1e-10 error 0.0775838
123*557cf195SMatthew G. Knepley   Interpolation trig (3, 0) tests FAIL for order 4 derivatives at tolerance 1e-10 error 2.36926
124*557cf195SMatthew G. Knepley   Interpolation trig (3, 1) tests FAIL for order 4 at tolerance 1e-10 error 0.0775464
125*557cf195SMatthew G. Knepley   Interpolation trig (3, 1) tests FAIL for order 4 derivatives at tolerance 1e-10 error 2.36929
126*557cf195SMatthew G. Knepley */
127*557cf195SMatthew G. Knepley 
128*557cf195SMatthew G. Knepley typedef struct {
129*557cf195SMatthew G. Knepley   /* Domain and mesh definition */
130*557cf195SMatthew G. Knepley   PetscInt  dim;               /* The topological mesh dimension */
131*557cf195SMatthew G. Knepley   PetscBool simplex;           /* Flag for simplex or tensor product mesh */
132*557cf195SMatthew G. Knepley   PetscBool interpolate;       /* Generate intermediate mesh elements */
133*557cf195SMatthew G. Knepley   PetscReal refinementLimit;   /* The largest allowable cell volume */
134*557cf195SMatthew G. Knepley   /* Element definition */
135*557cf195SMatthew G. Knepley   PetscInt  qorder;            /* Order of the quadrature */
136*557cf195SMatthew G. Knepley   PetscInt  Nc;                /* Number of field components */
137*557cf195SMatthew G. Knepley   PetscFE   fe;                /* The finite element */
138*557cf195SMatthew G. Knepley   /* Testing space */
139*557cf195SMatthew G. Knepley   PetscInt  porder;            /* Order of polynomials to test */
140*557cf195SMatthew G. Knepley   PetscReal constants[3];      /* Constant values for each dimension */
141*557cf195SMatthew G. Knepley   PetscInt  m;                 /* The frequency of sinusoids to use */
142*557cf195SMatthew G. Knepley   PetscInt  dir;               /* The direction of sinusoids to use */
143*557cf195SMatthew G. Knepley   /* Adaptation */
144*557cf195SMatthew G. Knepley   PetscInt  K;                 /* Number of coarse modes used for optimization */
145*557cf195SMatthew G. Knepley   PetscBool usePoly;           /* Use polynomials, or harmonics, to adapt interpolator */
146*557cf195SMatthew G. Knepley } AppCtx;
147*557cf195SMatthew G. Knepley 
148*557cf195SMatthew G. Knepley typedef enum {INTERPOLATION, RESTRICTION, INJECTION} InterpType;
149*557cf195SMatthew G. Knepley 
150*557cf195SMatthew G. Knepley /* u = 1 */
151*557cf195SMatthew G. Knepley PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
152*557cf195SMatthew G. Knepley {
153*557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
154*557cf195SMatthew G. Knepley   PetscInt d = user->dir;
155*557cf195SMatthew G. Knepley 
156*557cf195SMatthew G. Knepley   if (Nc > 1) {
157*557cf195SMatthew G. Knepley     for (d = 0; d < Nc; ++d) u[d] = user->constants[d];
158*557cf195SMatthew G. Knepley   } else {
159*557cf195SMatthew G. Knepley     u[0] = user->constants[d];
160*557cf195SMatthew G. Knepley   }
161*557cf195SMatthew G. Knepley   return 0;
162*557cf195SMatthew G. Knepley }
163*557cf195SMatthew G. Knepley PetscErrorCode constantDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nc, PetscScalar *u, void *ctx)
164*557cf195SMatthew G. Knepley {
165*557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
166*557cf195SMatthew G. Knepley   PetscInt d = user->dir;
167*557cf195SMatthew G. Knepley 
168*557cf195SMatthew G. Knepley   if (Nc > 1) {
169*557cf195SMatthew G. Knepley     for (d = 0; d < Nc; ++d) u[d] = 0.0;
170*557cf195SMatthew G. Knepley   } else {
171*557cf195SMatthew G. Knepley     u[0] = user->constants[d];
172*557cf195SMatthew G. Knepley   }
173*557cf195SMatthew G. Knepley   return 0;
174*557cf195SMatthew G. Knepley }
175*557cf195SMatthew G. Knepley 
176*557cf195SMatthew G. Knepley /* u = x */
177*557cf195SMatthew G. Knepley PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
178*557cf195SMatthew G. Knepley {
179*557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
180*557cf195SMatthew G. Knepley   PetscInt d = user->dir;
181*557cf195SMatthew G. Knepley 
182*557cf195SMatthew G. Knepley   if (Nc > 1) {
183*557cf195SMatthew G. Knepley     for (d = 0; d < Nc; ++d) u[d] = coords[d];
184*557cf195SMatthew G. Knepley   } else {
185*557cf195SMatthew G. Knepley     u[0] = coords[d];
186*557cf195SMatthew G. Knepley   }
187*557cf195SMatthew G. Knepley   return 0;
188*557cf195SMatthew G. Knepley }
189*557cf195SMatthew G. Knepley PetscErrorCode linearDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nc, PetscScalar *u, void *ctx)
190*557cf195SMatthew G. Knepley {
191*557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
192*557cf195SMatthew G. Knepley   PetscInt d = user->dir;
193*557cf195SMatthew G. Knepley 
194*557cf195SMatthew G. Knepley   if (Nc > 1) {
195*557cf195SMatthew G. Knepley     PetscInt e;
196*557cf195SMatthew G. Knepley     for (d = 0; d < Nc; ++d) {
197*557cf195SMatthew G. Knepley       u[d] = 0.0;
198*557cf195SMatthew G. Knepley       for (e = 0; e < dim; ++e) u[d] += (d == e ? 1.0 : 0.0) * n[e];
199*557cf195SMatthew G. Knepley     }
200*557cf195SMatthew G. Knepley   } else {
201*557cf195SMatthew G. Knepley     u[0] = n[d];
202*557cf195SMatthew G. Knepley   }
203*557cf195SMatthew G. Knepley   return 0;
204*557cf195SMatthew G. Knepley }
205*557cf195SMatthew G. Knepley 
206*557cf195SMatthew G. Knepley /* u = x^2 or u = (x^2, xy) or u = (xy, yz, zx) */
207*557cf195SMatthew G. Knepley PetscErrorCode quadratic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
208*557cf195SMatthew G. Knepley {
209*557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
210*557cf195SMatthew G. Knepley   PetscInt d = user->dir;
211*557cf195SMatthew G. Knepley 
212*557cf195SMatthew G. Knepley   if (Nc > 1) {
213*557cf195SMatthew G. Knepley     if (Nc > 2) {u[0] = coords[0]*coords[1]; u[1] = coords[1]*coords[2]; u[2] = coords[2]*coords[0];}
214*557cf195SMatthew G. Knepley     else        {u[0] = coords[0]*coords[0]; u[1] = coords[0]*coords[1];}
215*557cf195SMatthew G. Knepley   } else {
216*557cf195SMatthew G. Knepley     u[0] = coords[d]*coords[d];
217*557cf195SMatthew G. Knepley   }
218*557cf195SMatthew G. Knepley   return 0;
219*557cf195SMatthew G. Knepley }
220*557cf195SMatthew G. Knepley PetscErrorCode quadraticDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nc, PetscScalar *u, void *ctx)
221*557cf195SMatthew G. Knepley {
222*557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
223*557cf195SMatthew G. Knepley   PetscInt d = user->dir;
224*557cf195SMatthew G. Knepley 
225*557cf195SMatthew G. Knepley   if (Nc > 1) {
226*557cf195SMatthew G. Knepley     if (Nc > 2) {u[0] = coords[1]*n[0] + coords[0]*n[1]; u[1] = coords[2]*n[1] + coords[1]*n[2]; u[2] = coords[2]*n[0] + coords[0]*n[2];}
227*557cf195SMatthew G. Knepley     else        {u[0] = 2.0*coords[0]*n[0]; u[1] = coords[1]*n[0] + coords[0]*n[1];}
228*557cf195SMatthew G. Knepley   } else {
229*557cf195SMatthew G. Knepley     u[0] = 2.0*coords[d]*n[d];
230*557cf195SMatthew G. Knepley   }
231*557cf195SMatthew G. Knepley   return 0;
232*557cf195SMatthew G. Knepley }
233*557cf195SMatthew G. Knepley 
234*557cf195SMatthew G. Knepley /* u = x^3 or u = (x^3, x^2y) or u = (x^2y, y^2z, z^2x) */
235*557cf195SMatthew G. Knepley PetscErrorCode cubic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
236*557cf195SMatthew G. Knepley {
237*557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
238*557cf195SMatthew G. Knepley   PetscInt d = user->dir;
239*557cf195SMatthew G. Knepley 
240*557cf195SMatthew G. Knepley   if (Nc > 1) {
241*557cf195SMatthew G. Knepley     if (Nc > 2) {u[0] = coords[0]*coords[0]*coords[1]; u[1] = coords[1]*coords[1]*coords[2]; u[2] = coords[2]*coords[2]*coords[0];}
242*557cf195SMatthew G. Knepley     else        {u[0] = coords[0]*coords[0]*coords[0]; u[1] = coords[0]*coords[0]*coords[1];}
243*557cf195SMatthew G. Knepley   } else {
244*557cf195SMatthew G. Knepley     u[0] = coords[d]*coords[d]*coords[d];
245*557cf195SMatthew G. Knepley   }
246*557cf195SMatthew G. Knepley   return 0;
247*557cf195SMatthew G. Knepley }
248*557cf195SMatthew G. Knepley PetscErrorCode cubicDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nc, PetscScalar *u, void *ctx)
249*557cf195SMatthew G. Knepley {
250*557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
251*557cf195SMatthew G. Knepley   PetscInt d = user->dir;
252*557cf195SMatthew G. Knepley 
253*557cf195SMatthew G. Knepley   if (Nc > 1) {
254*557cf195SMatthew G. Knepley     if (Nc > 2) {u[0] = 2.0*coords[0]*coords[1]*n[0] + coords[0]*coords[0]*n[1]; u[1] = 2.0*coords[1]*coords[2]*n[1] + coords[1]*coords[1]*n[2]; u[2] = 2.0*coords[2]*coords[0]*n[2] + coords[2]*coords[2]*n[0];}
255*557cf195SMatthew G. Knepley     else        {u[0] = 3.0*coords[0]*coords[0]*n[0]; u[1] = 2.0*coords[0]*coords[1]*n[0] + coords[0]*coords[0]*n[1];}
256*557cf195SMatthew G. Knepley   } else {
257*557cf195SMatthew G. Knepley     u[0] = 3.0*coords[d]*coords[d]*n[d];
258*557cf195SMatthew G. Knepley   }
259*557cf195SMatthew G. Knepley   return 0;
260*557cf195SMatthew G. Knepley }
261*557cf195SMatthew G. Knepley 
262*557cf195SMatthew G. Knepley /* u = x^4 or u = (x^4, x^2y^2) or u = (x^2y^2, y^2z^2, z^2x^2) */
263*557cf195SMatthew G. Knepley PetscErrorCode quartic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
264*557cf195SMatthew G. Knepley {
265*557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
266*557cf195SMatthew G. Knepley   PetscInt d = user->dir;
267*557cf195SMatthew G. Knepley 
268*557cf195SMatthew G. Knepley   if (Nc > 1) {
269*557cf195SMatthew G. Knepley     if (Nc > 2) {u[0] = coords[0]*coords[0]*coords[1]*coords[1]; u[1] = coords[1]*coords[1]*coords[2]*coords[2]; u[2] = coords[2]*coords[2]*coords[0]*coords[0];}
270*557cf195SMatthew G. Knepley     else        {u[0] = coords[0]*coords[0]*coords[0]*coords[0]; u[1] = coords[0]*coords[0]*coords[1]*coords[1];}
271*557cf195SMatthew G. Knepley   } else {
272*557cf195SMatthew G. Knepley     u[0] = coords[d]*coords[d]*coords[d]*coords[d];
273*557cf195SMatthew G. Knepley   }
274*557cf195SMatthew G. Knepley   return 0;
275*557cf195SMatthew G. Knepley }
276*557cf195SMatthew G. Knepley PetscErrorCode quarticDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nc, PetscScalar *u, void *ctx)
277*557cf195SMatthew G. Knepley {
278*557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
279*557cf195SMatthew G. Knepley   PetscInt d = user->dir;
280*557cf195SMatthew G. Knepley 
281*557cf195SMatthew G. Knepley   if (Nc > 1) {
282*557cf195SMatthew G. Knepley     if (Nc > 2) {u[0] = 2.0*coords[0]*coords[1]*coords[1]*n[0] + 2.0*coords[0]*coords[0]*coords[1]*n[1];
283*557cf195SMatthew G. Knepley                  u[1] = 2.0*coords[1]*coords[2]*coords[2]*n[1] + 2.0*coords[1]*coords[1]*coords[2]*n[2];
284*557cf195SMatthew G. Knepley                  u[2] = 2.0*coords[2]*coords[0]*coords[0]*n[2] + 2.0*coords[2]*coords[2]*coords[0]*n[0];}
285*557cf195SMatthew G. Knepley     else        {u[0] = 4.0*coords[0]*coords[0]*coords[0]*n[0]; u[1] = 2.0*coords[0]*coords[1]*coords[1]*n[0] + 2.0*coords[0]*coords[0]*coords[1]*n[1];}
286*557cf195SMatthew G. Knepley   } else {
287*557cf195SMatthew G. Knepley     u[0] = 4.0*coords[d]*coords[d]*coords[d]*n[d];
288*557cf195SMatthew G. Knepley   }
289*557cf195SMatthew G. Knepley   return 0;
290*557cf195SMatthew G. Knepley }
291*557cf195SMatthew G. Knepley 
292*557cf195SMatthew G. Knepley PetscErrorCode mytanh(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
293*557cf195SMatthew G. Knepley {
294*557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
295*557cf195SMatthew G. Knepley   PetscInt d = user->dir;
296*557cf195SMatthew G. Knepley 
297*557cf195SMatthew G. Knepley   if (Nc > 1) {
298*557cf195SMatthew G. Knepley     for (d = 0; d < Nc; ++d) u[d] = PetscTanhReal(coords[d] - 0.5);
299*557cf195SMatthew G. Knepley   } else {
300*557cf195SMatthew G. Knepley     u[0] = PetscTanhReal(coords[d] - 0.5);
301*557cf195SMatthew G. Knepley   }
302*557cf195SMatthew G. Knepley   return 0;
303*557cf195SMatthew G. Knepley }
304*557cf195SMatthew G. Knepley PetscErrorCode mytanhDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nc, PetscScalar *u, void *ctx)
305*557cf195SMatthew G. Knepley {
306*557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
307*557cf195SMatthew G. Knepley   PetscInt d = user->dir;
308*557cf195SMatthew G. Knepley 
309*557cf195SMatthew G. Knepley   if (Nc > 1) {
310*557cf195SMatthew G. Knepley     for (d = 0; d < Nc; ++d) u[d] = 1.0/PetscSqr(PetscCoshReal(coords[d] - 0.5)) * n[d];
311*557cf195SMatthew G. Knepley   } else {
312*557cf195SMatthew G. Knepley     u[0] = 1.0/PetscSqr(PetscCoshReal(coords[d] - 0.5)) * n[d];
313*557cf195SMatthew G. Knepley   }
314*557cf195SMatthew G. Knepley   return 0;
315*557cf195SMatthew G. Knepley }
316*557cf195SMatthew G. Knepley 
317*557cf195SMatthew G. Knepley PetscErrorCode trig(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
318*557cf195SMatthew G. Knepley {
319*557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
320*557cf195SMatthew G. Knepley   PetscInt m = user->m, d = user->dir;
321*557cf195SMatthew G. Knepley 
322*557cf195SMatthew G. Knepley   if (Nc > 1) {
323*557cf195SMatthew G. Knepley     for (d = 0; d < Nc; ++d) u[d] = PetscSinReal(PETSC_PI*m*coords[d]);
324*557cf195SMatthew G. Knepley   } else {
325*557cf195SMatthew G. Knepley     u[0] = PetscSinReal(PETSC_PI*m*coords[d]);
326*557cf195SMatthew G. Knepley   }
327*557cf195SMatthew G. Knepley   return 0;
328*557cf195SMatthew G. Knepley }
329*557cf195SMatthew G. Knepley PetscErrorCode trigDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nc, PetscScalar *u, void *ctx)
330*557cf195SMatthew G. Knepley {
331*557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
332*557cf195SMatthew G. Knepley   PetscInt m = user->m, d = user->dir;
333*557cf195SMatthew G. Knepley 
334*557cf195SMatthew G. Knepley   if (Nc > 1) {
335*557cf195SMatthew G. Knepley     for (d = 0; d < Nc; ++d) u[d] = PETSC_PI*m*PetscCosReal(PETSC_PI*m*coords[d]) * n[d];
336*557cf195SMatthew G. Knepley   } else {
337*557cf195SMatthew G. Knepley     u[0] = PETSC_PI*m*PetscCosReal(PETSC_PI*m*coords[d]) * n[d];
338*557cf195SMatthew G. Knepley   }
339*557cf195SMatthew G. Knepley   return 0;
340*557cf195SMatthew G. Knepley }
341*557cf195SMatthew G. Knepley 
342*557cf195SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
343*557cf195SMatthew G. Knepley {
344*557cf195SMatthew G. Knepley   PetscErrorCode ierr;
345*557cf195SMatthew G. Knepley 
346*557cf195SMatthew G. Knepley   PetscFunctionBeginUser;
347*557cf195SMatthew G. Knepley   options->dim             = 2;
348*557cf195SMatthew G. Knepley   options->simplex         = PETSC_TRUE;
349*557cf195SMatthew G. Knepley   options->interpolate     = PETSC_TRUE;
350*557cf195SMatthew G. Knepley   options->refinementLimit = 0.0;
351*557cf195SMatthew G. Knepley   options->qorder          = 0;
352*557cf195SMatthew G. Knepley   options->Nc              = PETSC_DEFAULT;
353*557cf195SMatthew G. Knepley   options->porder          = 0;
354*557cf195SMatthew G. Knepley   options->m               = 1;
355*557cf195SMatthew G. Knepley   options->dir             = 0;
356*557cf195SMatthew G. Knepley   options->K               = 0;
357*557cf195SMatthew G. Knepley   options->usePoly         = PETSC_TRUE;
358*557cf195SMatthew G. Knepley 
359*557cf195SMatthew G. Knepley   ierr = PetscOptionsBegin(comm, "", "Projection Test Options", "DMPlex");CHKERRQ(ierr);
360*557cf195SMatthew G. Knepley   ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex8.c", options->dim, &options->dim, NULL);CHKERRQ(ierr);
361*557cf195SMatthew G. Knepley   ierr = PetscOptionsBool("-simplex", "Flag for simplices or hexhedra", "ex8.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr);
362*557cf195SMatthew G. Knepley   ierr = PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex8.c", options->interpolate, &options->interpolate, NULL);CHKERRQ(ierr);
363*557cf195SMatthew G. Knepley   ierr = PetscOptionsReal("-refinement_limit", "The largest allowable cell volume", "ex8.c", options->refinementLimit, &options->refinementLimit, NULL);CHKERRQ(ierr);
364*557cf195SMatthew G. Knepley   ierr = PetscOptionsInt("-qorder", "The quadrature order", "ex8.c", options->qorder, &options->qorder, NULL);CHKERRQ(ierr);
365*557cf195SMatthew G. Knepley   ierr = PetscOptionsInt("-num_comp", "The number of field components", "ex8.c", options->Nc, &options->Nc, NULL);CHKERRQ(ierr);
366*557cf195SMatthew G. Knepley   ierr = PetscOptionsInt("-porder", "The order of polynomials to test", "ex8.c", options->porder, &options->porder, NULL);CHKERRQ(ierr);
367*557cf195SMatthew G. Knepley   ierr = PetscOptionsInt("-K", "The number of coarse modes used in optimization", "ex8.c", options->K, &options->K, NULL);CHKERRQ(ierr);
368*557cf195SMatthew G. Knepley   ierr = PetscOptionsBool("-use_poly", "Use polynomials (or harmonics) to adapt interpolator", "ex8.c", options->usePoly, &options->usePoly, NULL);CHKERRQ(ierr);
369*557cf195SMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
370*557cf195SMatthew G. Knepley 
371*557cf195SMatthew G. Knepley   options->Nc = options->Nc < 0 ? options->dim : options->Nc;
372*557cf195SMatthew G. Knepley 
373*557cf195SMatthew G. Knepley   PetscFunctionReturn(0);
374*557cf195SMatthew G. Knepley }
375*557cf195SMatthew G. Knepley 
376*557cf195SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
377*557cf195SMatthew G. Knepley {
378*557cf195SMatthew G. Knepley   DM               pdm      = NULL;
379*557cf195SMatthew G. Knepley   PetscInt         cells[3] = {2, 2, 2};
380*557cf195SMatthew G. Knepley   PetscPartitioner part;
381*557cf195SMatthew G. Knepley   PetscErrorCode   ierr;
382*557cf195SMatthew G. Knepley 
383*557cf195SMatthew G. Knepley   PetscFunctionBeginUser;
384*557cf195SMatthew G. Knepley   ierr = PetscOptionsGetInt(NULL, NULL, "-da_grid_x", &cells[0], NULL);CHKERRQ(ierr);
385*557cf195SMatthew G. Knepley   ierr = PetscOptionsGetInt(NULL, NULL, "-da_grid_y", &cells[1], NULL);CHKERRQ(ierr);
386*557cf195SMatthew G. Knepley   ierr = PetscOptionsGetInt(NULL, NULL, "-da_grid_z", &cells[2], NULL);CHKERRQ(ierr);
387*557cf195SMatthew G. Knepley   ierr = DMPlexCreateBoxMesh(comm, user->dim, user->simplex, cells, NULL, NULL, NULL, user->interpolate, dm);CHKERRQ(ierr);
388*557cf195SMatthew G. Knepley   /* Refine mesh using a volume constraint */
389*557cf195SMatthew G. Knepley   if (user->simplex) {
390*557cf195SMatthew G. Knepley     DM rdm = NULL;
391*557cf195SMatthew G. Knepley 
392*557cf195SMatthew G. Knepley     ierr = DMPlexSetRefinementLimit(*dm, user->refinementLimit);CHKERRQ(ierr);
393*557cf195SMatthew G. Knepley     ierr = DMRefine(*dm, comm, &rdm);CHKERRQ(ierr);
394*557cf195SMatthew G. Knepley     if (rdm) {
395*557cf195SMatthew G. Knepley       ierr = DMDestroy(dm);CHKERRQ(ierr);
396*557cf195SMatthew G. Knepley       *dm  = rdm;
397*557cf195SMatthew G. Knepley     }
398*557cf195SMatthew G. Knepley   }
399*557cf195SMatthew G. Knepley   ierr = DMPlexSetRefinementUniform(*dm, PETSC_TRUE);CHKERRQ(ierr);
400*557cf195SMatthew G. Knepley   /* Distribute mesh over processes */
401*557cf195SMatthew G. Knepley   ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr);
402*557cf195SMatthew G. Knepley   ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
403*557cf195SMatthew G. Knepley   ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr);
404*557cf195SMatthew G. Knepley   if (pdm) {
405*557cf195SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
406*557cf195SMatthew G. Knepley     *dm  = pdm;
407*557cf195SMatthew G. Knepley   }
408*557cf195SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) *dm, user->simplex ? "Simplicial Mesh" : "Hexahedral Mesh");CHKERRQ(ierr);
409*557cf195SMatthew G. Knepley   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
410*557cf195SMatthew G. Knepley   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
411*557cf195SMatthew G. Knepley   PetscFunctionReturn(0);
412*557cf195SMatthew G. Knepley }
413*557cf195SMatthew G. Knepley 
414*557cf195SMatthew G. Knepley /* Setup functions to approximate */
415*557cf195SMatthew G. Knepley static PetscErrorCode SetupFunctions(DM dm, PetscBool usePoly, PetscInt order, PetscInt dir, PetscErrorCode (**exactFuncs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *),
416*557cf195SMatthew G. Knepley                                      PetscErrorCode (**exactFuncDers)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), AppCtx *user)
417*557cf195SMatthew G. Knepley {
418*557cf195SMatthew G. Knepley   PetscInt       dim;
419*557cf195SMatthew G. Knepley   PetscErrorCode ierr;
420*557cf195SMatthew G. Knepley 
421*557cf195SMatthew G. Knepley   PetscFunctionBeginUser;
422*557cf195SMatthew G. Knepley   user->dir = dir;
423*557cf195SMatthew G. Knepley   if (usePoly) {
424*557cf195SMatthew G. Knepley     switch (order) {
425*557cf195SMatthew G. Knepley     case 0:
426*557cf195SMatthew G. Knepley       exactFuncs[0]    = constant;
427*557cf195SMatthew G. Knepley       exactFuncDers[0] = constantDer;
428*557cf195SMatthew G. Knepley       break;
429*557cf195SMatthew G. Knepley     case 1:
430*557cf195SMatthew G. Knepley       exactFuncs[0]    = linear;
431*557cf195SMatthew G. Knepley       exactFuncDers[0] = linearDer;
432*557cf195SMatthew G. Knepley       break;
433*557cf195SMatthew G. Knepley     case 2:
434*557cf195SMatthew G. Knepley       exactFuncs[0]    = quadratic;
435*557cf195SMatthew G. Knepley       exactFuncDers[0] = quadraticDer;
436*557cf195SMatthew G. Knepley       break;
437*557cf195SMatthew G. Knepley     case 3:
438*557cf195SMatthew G. Knepley       exactFuncs[0]    = cubic;
439*557cf195SMatthew G. Knepley       exactFuncDers[0] = cubicDer;
440*557cf195SMatthew G. Knepley       break;
441*557cf195SMatthew G. Knepley     case 4:
442*557cf195SMatthew G. Knepley       exactFuncs[0]    = quartic;
443*557cf195SMatthew G. Knepley       exactFuncDers[0] = quarticDer;
444*557cf195SMatthew G. Knepley       break;
445*557cf195SMatthew G. Knepley     default:
446*557cf195SMatthew G. Knepley       ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
447*557cf195SMatthew G. Knepley       SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Could not determine functions to test for dimension %d order %d", dim, order);
448*557cf195SMatthew G. Knepley     }
449*557cf195SMatthew G. Knepley   } else {
450*557cf195SMatthew G. Knepley     user->m          = order;
451*557cf195SMatthew G. Knepley     exactFuncs[0]    = trig;
452*557cf195SMatthew G. Knepley     exactFuncDers[0] = trigDer;
453*557cf195SMatthew G. Knepley   }
454*557cf195SMatthew G. Knepley   PetscFunctionReturn(0);
455*557cf195SMatthew G. Knepley }
456*557cf195SMatthew G. Knepley 
457*557cf195SMatthew G. Knepley static PetscErrorCode ComputeError(DM dm, PetscErrorCode (**exactFuncs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *),
458*557cf195SMatthew G. Knepley                                    PetscErrorCode (**exactFuncDers)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *),
459*557cf195SMatthew G. Knepley                                    void **exactCtxs, PetscReal *error, PetscReal *errorDer, AppCtx *user)
460*557cf195SMatthew G. Knepley {
461*557cf195SMatthew G. Knepley   Vec            u;
462*557cf195SMatthew G. Knepley   PetscReal      n[3] = {1.0, 1.0, 1.0};
463*557cf195SMatthew G. Knepley   PetscErrorCode ierr;
464*557cf195SMatthew G. Knepley 
465*557cf195SMatthew G. Knepley   PetscFunctionBeginUser;
466*557cf195SMatthew G. Knepley   ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr);
467*557cf195SMatthew G. Knepley   /* Project function into FE function space */
468*557cf195SMatthew G. Knepley   ierr = DMProjectFunction(dm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, u);CHKERRQ(ierr);
469*557cf195SMatthew G. Knepley   ierr = VecViewFromOptions(u, NULL, "-projection_view");CHKERRQ(ierr);
470*557cf195SMatthew G. Knepley   /* Compare approximation to exact in L_2 */
471*557cf195SMatthew G. Knepley   ierr = DMComputeL2Diff(dm, 0.0, exactFuncs, exactCtxs, u, error);CHKERRQ(ierr);
472*557cf195SMatthew G. Knepley   ierr = DMComputeL2GradientDiff(dm, 0.0, exactFuncDers, exactCtxs, u, n, errorDer);CHKERRQ(ierr);
473*557cf195SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr);
474*557cf195SMatthew G. Knepley   PetscFunctionReturn(0);
475*557cf195SMatthew G. Knepley }
476*557cf195SMatthew G. Knepley 
477*557cf195SMatthew G. Knepley static PetscErrorCode CheckFunctions(DM dm, PetscInt order, AppCtx *user)
478*557cf195SMatthew G. Knepley {
479*557cf195SMatthew G. Knepley   PetscErrorCode (*exactFuncs[1]) (PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
480*557cf195SMatthew G. Knepley   PetscErrorCode (*exactFuncDers[1]) (PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx);
481*557cf195SMatthew G. Knepley   void            *exactCtxs[3];
482*557cf195SMatthew G. Knepley   MPI_Comm         comm;
483*557cf195SMatthew G. Knepley   PetscReal        error, errorDer, tol = PETSC_SMALL;
484*557cf195SMatthew G. Knepley   PetscErrorCode   ierr;
485*557cf195SMatthew G. Knepley 
486*557cf195SMatthew G. Knepley   PetscFunctionBeginUser;
487*557cf195SMatthew G. Knepley   exactCtxs[0]       = user;
488*557cf195SMatthew G. Knepley   exactCtxs[1]       = user;
489*557cf195SMatthew G. Knepley   exactCtxs[2]       = user;
490*557cf195SMatthew G. Knepley   user->constants[0] = 1.0;
491*557cf195SMatthew G. Knepley   user->constants[1] = 2.0;
492*557cf195SMatthew G. Knepley   user->constants[2] = 3.0;
493*557cf195SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
494*557cf195SMatthew G. Knepley   ierr = SetupFunctions(dm, PETSC_TRUE, order, 0, exactFuncs, exactFuncDers, user);CHKERRQ(ierr);
495*557cf195SMatthew G. Knepley   ierr = ComputeError(dm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user);CHKERRQ(ierr);
496*557cf195SMatthew G. Knepley   /* Report result */
497*557cf195SMatthew G. Knepley   if (error > tol)    {ierr = PetscPrintf(comm, "Function tests FAIL for order %D at tolerance %g error %g\n", order, (double)tol,(double) error);CHKERRQ(ierr);}
498*557cf195SMatthew G. Knepley   else                {ierr = PetscPrintf(comm, "Function tests pass for order %D at tolerance %g\n", order, (double)tol);CHKERRQ(ierr);}
499*557cf195SMatthew G. Knepley   if (errorDer > tol) {ierr = PetscPrintf(comm, "Function tests FAIL for order %D derivatives at tolerance %g error %g\n", order, (double)tol, (double)errorDer);CHKERRQ(ierr);}
500*557cf195SMatthew G. Knepley   else                {ierr = PetscPrintf(comm, "Function tests pass for order %D derivatives at tolerance %g\n", order, (double)tol);CHKERRQ(ierr);}
501*557cf195SMatthew G. Knepley   PetscFunctionReturn(0);
502*557cf195SMatthew G. Knepley }
503*557cf195SMatthew G. Knepley 
504*557cf195SMatthew G. Knepley /* Compare approximation to exact in L_2 */
505*557cf195SMatthew G. Knepley static PetscErrorCode CheckTransferError(DM fdm, PetscBool usePoly, PetscInt order, PetscInt dir, const char *testname, Vec fu, AppCtx *user)
506*557cf195SMatthew G. Knepley {
507*557cf195SMatthew G. Knepley   PetscErrorCode (*exactFuncs[1]) (PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
508*557cf195SMatthew G. Knepley   PetscErrorCode (*exactFuncDers[1]) (PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx);
509*557cf195SMatthew G. Knepley   PetscReal        n[3] = {1.0, 1.0, 1.0};
510*557cf195SMatthew G. Knepley   void            *exactCtxs[3];
511*557cf195SMatthew G. Knepley   MPI_Comm         comm;
512*557cf195SMatthew G. Knepley   PetscReal        error, errorDer, tol = PETSC_SMALL;
513*557cf195SMatthew G. Knepley   PetscErrorCode   ierr;
514*557cf195SMatthew G. Knepley 
515*557cf195SMatthew G. Knepley   PetscFunctionBeginUser;
516*557cf195SMatthew G. Knepley   exactCtxs[0]       = user;
517*557cf195SMatthew G. Knepley   exactCtxs[1]       = user;
518*557cf195SMatthew G. Knepley   exactCtxs[2]       = user;
519*557cf195SMatthew G. Knepley   user->constants[0] = 1.0;
520*557cf195SMatthew G. Knepley   user->constants[1] = 2.0;
521*557cf195SMatthew G. Knepley   user->constants[2] = 3.0;
522*557cf195SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) fdm, &comm);CHKERRQ(ierr);
523*557cf195SMatthew G. Knepley   ierr = SetupFunctions(fdm, usePoly, order, dir, exactFuncs, exactFuncDers, user);CHKERRQ(ierr);
524*557cf195SMatthew G. Knepley   ierr = DMComputeL2Diff(fdm, 0.0, exactFuncs, exactCtxs, fu, &error);CHKERRQ(ierr);
525*557cf195SMatthew G. Knepley   ierr = DMComputeL2GradientDiff(fdm, 0.0, exactFuncDers, exactCtxs, fu, n, &errorDer);CHKERRQ(ierr);
526*557cf195SMatthew G. Knepley   /* Report result */
527*557cf195SMatthew G. Knepley   if (error > tol)    {ierr = PetscPrintf(comm, "%s tests FAIL for order %D at tolerance %g error %g\n", testname, order, (double)tol, (double)error);CHKERRQ(ierr);}
528*557cf195SMatthew G. Knepley   else                {ierr = PetscPrintf(comm, "%s tests pass for order %D at tolerance %g\n", testname, order, (double)tol);CHKERRQ(ierr);}
529*557cf195SMatthew G. Knepley   if (errorDer > tol) {ierr = PetscPrintf(comm, "%s tests FAIL for order %D derivatives at tolerance %g error %g\n", testname, order, (double)tol, (double)errorDer);CHKERRQ(ierr);}
530*557cf195SMatthew G. Knepley   else                {ierr = PetscPrintf(comm, "%s tests pass for order %D derivatives at tolerance %g\n", testname, order, (double)tol);CHKERRQ(ierr);}
531*557cf195SMatthew G. Knepley   PetscFunctionReturn(0);
532*557cf195SMatthew G. Knepley }
533*557cf195SMatthew G. Knepley 
534*557cf195SMatthew G. Knepley static PetscErrorCode CheckTransfer(DM dm, InterpType inType, PetscInt order, AppCtx *user)
535*557cf195SMatthew G. Knepley {
536*557cf195SMatthew G. Knepley   PetscErrorCode (*exactFuncs[1]) (PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
537*557cf195SMatthew G. Knepley   PetscErrorCode (*exactFuncDers[1]) (PetscInt, PetscReal, const PetscReal x[], const PetscReal n[], PetscInt, PetscScalar *u, void *ctx);
538*557cf195SMatthew G. Knepley   void           *exactCtxs[3] = {user, user, user};
539*557cf195SMatthew G. Knepley   DM              rdm, idm, fdm;
540*557cf195SMatthew G. Knepley   Mat             Interp, InterpAdapt = NULL;
541*557cf195SMatthew G. Knepley   Vec             iu, fu, scaling;
542*557cf195SMatthew G. Knepley   MPI_Comm        comm;
543*557cf195SMatthew G. Knepley   const char     *testname;
544*557cf195SMatthew G. Knepley   char            checkname[PETSC_MAX_PATH_LEN];
545*557cf195SMatthew G. Knepley   PetscErrorCode  ierr;
546*557cf195SMatthew G. Knepley 
547*557cf195SMatthew G. Knepley   PetscFunctionBeginUser;
548*557cf195SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
549*557cf195SMatthew G. Knepley   ierr = DMRefine(dm, comm, &rdm);CHKERRQ(ierr);
550*557cf195SMatthew G. Knepley   ierr = DMSetCoarseDM(rdm, dm);CHKERRQ(ierr);
551*557cf195SMatthew G. Knepley   ierr = DMCopyDisc(dm, rdm);CHKERRQ(ierr);
552*557cf195SMatthew G. Knepley   switch (inType) {
553*557cf195SMatthew G. Knepley   case INTERPOLATION:
554*557cf195SMatthew G. Knepley     testname = "Interpolation";
555*557cf195SMatthew G. Knepley     idm = dm;
556*557cf195SMatthew G. Knepley     fdm = rdm;
557*557cf195SMatthew G. Knepley     break;
558*557cf195SMatthew G. Knepley   case RESTRICTION:
559*557cf195SMatthew G. Knepley     testname = "Restriction";
560*557cf195SMatthew G. Knepley     idm = rdm;
561*557cf195SMatthew G. Knepley     fdm = dm;
562*557cf195SMatthew G. Knepley     break;
563*557cf195SMatthew G. Knepley   case INJECTION:
564*557cf195SMatthew G. Knepley     testname = "Injection";
565*557cf195SMatthew G. Knepley     idm = rdm;
566*557cf195SMatthew G. Knepley     fdm = dm;
567*557cf195SMatthew G. Knepley     break;
568*557cf195SMatthew G. Knepley   }
569*557cf195SMatthew G. Knepley   ierr = DMGetGlobalVector(idm, &iu);CHKERRQ(ierr);
570*557cf195SMatthew G. Knepley   ierr = DMGetGlobalVector(fdm, &fu);CHKERRQ(ierr);
571*557cf195SMatthew G. Knepley   ierr = DMSetApplicationContext(dm, user);CHKERRQ(ierr);
572*557cf195SMatthew G. Knepley   ierr = DMSetApplicationContext(rdm, user);CHKERRQ(ierr);
573*557cf195SMatthew G. Knepley   /* Project function into initial FE function space */
574*557cf195SMatthew G. Knepley   ierr = SetupFunctions(dm, PETSC_TRUE, order, 0, exactFuncs, exactFuncDers, user);CHKERRQ(ierr);
575*557cf195SMatthew G. Knepley   ierr = DMProjectFunction(idm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, iu);CHKERRQ(ierr);
576*557cf195SMatthew G. Knepley   /* Interpolate function into final FE function space */
577*557cf195SMatthew G. Knepley   switch (inType) {
578*557cf195SMatthew G. Knepley   case INTERPOLATION:
579*557cf195SMatthew G. Knepley     ierr = DMCreateInterpolation(dm, rdm, &Interp, &scaling);CHKERRQ(ierr);
580*557cf195SMatthew G. Knepley     ierr = MatInterpolate(Interp, iu, fu);CHKERRQ(ierr);
581*557cf195SMatthew G. Knepley     break;
582*557cf195SMatthew G. Knepley   case RESTRICTION:
583*557cf195SMatthew G. Knepley     ierr = DMCreateInterpolation(dm, rdm, &Interp, &scaling);CHKERRQ(ierr);
584*557cf195SMatthew G. Knepley     ierr = MatRestrict(Interp, iu, fu);CHKERRQ(ierr);
585*557cf195SMatthew G. Knepley     ierr = VecPointwiseMult(fu, scaling, fu);CHKERRQ(ierr);
586*557cf195SMatthew G. Knepley     break;
587*557cf195SMatthew G. Knepley   case INJECTION:
588*557cf195SMatthew G. Knepley     ierr = DMCreateInjection(dm, rdm, &Interp);CHKERRQ(ierr);
589*557cf195SMatthew G. Knepley     ierr = MatRestrict(Interp, iu, fu);CHKERRQ(ierr);
590*557cf195SMatthew G. Knepley     break;
591*557cf195SMatthew G. Knepley   }
592*557cf195SMatthew G. Knepley   ierr = CheckTransferError(fdm, PETSC_TRUE, order, 0, testname, fu, user);CHKERRQ(ierr);
593*557cf195SMatthew G. Knepley   if (user->K && (inType == INTERPOLATION)) {
594*557cf195SMatthew G. Knepley     KSP      smoother;
595*557cf195SMatthew G. Knepley     Mat      A;
596*557cf195SMatthew G. Knepley     Vec     *iV, *fV;
597*557cf195SMatthew G. Knepley     PetscInt k, dim, d;
598*557cf195SMatthew G. Knepley 
599*557cf195SMatthew G. Knepley     ierr = PetscPrintf(comm, " Adapting interpolator using %s\n", user->usePoly ? "polynomials" : "harmonics");CHKERRQ(ierr);
600*557cf195SMatthew G. Knepley     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
601*557cf195SMatthew G. Knepley     ierr = PetscMalloc2(user->K*dim, &iV, user->K*dim, &fV);CHKERRQ(ierr);
602*557cf195SMatthew G. Knepley     /* Project coarse modes into initial and final FE function space */
603*557cf195SMatthew G. Knepley     for (k = 0; k < user->K; ++k) {
604*557cf195SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
605*557cf195SMatthew G. Knepley         ierr = DMGetGlobalVector(idm, &iV[k*dim+d]);CHKERRQ(ierr);
606*557cf195SMatthew G. Knepley         ierr = DMGetGlobalVector(fdm, &fV[k*dim+d]);CHKERRQ(ierr);
607*557cf195SMatthew G. Knepley         ierr = SetupFunctions(idm, user->usePoly, user->usePoly ? k : k+1, d, exactFuncs, exactFuncDers, user);CHKERRQ(ierr);
608*557cf195SMatthew G. Knepley         ierr = DMProjectFunction(idm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, iV[k*dim+d]);CHKERRQ(ierr);
609*557cf195SMatthew G. Knepley         ierr = DMProjectFunction(fdm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, fV[k*dim+d]);CHKERRQ(ierr);
610*557cf195SMatthew G. Knepley       }
611*557cf195SMatthew G. Knepley     }
612*557cf195SMatthew G. Knepley     /* Adapt interpolator */
613*557cf195SMatthew G. Knepley     ierr = DMCreateMatrix(rdm, &A);CHKERRQ(ierr);
614*557cf195SMatthew G. Knepley     ierr = MatShift(A, 1.0);CHKERRQ(ierr);
615*557cf195SMatthew G. Knepley     ierr = KSPCreate(comm, &smoother);CHKERRQ(ierr);
616*557cf195SMatthew G. Knepley     ierr = KSPSetFromOptions(smoother);CHKERRQ(ierr);
617*557cf195SMatthew G. Knepley     ierr = KSPSetOperators(smoother, A, A);CHKERRQ(ierr);
618*557cf195SMatthew G. Knepley     ierr = DMAdaptInterpolator(dm, rdm, Interp, smoother, user->K*dim, fV, iV, &InterpAdapt, user);CHKERRQ(ierr);
619*557cf195SMatthew G. Knepley     /* Interpolate function into final FE function space */
620*557cf195SMatthew G. Knepley     ierr = PetscSNPrintf(checkname, PETSC_MAX_PATH_LEN, "  %s poly", testname);CHKERRQ(ierr);
621*557cf195SMatthew G. Knepley     ierr = MatInterpolate(InterpAdapt, iu, fu);CHKERRQ(ierr);
622*557cf195SMatthew G. Knepley     ierr = CheckTransferError(fdm, PETSC_TRUE, order, 0, checkname, fu, user);CHKERRQ(ierr);
623*557cf195SMatthew G. Knepley     for (k = 0; k < user->K; ++k) {
624*557cf195SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
625*557cf195SMatthew G. Knepley         ierr = PetscSNPrintf(checkname, PETSC_MAX_PATH_LEN, "  %s trig (%D, %D)", testname, k, d);CHKERRQ(ierr);
626*557cf195SMatthew G. Knepley         ierr = MatInterpolate(InterpAdapt, iV[k*dim+d], fV[k*dim+d]);CHKERRQ(ierr);
627*557cf195SMatthew G. Knepley         ierr = CheckTransferError(fdm, PETSC_FALSE, k+1, d, checkname, fV[k*dim+d], user);CHKERRQ(ierr);
628*557cf195SMatthew G. Knepley       }
629*557cf195SMatthew G. Knepley     }
630*557cf195SMatthew G. Knepley     /* Cleanup */
631*557cf195SMatthew G. Knepley     ierr = KSPDestroy(&smoother);CHKERRQ(ierr);
632*557cf195SMatthew G. Knepley     ierr = MatDestroy(&A);CHKERRQ(ierr);
633*557cf195SMatthew G. Knepley     for (k = 0; k < user->K; ++k) {
634*557cf195SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
635*557cf195SMatthew G. Knepley         ierr = DMRestoreGlobalVector(idm, &iV[k*dim+d]);CHKERRQ(ierr);
636*557cf195SMatthew G. Knepley         ierr = DMRestoreGlobalVector(fdm, &fV[k*dim+d]);CHKERRQ(ierr);
637*557cf195SMatthew G. Knepley       }
638*557cf195SMatthew G. Knepley     }
639*557cf195SMatthew G. Knepley     ierr = PetscFree2(iV, fV);CHKERRQ(ierr);
640*557cf195SMatthew G. Knepley     ierr = MatDestroy(&InterpAdapt);CHKERRQ(ierr);
641*557cf195SMatthew G. Knepley   }
642*557cf195SMatthew G. Knepley   ierr = DMRestoreGlobalVector(idm, &iu);CHKERRQ(ierr);
643*557cf195SMatthew G. Knepley   ierr = DMRestoreGlobalVector(fdm, &fu);CHKERRQ(ierr);
644*557cf195SMatthew G. Knepley   ierr = MatDestroy(&Interp);CHKERRQ(ierr);
645*557cf195SMatthew G. Knepley   ierr = VecDestroy(&scaling);CHKERRQ(ierr);
646*557cf195SMatthew G. Knepley   ierr = DMDestroy(&rdm);CHKERRQ(ierr);
647*557cf195SMatthew G. Knepley   PetscFunctionReturn(0);
648*557cf195SMatthew G. Knepley }
649*557cf195SMatthew G. Knepley 
650*557cf195SMatthew G. Knepley int main(int argc, char **argv)
651*557cf195SMatthew G. Knepley {
652*557cf195SMatthew G. Knepley   DM             dm;
653*557cf195SMatthew G. Knepley   AppCtx         user;                 /* user-defined work context */
654*557cf195SMatthew G. Knepley   PetscErrorCode ierr;
655*557cf195SMatthew G. Knepley 
656*557cf195SMatthew G. Knepley   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
657*557cf195SMatthew G. Knepley   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
658*557cf195SMatthew G. Knepley   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
659*557cf195SMatthew G. Knepley   ierr = PetscFECreateDefault(PETSC_COMM_WORLD, user.dim, user.Nc, user.simplex, NULL, user.qorder, &user.fe);CHKERRQ(ierr);
660*557cf195SMatthew G. Knepley   ierr = DMSetField(dm, 0, NULL, (PetscObject) user.fe);CHKERRQ(ierr);
661*557cf195SMatthew G. Knepley   ierr = DMCreateDS(dm);CHKERRQ(ierr);
662*557cf195SMatthew G. Knepley   ierr = CheckFunctions(dm, user.porder, &user);CHKERRQ(ierr);
663*557cf195SMatthew G. Knepley   ierr = CheckTransfer(dm, INTERPOLATION, user.porder, &user);CHKERRQ(ierr);
664*557cf195SMatthew G. Knepley   ierr = CheckTransfer(dm, INJECTION,  user.porder, &user);CHKERRQ(ierr);
665*557cf195SMatthew G. Knepley   ierr = PetscFEDestroy(&user.fe);CHKERRQ(ierr);
666*557cf195SMatthew G. Knepley   ierr = DMDestroy(&dm);CHKERRQ(ierr);
667*557cf195SMatthew G. Knepley   ierr = PetscFinalize();
668*557cf195SMatthew G. Knepley   return ierr;
669*557cf195SMatthew G. Knepley }
670*557cf195SMatthew G. Knepley 
671*557cf195SMatthew G. Knepley /*TEST
672*557cf195SMatthew G. Knepley 
673*557cf195SMatthew G. Knepley   # TODO dim 3 will not work until I get composite elements in 3D (see plexrefine.c:34)
674*557cf195SMatthew G. Knepley   # 2D/3D P_1 on a simplex
675*557cf195SMatthew G. Knepley   test:
676*557cf195SMatthew G. Knepley     suffix: p1
677*557cf195SMatthew G. Knepley     requires: triangle ctetgen
678*557cf195SMatthew G. Knepley     args: -dim {{2}separate output} -petscspace_degree 1 -num_comp 1 -qorder 1 -porder {{1}separate output}
679*557cf195SMatthew G. Knepley   test:
680*557cf195SMatthew G. Knepley     suffix: p1_pragmatic
681*557cf195SMatthew G. Knepley     requires: triangle ctetgen pragmatic
682*557cf195SMatthew G. Knepley     args: -dim {{2}separate output} -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder {{1 2}separate output}
683*557cf195SMatthew G. Knepley   test:
684*557cf195SMatthew G. Knepley     suffix: p1_adapt
685*557cf195SMatthew G. Knepley     requires: triangle ctetgen
686*557cf195SMatthew G. Knepley     args: -dim {{2}separate output} -dm_refine 3 -petscspace_degree 1 -qorder 1 -porder {{1 2}separate output}
687*557cf195SMatthew G. Knepley 
688*557cf195SMatthew G. Knepley   # TODO dim 3 will not work until I get composite elements in 3D (see plexrefine.c:34)
689*557cf195SMatthew G. Knepley   # 2D/3D P_2 on a simplex
690*557cf195SMatthew G. Knepley   test:
691*557cf195SMatthew G. Knepley     suffix: p2
692*557cf195SMatthew G. Knepley     requires: triangle ctetgen
693*557cf195SMatthew G. Knepley     args: -dim {{2}separate output} -petscspace_degree 2 -qorder 2 -porder {{1 2 3}separate output}
694*557cf195SMatthew G. Knepley   test:
695*557cf195SMatthew G. Knepley     suffix: p2_pragmatic
696*557cf195SMatthew G. Knepley     requires: triangle ctetgen pragmatic
697*557cf195SMatthew G. Knepley     args: -dim {{2}separate output} -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder {{1 2 3}separate output}
698*557cf195SMatthew G. Knepley 
699*557cf195SMatthew G. Knepley   # TODO dim 3 will not work until I get composite elements in 3D (see plexrefine.c:34)
700*557cf195SMatthew G. Knepley   # TODO This is broken. Check ex3 which worked
701*557cf195SMatthew G. Knepley   # 2D/3D P_3 on a simplex
702*557cf195SMatthew G. Knepley   test:
703*557cf195SMatthew G. Knepley     suffix: p3
704*557cf195SMatthew G. Knepley     requires: triangle ctetgen !single
705*557cf195SMatthew G. Knepley     args: -dim {{2}separate output} -petscspace_degree 3 -qorder 3 -porder {{1 2 3 4}separate output}
706*557cf195SMatthew G. Knepley   test:
707*557cf195SMatthew G. Knepley     suffix: p3_pragmatic
708*557cf195SMatthew G. Knepley     requires: triangle ctetgen pragmatic !single
709*557cf195SMatthew G. Knepley     args: -dim {{2}separate output} -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -porder {{1 2 3 4}separate output}
710*557cf195SMatthew G. Knepley 
711*557cf195SMatthew G. Knepley   # 2D/3D Q_1 on a tensor cell
712*557cf195SMatthew G. Knepley   test:
713*557cf195SMatthew G. Knepley     suffix: q1
714*557cf195SMatthew G. Knepley     requires: mpi_type_get_envelope
715*557cf195SMatthew G. Knepley     args: -dim {{2 3}separate output} -simplex 0 -petscspace_degree 1 -qorder 1 -porder {{1 2}separate output}
716*557cf195SMatthew G. Knepley 
717*557cf195SMatthew G. Knepley   # 2D/3D Q_2 on a tensor cell
718*557cf195SMatthew G. Knepley   test:
719*557cf195SMatthew G. Knepley     suffix: q2
720*557cf195SMatthew G. Knepley     requires: mpi_type_get_envelope
721*557cf195SMatthew G. Knepley     args: -dim {{2 3}separate output} -simplex 0 -petscspace_degree 2 -qorder 2 -porder {{1 2 3}separate output}
722*557cf195SMatthew G. Knepley 
723*557cf195SMatthew G. Knepley   # 2D/3D Q_3 on a tensor cell
724*557cf195SMatthew G. Knepley   test:
725*557cf195SMatthew G. Knepley     suffix: q3
726*557cf195SMatthew G. Knepley     requires: mpi_type_get_envelope !single
727*557cf195SMatthew G. Knepley     args: -dim {{2 3}separate output} -simplex 0 -petscspace_degree 3 -qorder 3 -porder {{1 2 3 4}separate output}
728*557cf195SMatthew G. Knepley 
729*557cf195SMatthew G. Knepley   # 2D/3D P_1disc on a triangle/quadrilateral
730*557cf195SMatthew G. Knepley   # TODO Missing injection functional for simplices
731*557cf195SMatthew G. Knepley   test:
732*557cf195SMatthew G. Knepley     suffix: p1d
733*557cf195SMatthew G. Knepley     requires: triangle ctetgen
734*557cf195SMatthew G. Knepley     args: -dim {{2}separate output} -simplex {{0}separate output} -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder {{1 2}separate output}
735*557cf195SMatthew G. Knepley 
736*557cf195SMatthew G. Knepley TEST*/
737