xref: /petsc/src/snes/tests/ex8.c (revision d6837840dcd6a76b85f9d85f7ce6c5e3b767934f)
1557cf195SMatthew G. Knepley static char help[] = "Test adaptive interpolation of functions of a given polynomial order\n\n";
2557cf195SMatthew G. Knepley 
3557cf195SMatthew G. Knepley #include <petscdmplex.h>
4557cf195SMatthew G. Knepley #include <petscsnes.h>
5557cf195SMatthew G. Knepley 
6557cf195SMatthew G. Knepley /*
7557cf195SMatthew G. Knepley   What properties does the adapted interpolator have?
8557cf195SMatthew G. Knepley 
9557cf195SMatthew G. Knepley 1) If we adapt to quadratics, we can get lower interpolation error for quadratics (than local interpolation) when using a linear basis
10557cf195SMatthew G. Knepley 
11557cf195SMatthew G. Knepley $ ./ex8 -dm_refine 2 -petscspace_degree 1 -qorder 1 -dim 2 -porder 2 -K 2 -num_comp 1 -use_poly 1
12557cf195SMatthew G. Knepley Function tests FAIL for order 2 at tolerance 1e-10 error 0.00273757
13557cf195SMatthew G. Knepley Function tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.0721688
14557cf195SMatthew G. Knepley Interpolation tests FAIL for order 2 at tolerance 1e-10 error 0.00284555
15557cf195SMatthew G. Knepley Interpolation tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.0721688
16557cf195SMatthew G. Knepley  Adapting interpolator using polynomials
17557cf195SMatthew G. Knepley The number of input vectors 4 < 7 the maximum number of column entries
18557cf195SMatthew G. Knepley   Interpolation poly tests FAIL for order 2 at tolerance 1e-10 error 0.00659864
19557cf195SMatthew G. Knepley   Interpolation poly tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.0836582
20557cf195SMatthew G. Knepley   Interpolation trig (0, 0) tests FAIL for order 1 at tolerance 1e-10 error 0.476194
21557cf195SMatthew G. Knepley   Interpolation trig (0, 0) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22144
22557cf195SMatthew G. Knepley   Interpolation trig (0, 1) tests FAIL for order 1 at tolerance 1e-10 error 1.39768
23557cf195SMatthew G. Knepley   Interpolation trig (0, 1) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22144
24557cf195SMatthew G. Knepley   Interpolation trig (1, 0) tests FAIL for order 2 at tolerance 1e-10 error 1.07315
25557cf195SMatthew G. Knepley   Interpolation trig (1, 0) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55403
26557cf195SMatthew G. Knepley   Interpolation trig (1, 1) tests FAIL for order 2 at tolerance 1e-10 error 1.07315
27557cf195SMatthew G. Knepley   Interpolation trig (1, 1) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55403
28557cf195SMatthew G. Knepley 
29557cf195SMatthew G. Knepley $ ./ex8 -dm_refine 2 -petscspace_degree 1 -qorder 1 -dim 2 -porder 2 -K 3 -num_comp 1 -use_poly 1
30557cf195SMatthew G. Knepley Function tests FAIL for order 2 at tolerance 1e-10 error 0.00273757
31557cf195SMatthew G. Knepley Function tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.0721688
32557cf195SMatthew G. Knepley Interpolation tests FAIL for order 2 at tolerance 1e-10 error 0.00284555
33557cf195SMatthew G. Knepley Interpolation tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.0721688
34557cf195SMatthew G. Knepley  Adapting interpolator using polynomials
35557cf195SMatthew G. Knepley The number of input vectors 6 < 7 the maximum number of column entries
36557cf195SMatthew G. Knepley   Interpolation poly tests FAIL for order 2 at tolerance 1e-10 error 0.00194055
37557cf195SMatthew G. Knepley   Interpolation poly tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.0525591
38557cf195SMatthew G. Knepley   Interpolation trig (0, 0) tests FAIL for order 1 at tolerance 1e-10 error 0.476255
39557cf195SMatthew G. Knepley   Interpolation trig (0, 0) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22132
40557cf195SMatthew G. Knepley   Interpolation trig (0, 1) tests FAIL for order 1 at tolerance 1e-10 error 1.39785
41557cf195SMatthew G. Knepley   Interpolation trig (0, 1) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22119
42557cf195SMatthew G. Knepley   Interpolation trig (1, 0) tests FAIL for order 2 at tolerance 1e-10 error 1.0727
43557cf195SMatthew G. Knepley   Interpolation trig (1, 0) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55364
44557cf195SMatthew G. Knepley   Interpolation trig (1, 1) tests FAIL for order 2 at tolerance 1e-10 error 1.0727
45557cf195SMatthew G. Knepley   Interpolation trig (1, 1) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55364
46557cf195SMatthew G. Knepley   Interpolation trig (2, 0) tests FAIL for order 3 at tolerance 1e-10 error 0.705258
47557cf195SMatthew G. Knepley   Interpolation trig (2, 0) tests FAIL for order 3 derivatives at tolerance 1e-10 error 6.82037
48557cf195SMatthew G. Knepley   Interpolation trig (2, 1) tests FAIL for order 3 at tolerance 1e-10 error 0.705258
49557cf195SMatthew G. Knepley   Interpolation trig (2, 1) tests FAIL for order 3 derivatives at tolerance 1e-10 error 6.82037
50557cf195SMatthew G. Knepley 
51557cf195SMatthew G. Knepley 2) We can more accurately capture low harmonics
52557cf195SMatthew G. Knepley 
53557cf195SMatthew G. Knepley If we adapt polynomials, we can be exact
54557cf195SMatthew G. Knepley 
55557cf195SMatthew G. Knepley $ ./ex8 -dm_refine 2 -petscspace_degree 1 -qorder 1 -dim 2 -porder 1 -K 2 -num_comp 1 -use_poly 1
56557cf195SMatthew G. Knepley Function tests pass for order 1 at tolerance 1e-10
57557cf195SMatthew G. Knepley Function tests pass for order 1 derivatives at tolerance 1e-10
58557cf195SMatthew G. Knepley Interpolation tests pass for order 1 at tolerance 1e-10
59557cf195SMatthew G. Knepley Interpolation tests pass for order 1 derivatives at tolerance 1e-10
60557cf195SMatthew G. Knepley  Adapting interpolator using polynomials
61557cf195SMatthew G. Knepley The number of input vectors 4 < 7 the maximum number of column entries
62557cf195SMatthew G. Knepley   Interpolation poly tests pass for order 1 at tolerance 1e-10
63557cf195SMatthew G. Knepley   Interpolation poly tests pass for order 1 derivatives at tolerance 1e-10
64557cf195SMatthew G. Knepley   Interpolation trig (0, 0) tests FAIL for order 1 at tolerance 1e-10 error 0.476194
65557cf195SMatthew G. Knepley   Interpolation trig (0, 0) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22144
66557cf195SMatthew G. Knepley   Interpolation trig (0, 1) tests FAIL for order 1 at tolerance 1e-10 error 1.39768
67557cf195SMatthew G. Knepley   Interpolation trig (0, 1) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22144
68557cf195SMatthew G. Knepley   Interpolation trig (1, 0) tests FAIL for order 2 at tolerance 1e-10 error 1.07315
69557cf195SMatthew G. Knepley   Interpolation trig (1, 0) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55403
70557cf195SMatthew G. Knepley   Interpolation trig (1, 1) tests FAIL for order 2 at tolerance 1e-10 error 1.07315
71557cf195SMatthew G. Knepley   Interpolation trig (1, 1) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55403
72557cf195SMatthew G. Knepley 
73557cf195SMatthew G. Knepley and least for small K,
74557cf195SMatthew G. Knepley 
75557cf195SMatthew G. Knepley $ ./ex8 -dm_refine 2 -petscspace_degree 1 -qorder 1 -dim 2 -porder 1 -K 4 -num_comp 1 -use_poly 1
76557cf195SMatthew G. Knepley Function tests pass for order 1 at tolerance 1e-10
77557cf195SMatthew G. Knepley Function tests pass for order 1 derivatives at tolerance 1e-10
78557cf195SMatthew G. Knepley Interpolation tests pass for order 1 at tolerance 1e-10
79557cf195SMatthew G. Knepley Interpolation tests pass for order 1 derivatives at tolerance 1e-10
80557cf195SMatthew G. Knepley  Adapting interpolator using polynomials
81557cf195SMatthew G. Knepley   Interpolation poly tests FAIL for order 1 at tolerance 1e-10 error 0.0015351
82557cf195SMatthew G. Knepley   Interpolation poly tests FAIL for order 1 derivatives at tolerance 1e-10 error 0.0427369
83557cf195SMatthew G. Knepley   Interpolation trig (0, 0) tests FAIL for order 1 at tolerance 1e-10 error 0.476359
84557cf195SMatthew G. Knepley   Interpolation trig (0, 0) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22115
85557cf195SMatthew G. Knepley   Interpolation trig (0, 1) tests FAIL for order 1 at tolerance 1e-10 error 1.3981
86557cf195SMatthew G. Knepley   Interpolation trig (0, 1) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22087
87557cf195SMatthew G. Knepley   Interpolation trig (1, 0) tests FAIL for order 2 at tolerance 1e-10 error 1.07228
88557cf195SMatthew G. Knepley   Interpolation trig (1, 0) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55238
89557cf195SMatthew G. Knepley   Interpolation trig (1, 1) tests FAIL for order 2 at tolerance 1e-10 error 1.07228
90557cf195SMatthew G. Knepley   Interpolation trig (1, 1) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55238
91557cf195SMatthew G. Knepley   Interpolation trig (2, 0) tests FAIL for order 3 at tolerance 1e-10 error 0.704947
92557cf195SMatthew G. Knepley   Interpolation trig (2, 0) tests FAIL for order 3 derivatives at tolerance 1e-10 error 6.82254
93557cf195SMatthew G. Knepley   Interpolation trig (2, 1) tests FAIL for order 3 at tolerance 1e-10 error 0.704948
94557cf195SMatthew G. Knepley   Interpolation trig (2, 1) tests FAIL for order 3 derivatives at tolerance 1e-10 error 6.82254
95557cf195SMatthew G. Knepley   Interpolation trig (3, 0) tests FAIL for order 4 at tolerance 1e-10 error 0.893279
96557cf195SMatthew G. Knepley   Interpolation trig (3, 0) tests FAIL for order 4 derivatives at tolerance 1e-10 error 8.93718
97557cf195SMatthew G. Knepley   Interpolation trig (3, 1) tests FAIL for order 4 at tolerance 1e-10 error 0.89328
98557cf195SMatthew G. Knepley   Interpolation trig (3, 1) tests FAIL for order 4 derivatives at tolerance 1e-10 error 8.93717
99557cf195SMatthew G. Knepley 
100557cf195SMatthew G. Knepley but adapting to harmonics gives alright polynomials errors and much better harmonics errors.
101557cf195SMatthew G. Knepley 
102557cf195SMatthew G. Knepley $ ./ex8 -dm_refine 2 -petscspace_degree 1 -qorder 1 -dim 2 -porder 1 -K 4 -num_comp 1 -use_poly 0
103557cf195SMatthew G. Knepley Function tests pass for order 1 at tolerance 1e-10
104557cf195SMatthew G. Knepley Function tests pass for order 1 derivatives at tolerance 1e-10
105557cf195SMatthew G. Knepley Interpolation tests pass for order 1 at tolerance 1e-10
106557cf195SMatthew G. Knepley Interpolation tests pass for order 1 derivatives at tolerance 1e-10
107557cf195SMatthew G. Knepley  Adapting interpolator using harmonics
108557cf195SMatthew G. Knepley   Interpolation poly tests FAIL for order 1 at tolerance 1e-10 error 0.0720606
109557cf195SMatthew G. Knepley   Interpolation poly tests FAIL for order 1 derivatives at tolerance 1e-10 error 1.97779
110557cf195SMatthew G. Knepley   Interpolation trig (0, 0) tests FAIL for order 1 at tolerance 1e-10 error 0.0398055
111557cf195SMatthew G. Knepley   Interpolation trig (0, 0) tests FAIL for order 1 derivatives at tolerance 1e-10 error 0.995963
112557cf195SMatthew G. Knepley   Interpolation trig (0, 1) tests FAIL for order 1 at tolerance 1e-10 error 0.0398051
113557cf195SMatthew G. Knepley   Interpolation trig (0, 1) tests FAIL for order 1 derivatives at tolerance 1e-10 error 0.995964
114557cf195SMatthew G. Knepley   Interpolation trig (1, 0) tests FAIL for order 2 at tolerance 1e-10 error 0.0238441
115557cf195SMatthew G. Knepley   Interpolation trig (1, 0) tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.888611
116557cf195SMatthew G. Knepley   Interpolation trig (1, 1) tests FAIL for order 2 at tolerance 1e-10 error 0.0238346
117557cf195SMatthew G. Knepley   Interpolation trig (1, 1) tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.888612
118557cf195SMatthew G. Knepley   Interpolation trig (2, 0) tests FAIL for order 3 at tolerance 1e-10 error 0.0537968
119557cf195SMatthew G. Knepley   Interpolation trig (2, 0) tests FAIL for order 3 derivatives at tolerance 1e-10 error 1.57665
120557cf195SMatthew G. Knepley   Interpolation trig (2, 1) tests FAIL for order 3 at tolerance 1e-10 error 0.0537779
121557cf195SMatthew G. Knepley   Interpolation trig (2, 1) tests FAIL for order 3 derivatives at tolerance 1e-10 error 1.57666
122557cf195SMatthew G. Knepley   Interpolation trig (3, 0) tests FAIL for order 4 at tolerance 1e-10 error 0.0775838
123557cf195SMatthew G. Knepley   Interpolation trig (3, 0) tests FAIL for order 4 derivatives at tolerance 1e-10 error 2.36926
124557cf195SMatthew G. Knepley   Interpolation trig (3, 1) tests FAIL for order 4 at tolerance 1e-10 error 0.0775464
125557cf195SMatthew G. Knepley   Interpolation trig (3, 1) tests FAIL for order 4 derivatives at tolerance 1e-10 error 2.36929
126557cf195SMatthew G. Knepley */
127557cf195SMatthew G. Knepley 
128557cf195SMatthew G. Knepley typedef struct {
129557cf195SMatthew G. Knepley   /* Domain and mesh definition */
130557cf195SMatthew G. Knepley   PetscInt  dim;               /* The topological mesh dimension */
131557cf195SMatthew G. Knepley   PetscBool simplex;           /* Flag for simplex or tensor product mesh */
132557cf195SMatthew G. Knepley   PetscBool interpolate;       /* Generate intermediate mesh elements */
133557cf195SMatthew G. Knepley   PetscReal refinementLimit;   /* The largest allowable cell volume */
134557cf195SMatthew G. Knepley   /* Element definition */
135557cf195SMatthew G. Knepley   PetscInt  qorder;            /* Order of the quadrature */
136557cf195SMatthew G. Knepley   PetscInt  Nc;                /* Number of field components */
137557cf195SMatthew G. Knepley   PetscFE   fe;                /* The finite element */
138557cf195SMatthew G. Knepley   /* Testing space */
139557cf195SMatthew G. Knepley   PetscInt  porder;            /* Order of polynomials to test */
140557cf195SMatthew G. Knepley   PetscReal constants[3];      /* Constant values for each dimension */
141557cf195SMatthew G. Knepley   PetscInt  m;                 /* The frequency of sinusoids to use */
142557cf195SMatthew G. Knepley   PetscInt  dir;               /* The direction of sinusoids to use */
143557cf195SMatthew G. Knepley   /* Adaptation */
144557cf195SMatthew G. Knepley   PetscInt  K;                 /* Number of coarse modes used for optimization */
145557cf195SMatthew G. Knepley   PetscBool usePoly;           /* Use polynomials, or harmonics, to adapt interpolator */
146557cf195SMatthew G. Knepley } AppCtx;
147557cf195SMatthew G. Knepley 
148557cf195SMatthew G. Knepley typedef enum {INTERPOLATION, RESTRICTION, INJECTION} InterpType;
149557cf195SMatthew G. Knepley 
150557cf195SMatthew G. Knepley /* u = 1 */
151557cf195SMatthew G. Knepley PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
152557cf195SMatthew G. Knepley {
153557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
154557cf195SMatthew G. Knepley   PetscInt d = user->dir;
155557cf195SMatthew G. Knepley 
156557cf195SMatthew G. Knepley   if (Nc > 1) {
157557cf195SMatthew G. Knepley     for (d = 0; d < Nc; ++d) u[d] = user->constants[d];
158557cf195SMatthew G. Knepley   } else {
159557cf195SMatthew G. Knepley     u[0] = user->constants[d];
160557cf195SMatthew G. Knepley   }
161557cf195SMatthew G. Knepley   return 0;
162557cf195SMatthew G. Knepley }
163557cf195SMatthew G. Knepley PetscErrorCode constantDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nc, PetscScalar *u, void *ctx)
164557cf195SMatthew G. Knepley {
165557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
166557cf195SMatthew G. Knepley   PetscInt d = user->dir;
167557cf195SMatthew G. Knepley 
168557cf195SMatthew G. Knepley   if (Nc > 1) {
169557cf195SMatthew G. Knepley     for (d = 0; d < Nc; ++d) u[d] = 0.0;
170557cf195SMatthew G. Knepley   } else {
171557cf195SMatthew G. Knepley     u[0] = user->constants[d];
172557cf195SMatthew G. Knepley   }
173557cf195SMatthew G. Knepley   return 0;
174557cf195SMatthew G. Knepley }
175557cf195SMatthew G. Knepley 
176557cf195SMatthew G. Knepley /* u = x */
177557cf195SMatthew G. Knepley PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
178557cf195SMatthew G. Knepley {
179557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
180557cf195SMatthew G. Knepley   PetscInt d = user->dir;
181557cf195SMatthew G. Knepley 
182557cf195SMatthew G. Knepley   if (Nc > 1) {
183557cf195SMatthew G. Knepley     for (d = 0; d < Nc; ++d) u[d] = coords[d];
184557cf195SMatthew G. Knepley   } else {
185557cf195SMatthew G. Knepley     u[0] = coords[d];
186557cf195SMatthew G. Knepley   }
187557cf195SMatthew G. Knepley   return 0;
188557cf195SMatthew G. Knepley }
189557cf195SMatthew G. Knepley PetscErrorCode linearDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nc, PetscScalar *u, void *ctx)
190557cf195SMatthew G. Knepley {
191557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
192557cf195SMatthew G. Knepley   PetscInt d = user->dir;
193557cf195SMatthew G. Knepley 
194557cf195SMatthew G. Knepley   if (Nc > 1) {
195557cf195SMatthew G. Knepley     PetscInt e;
196557cf195SMatthew G. Knepley     for (d = 0; d < Nc; ++d) {
197557cf195SMatthew G. Knepley       u[d] = 0.0;
198557cf195SMatthew G. Knepley       for (e = 0; e < dim; ++e) u[d] += (d == e ? 1.0 : 0.0) * n[e];
199557cf195SMatthew G. Knepley     }
200557cf195SMatthew G. Knepley   } else {
201557cf195SMatthew G. Knepley     u[0] = n[d];
202557cf195SMatthew G. Knepley   }
203557cf195SMatthew G. Knepley   return 0;
204557cf195SMatthew G. Knepley }
205557cf195SMatthew G. Knepley 
206557cf195SMatthew G. Knepley /* u = x^2 or u = (x^2, xy) or u = (xy, yz, zx) */
207557cf195SMatthew G. Knepley PetscErrorCode quadratic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
208557cf195SMatthew G. Knepley {
209557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
210557cf195SMatthew G. Knepley   PetscInt d = user->dir;
211557cf195SMatthew G. Knepley 
212557cf195SMatthew G. Knepley   if (Nc > 1) {
213557cf195SMatthew G. Knepley     if (Nc > 2) {u[0] = coords[0]*coords[1]; u[1] = coords[1]*coords[2]; u[2] = coords[2]*coords[0];}
214557cf195SMatthew G. Knepley     else        {u[0] = coords[0]*coords[0]; u[1] = coords[0]*coords[1];}
215557cf195SMatthew G. Knepley   } else {
216557cf195SMatthew G. Knepley     u[0] = coords[d]*coords[d];
217557cf195SMatthew G. Knepley   }
218557cf195SMatthew G. Knepley   return 0;
219557cf195SMatthew G. Knepley }
220557cf195SMatthew G. Knepley PetscErrorCode quadraticDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nc, PetscScalar *u, void *ctx)
221557cf195SMatthew G. Knepley {
222557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
223557cf195SMatthew G. Knepley   PetscInt d = user->dir;
224557cf195SMatthew G. Knepley 
225557cf195SMatthew G. Knepley   if (Nc > 1) {
226557cf195SMatthew 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];}
227557cf195SMatthew G. Knepley     else        {u[0] = 2.0*coords[0]*n[0]; u[1] = coords[1]*n[0] + coords[0]*n[1];}
228557cf195SMatthew G. Knepley   } else {
229557cf195SMatthew G. Knepley     u[0] = 2.0*coords[d]*n[d];
230557cf195SMatthew G. Knepley   }
231557cf195SMatthew G. Knepley   return 0;
232557cf195SMatthew G. Knepley }
233557cf195SMatthew G. Knepley 
234557cf195SMatthew G. Knepley /* u = x^3 or u = (x^3, x^2y) or u = (x^2y, y^2z, z^2x) */
235557cf195SMatthew G. Knepley PetscErrorCode cubic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
236557cf195SMatthew G. Knepley {
237557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
238557cf195SMatthew G. Knepley   PetscInt d = user->dir;
239557cf195SMatthew G. Knepley 
240557cf195SMatthew G. Knepley   if (Nc > 1) {
241557cf195SMatthew 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];}
242557cf195SMatthew G. Knepley     else        {u[0] = coords[0]*coords[0]*coords[0]; u[1] = coords[0]*coords[0]*coords[1];}
243557cf195SMatthew G. Knepley   } else {
244557cf195SMatthew G. Knepley     u[0] = coords[d]*coords[d]*coords[d];
245557cf195SMatthew G. Knepley   }
246557cf195SMatthew G. Knepley   return 0;
247557cf195SMatthew G. Knepley }
248557cf195SMatthew G. Knepley PetscErrorCode cubicDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nc, PetscScalar *u, void *ctx)
249557cf195SMatthew G. Knepley {
250557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
251557cf195SMatthew G. Knepley   PetscInt d = user->dir;
252557cf195SMatthew G. Knepley 
253557cf195SMatthew G. Knepley   if (Nc > 1) {
254557cf195SMatthew 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];}
255557cf195SMatthew 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];}
256557cf195SMatthew G. Knepley   } else {
257557cf195SMatthew G. Knepley     u[0] = 3.0*coords[d]*coords[d]*n[d];
258557cf195SMatthew G. Knepley   }
259557cf195SMatthew G. Knepley   return 0;
260557cf195SMatthew G. Knepley }
261557cf195SMatthew G. Knepley 
262557cf195SMatthew G. Knepley /* u = x^4 or u = (x^4, x^2y^2) or u = (x^2y^2, y^2z^2, z^2x^2) */
263557cf195SMatthew G. Knepley PetscErrorCode quartic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
264557cf195SMatthew G. Knepley {
265557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
266557cf195SMatthew G. Knepley   PetscInt d = user->dir;
267557cf195SMatthew G. Knepley 
268557cf195SMatthew G. Knepley   if (Nc > 1) {
269557cf195SMatthew 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];}
270557cf195SMatthew G. Knepley     else        {u[0] = coords[0]*coords[0]*coords[0]*coords[0]; u[1] = coords[0]*coords[0]*coords[1]*coords[1];}
271557cf195SMatthew G. Knepley   } else {
272557cf195SMatthew G. Knepley     u[0] = coords[d]*coords[d]*coords[d]*coords[d];
273557cf195SMatthew G. Knepley   }
274557cf195SMatthew G. Knepley   return 0;
275557cf195SMatthew G. Knepley }
276557cf195SMatthew G. Knepley PetscErrorCode quarticDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nc, PetscScalar *u, void *ctx)
277557cf195SMatthew G. Knepley {
278557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
279557cf195SMatthew G. Knepley   PetscInt d = user->dir;
280557cf195SMatthew G. Knepley 
281557cf195SMatthew G. Knepley   if (Nc > 1) {
282557cf195SMatthew 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];
283557cf195SMatthew G. Knepley                  u[1] = 2.0*coords[1]*coords[2]*coords[2]*n[1] + 2.0*coords[1]*coords[1]*coords[2]*n[2];
284557cf195SMatthew G. Knepley                  u[2] = 2.0*coords[2]*coords[0]*coords[0]*n[2] + 2.0*coords[2]*coords[2]*coords[0]*n[0];}
285557cf195SMatthew 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];}
286557cf195SMatthew G. Knepley   } else {
287557cf195SMatthew G. Knepley     u[0] = 4.0*coords[d]*coords[d]*coords[d]*n[d];
288557cf195SMatthew G. Knepley   }
289557cf195SMatthew G. Knepley   return 0;
290557cf195SMatthew G. Knepley }
291557cf195SMatthew G. Knepley 
292557cf195SMatthew G. Knepley PetscErrorCode mytanh(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
293557cf195SMatthew G. Knepley {
294557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
295557cf195SMatthew G. Knepley   PetscInt d = user->dir;
296557cf195SMatthew G. Knepley 
297557cf195SMatthew G. Knepley   if (Nc > 1) {
298557cf195SMatthew G. Knepley     for (d = 0; d < Nc; ++d) u[d] = PetscTanhReal(coords[d] - 0.5);
299557cf195SMatthew G. Knepley   } else {
300557cf195SMatthew G. Knepley     u[0] = PetscTanhReal(coords[d] - 0.5);
301557cf195SMatthew G. Knepley   }
302557cf195SMatthew G. Knepley   return 0;
303557cf195SMatthew G. Knepley }
304557cf195SMatthew G. Knepley PetscErrorCode mytanhDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nc, PetscScalar *u, void *ctx)
305557cf195SMatthew G. Knepley {
306557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
307557cf195SMatthew G. Knepley   PetscInt d = user->dir;
308557cf195SMatthew G. Knepley 
309557cf195SMatthew G. Knepley   if (Nc > 1) {
310557cf195SMatthew G. Knepley     for (d = 0; d < Nc; ++d) u[d] = 1.0/PetscSqr(PetscCoshReal(coords[d] - 0.5)) * n[d];
311557cf195SMatthew G. Knepley   } else {
312557cf195SMatthew G. Knepley     u[0] = 1.0/PetscSqr(PetscCoshReal(coords[d] - 0.5)) * n[d];
313557cf195SMatthew G. Knepley   }
314557cf195SMatthew G. Knepley   return 0;
315557cf195SMatthew G. Knepley }
316557cf195SMatthew G. Knepley 
317557cf195SMatthew G. Knepley PetscErrorCode trig(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
318557cf195SMatthew G. Knepley {
319557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
320557cf195SMatthew G. Knepley   PetscInt m = user->m, d = user->dir;
321557cf195SMatthew G. Knepley 
322557cf195SMatthew G. Knepley   if (Nc > 1) {
323557cf195SMatthew G. Knepley     for (d = 0; d < Nc; ++d) u[d] = PetscSinReal(PETSC_PI*m*coords[d]);
324557cf195SMatthew G. Knepley   } else {
325557cf195SMatthew G. Knepley     u[0] = PetscSinReal(PETSC_PI*m*coords[d]);
326557cf195SMatthew G. Knepley   }
327557cf195SMatthew G. Knepley   return 0;
328557cf195SMatthew G. Knepley }
329557cf195SMatthew G. Knepley PetscErrorCode trigDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nc, PetscScalar *u, void *ctx)
330557cf195SMatthew G. Knepley {
331557cf195SMatthew G. Knepley   AppCtx  *user = (AppCtx *) ctx;
332557cf195SMatthew G. Knepley   PetscInt m = user->m, d = user->dir;
333557cf195SMatthew G. Knepley 
334557cf195SMatthew G. Knepley   if (Nc > 1) {
335557cf195SMatthew G. Knepley     for (d = 0; d < Nc; ++d) u[d] = PETSC_PI*m*PetscCosReal(PETSC_PI*m*coords[d]) * n[d];
336557cf195SMatthew G. Knepley   } else {
337557cf195SMatthew G. Knepley     u[0] = PETSC_PI*m*PetscCosReal(PETSC_PI*m*coords[d]) * n[d];
338557cf195SMatthew G. Knepley   }
339557cf195SMatthew G. Knepley   return 0;
340557cf195SMatthew G. Knepley }
341557cf195SMatthew G. Knepley 
342557cf195SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
343557cf195SMatthew G. Knepley {
344557cf195SMatthew G. Knepley   PetscErrorCode ierr;
345557cf195SMatthew G. Knepley 
346557cf195SMatthew G. Knepley   PetscFunctionBeginUser;
347557cf195SMatthew G. Knepley   options->dim             = 2;
348557cf195SMatthew G. Knepley   options->simplex         = PETSC_TRUE;
349557cf195SMatthew G. Knepley   options->interpolate     = PETSC_TRUE;
350557cf195SMatthew G. Knepley   options->refinementLimit = 0.0;
351557cf195SMatthew G. Knepley   options->qorder          = 0;
352557cf195SMatthew G. Knepley   options->Nc              = PETSC_DEFAULT;
353557cf195SMatthew G. Knepley   options->porder          = 0;
354557cf195SMatthew G. Knepley   options->m               = 1;
355557cf195SMatthew G. Knepley   options->dir             = 0;
356557cf195SMatthew G. Knepley   options->K               = 0;
357557cf195SMatthew G. Knepley   options->usePoly         = PETSC_TRUE;
358557cf195SMatthew G. Knepley 
359557cf195SMatthew G. Knepley   ierr = PetscOptionsBegin(comm, "", "Projection Test Options", "DMPlex");CHKERRQ(ierr);
360557cf195SMatthew G. Knepley   ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex8.c", options->dim, &options->dim, NULL);CHKERRQ(ierr);
361*d6837840SMatthew G. Knepley   ierr = PetscOptionsBool("-simplex", "Flag for simplices or hexahedra", "ex8.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr);
362557cf195SMatthew G. Knepley   ierr = PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex8.c", options->interpolate, &options->interpolate, NULL);CHKERRQ(ierr);
363557cf195SMatthew G. Knepley   ierr = PetscOptionsReal("-refinement_limit", "The largest allowable cell volume", "ex8.c", options->refinementLimit, &options->refinementLimit, NULL);CHKERRQ(ierr);
364557cf195SMatthew G. Knepley   ierr = PetscOptionsInt("-qorder", "The quadrature order", "ex8.c", options->qorder, &options->qorder, NULL);CHKERRQ(ierr);
365557cf195SMatthew G. Knepley   ierr = PetscOptionsInt("-num_comp", "The number of field components", "ex8.c", options->Nc, &options->Nc, NULL);CHKERRQ(ierr);
366557cf195SMatthew G. Knepley   ierr = PetscOptionsInt("-porder", "The order of polynomials to test", "ex8.c", options->porder, &options->porder, NULL);CHKERRQ(ierr);
367557cf195SMatthew G. Knepley   ierr = PetscOptionsInt("-K", "The number of coarse modes used in optimization", "ex8.c", options->K, &options->K, NULL);CHKERRQ(ierr);
368557cf195SMatthew G. Knepley   ierr = PetscOptionsBool("-use_poly", "Use polynomials (or harmonics) to adapt interpolator", "ex8.c", options->usePoly, &options->usePoly, NULL);CHKERRQ(ierr);
369557cf195SMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
370557cf195SMatthew G. Knepley 
371557cf195SMatthew G. Knepley   options->Nc = options->Nc < 0 ? options->dim : options->Nc;
372557cf195SMatthew G. Knepley 
373557cf195SMatthew G. Knepley   PetscFunctionReturn(0);
374557cf195SMatthew G. Knepley }
375557cf195SMatthew G. Knepley 
376557cf195SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
377557cf195SMatthew G. Knepley {
378557cf195SMatthew G. Knepley   DM               pdm      = NULL;
379557cf195SMatthew G. Knepley   PetscInt         cells[3] = {2, 2, 2};
380557cf195SMatthew G. Knepley   PetscPartitioner part;
381557cf195SMatthew G. Knepley   PetscErrorCode   ierr;
382557cf195SMatthew G. Knepley 
383557cf195SMatthew G. Knepley   PetscFunctionBeginUser;
384557cf195SMatthew G. Knepley   ierr = PetscOptionsGetInt(NULL, NULL, "-da_grid_x", &cells[0], NULL);CHKERRQ(ierr);
385557cf195SMatthew G. Knepley   ierr = PetscOptionsGetInt(NULL, NULL, "-da_grid_y", &cells[1], NULL);CHKERRQ(ierr);
386557cf195SMatthew G. Knepley   ierr = PetscOptionsGetInt(NULL, NULL, "-da_grid_z", &cells[2], NULL);CHKERRQ(ierr);
387557cf195SMatthew G. Knepley   ierr = DMPlexCreateBoxMesh(comm, user->dim, user->simplex, cells, NULL, NULL, NULL, user->interpolate, dm);CHKERRQ(ierr);
388557cf195SMatthew G. Knepley   /* Refine mesh using a volume constraint */
389557cf195SMatthew G. Knepley   if (user->simplex) {
390557cf195SMatthew G. Knepley     DM rdm = NULL;
391557cf195SMatthew G. Knepley 
392557cf195SMatthew G. Knepley     ierr = DMPlexSetRefinementLimit(*dm, user->refinementLimit);CHKERRQ(ierr);
393557cf195SMatthew G. Knepley     ierr = DMRefine(*dm, comm, &rdm);CHKERRQ(ierr);
394557cf195SMatthew G. Knepley     if (rdm) {
395557cf195SMatthew G. Knepley       ierr = DMDestroy(dm);CHKERRQ(ierr);
396557cf195SMatthew G. Knepley       *dm  = rdm;
397557cf195SMatthew G. Knepley     }
398557cf195SMatthew G. Knepley   }
399557cf195SMatthew G. Knepley   ierr = DMPlexSetRefinementUniform(*dm, PETSC_TRUE);CHKERRQ(ierr);
400557cf195SMatthew G. Knepley   /* Distribute mesh over processes */
401557cf195SMatthew G. Knepley   ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr);
402557cf195SMatthew G. Knepley   ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
403557cf195SMatthew G. Knepley   ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr);
404557cf195SMatthew G. Knepley   if (pdm) {
405557cf195SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
406557cf195SMatthew G. Knepley     *dm  = pdm;
407557cf195SMatthew G. Knepley   }
408557cf195SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) *dm, user->simplex ? "Simplicial Mesh" : "Hexahedral Mesh");CHKERRQ(ierr);
409557cf195SMatthew G. Knepley   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
410557cf195SMatthew G. Knepley   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
411557cf195SMatthew G. Knepley   PetscFunctionReturn(0);
412557cf195SMatthew G. Knepley }
413557cf195SMatthew G. Knepley 
414557cf195SMatthew G. Knepley /* Setup functions to approximate */
415557cf195SMatthew G. Knepley static PetscErrorCode SetupFunctions(DM dm, PetscBool usePoly, PetscInt order, PetscInt dir, PetscErrorCode (**exactFuncs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *),
416557cf195SMatthew G. Knepley                                      PetscErrorCode (**exactFuncDers)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), AppCtx *user)
417557cf195SMatthew G. Knepley {
418557cf195SMatthew G. Knepley   PetscInt       dim;
419557cf195SMatthew G. Knepley   PetscErrorCode ierr;
420557cf195SMatthew G. Knepley 
421557cf195SMatthew G. Knepley   PetscFunctionBeginUser;
422557cf195SMatthew G. Knepley   user->dir = dir;
423557cf195SMatthew G. Knepley   if (usePoly) {
424557cf195SMatthew G. Knepley     switch (order) {
425557cf195SMatthew G. Knepley     case 0:
426557cf195SMatthew G. Knepley       exactFuncs[0]    = constant;
427557cf195SMatthew G. Knepley       exactFuncDers[0] = constantDer;
428557cf195SMatthew G. Knepley       break;
429557cf195SMatthew G. Knepley     case 1:
430557cf195SMatthew G. Knepley       exactFuncs[0]    = linear;
431557cf195SMatthew G. Knepley       exactFuncDers[0] = linearDer;
432557cf195SMatthew G. Knepley       break;
433557cf195SMatthew G. Knepley     case 2:
434557cf195SMatthew G. Knepley       exactFuncs[0]    = quadratic;
435557cf195SMatthew G. Knepley       exactFuncDers[0] = quadraticDer;
436557cf195SMatthew G. Knepley       break;
437557cf195SMatthew G. Knepley     case 3:
438557cf195SMatthew G. Knepley       exactFuncs[0]    = cubic;
439557cf195SMatthew G. Knepley       exactFuncDers[0] = cubicDer;
440557cf195SMatthew G. Knepley       break;
441557cf195SMatthew G. Knepley     case 4:
442557cf195SMatthew G. Knepley       exactFuncs[0]    = quartic;
443557cf195SMatthew G. Knepley       exactFuncDers[0] = quarticDer;
444557cf195SMatthew G. Knepley       break;
445557cf195SMatthew G. Knepley     default:
446557cf195SMatthew G. Knepley       ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
447557cf195SMatthew G. Knepley       SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Could not determine functions to test for dimension %d order %d", dim, order);
448557cf195SMatthew G. Knepley     }
449557cf195SMatthew G. Knepley   } else {
450557cf195SMatthew G. Knepley     user->m          = order;
451557cf195SMatthew G. Knepley     exactFuncs[0]    = trig;
452557cf195SMatthew G. Knepley     exactFuncDers[0] = trigDer;
453557cf195SMatthew G. Knepley   }
454557cf195SMatthew G. Knepley   PetscFunctionReturn(0);
455557cf195SMatthew G. Knepley }
456557cf195SMatthew G. Knepley 
457557cf195SMatthew G. Knepley static PetscErrorCode ComputeError(DM dm, PetscErrorCode (**exactFuncs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *),
458557cf195SMatthew G. Knepley                                    PetscErrorCode (**exactFuncDers)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *),
459557cf195SMatthew G. Knepley                                    void **exactCtxs, PetscReal *error, PetscReal *errorDer, AppCtx *user)
460557cf195SMatthew G. Knepley {
461557cf195SMatthew G. Knepley   Vec            u;
462557cf195SMatthew G. Knepley   PetscReal      n[3] = {1.0, 1.0, 1.0};
463557cf195SMatthew G. Knepley   PetscErrorCode ierr;
464557cf195SMatthew G. Knepley 
465557cf195SMatthew G. Knepley   PetscFunctionBeginUser;
466557cf195SMatthew G. Knepley   ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr);
467557cf195SMatthew G. Knepley   /* Project function into FE function space */
468557cf195SMatthew G. Knepley   ierr = DMProjectFunction(dm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, u);CHKERRQ(ierr);
469557cf195SMatthew G. Knepley   ierr = VecViewFromOptions(u, NULL, "-projection_view");CHKERRQ(ierr);
470557cf195SMatthew G. Knepley   /* Compare approximation to exact in L_2 */
471557cf195SMatthew G. Knepley   ierr = DMComputeL2Diff(dm, 0.0, exactFuncs, exactCtxs, u, error);CHKERRQ(ierr);
472557cf195SMatthew G. Knepley   ierr = DMComputeL2GradientDiff(dm, 0.0, exactFuncDers, exactCtxs, u, n, errorDer);CHKERRQ(ierr);
473557cf195SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr);
474557cf195SMatthew G. Knepley   PetscFunctionReturn(0);
475557cf195SMatthew G. Knepley }
476557cf195SMatthew G. Knepley 
477557cf195SMatthew G. Knepley static PetscErrorCode CheckFunctions(DM dm, PetscInt order, AppCtx *user)
478557cf195SMatthew G. Knepley {
479557cf195SMatthew G. Knepley   PetscErrorCode (*exactFuncs[1]) (PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
480557cf195SMatthew G. Knepley   PetscErrorCode (*exactFuncDers[1]) (PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx);
481557cf195SMatthew G. Knepley   void            *exactCtxs[3];
482557cf195SMatthew G. Knepley   MPI_Comm         comm;
483557cf195SMatthew G. Knepley   PetscReal        error, errorDer, tol = PETSC_SMALL;
484557cf195SMatthew G. Knepley   PetscErrorCode   ierr;
485557cf195SMatthew G. Knepley 
486557cf195SMatthew G. Knepley   PetscFunctionBeginUser;
487557cf195SMatthew G. Knepley   exactCtxs[0]       = user;
488557cf195SMatthew G. Knepley   exactCtxs[1]       = user;
489557cf195SMatthew G. Knepley   exactCtxs[2]       = user;
490557cf195SMatthew G. Knepley   user->constants[0] = 1.0;
491557cf195SMatthew G. Knepley   user->constants[1] = 2.0;
492557cf195SMatthew G. Knepley   user->constants[2] = 3.0;
493557cf195SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
494557cf195SMatthew G. Knepley   ierr = SetupFunctions(dm, PETSC_TRUE, order, 0, exactFuncs, exactFuncDers, user);CHKERRQ(ierr);
495557cf195SMatthew G. Knepley   ierr = ComputeError(dm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user);CHKERRQ(ierr);
496557cf195SMatthew G. Knepley   /* Report result */
497557cf195SMatthew 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);}
498557cf195SMatthew G. Knepley   else                {ierr = PetscPrintf(comm, "Function tests pass for order %D at tolerance %g\n", order, (double)tol);CHKERRQ(ierr);}
499557cf195SMatthew 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);}
500557cf195SMatthew G. Knepley   else                {ierr = PetscPrintf(comm, "Function tests pass for order %D derivatives at tolerance %g\n", order, (double)tol);CHKERRQ(ierr);}
501557cf195SMatthew G. Knepley   PetscFunctionReturn(0);
502557cf195SMatthew G. Knepley }
503557cf195SMatthew G. Knepley 
504557cf195SMatthew G. Knepley /* Compare approximation to exact in L_2 */
505557cf195SMatthew G. Knepley static PetscErrorCode CheckTransferError(DM fdm, PetscBool usePoly, PetscInt order, PetscInt dir, const char *testname, Vec fu, AppCtx *user)
506557cf195SMatthew G. Knepley {
507557cf195SMatthew G. Knepley   PetscErrorCode (*exactFuncs[1]) (PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
508557cf195SMatthew G. Knepley   PetscErrorCode (*exactFuncDers[1]) (PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx);
509557cf195SMatthew G. Knepley   PetscReal        n[3] = {1.0, 1.0, 1.0};
510557cf195SMatthew G. Knepley   void            *exactCtxs[3];
511557cf195SMatthew G. Knepley   MPI_Comm         comm;
512557cf195SMatthew G. Knepley   PetscReal        error, errorDer, tol = PETSC_SMALL;
513557cf195SMatthew G. Knepley   PetscErrorCode   ierr;
514557cf195SMatthew G. Knepley 
515557cf195SMatthew G. Knepley   PetscFunctionBeginUser;
516557cf195SMatthew G. Knepley   exactCtxs[0]       = user;
517557cf195SMatthew G. Knepley   exactCtxs[1]       = user;
518557cf195SMatthew G. Knepley   exactCtxs[2]       = user;
519557cf195SMatthew G. Knepley   user->constants[0] = 1.0;
520557cf195SMatthew G. Knepley   user->constants[1] = 2.0;
521557cf195SMatthew G. Knepley   user->constants[2] = 3.0;
522557cf195SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) fdm, &comm);CHKERRQ(ierr);
523557cf195SMatthew G. Knepley   ierr = SetupFunctions(fdm, usePoly, order, dir, exactFuncs, exactFuncDers, user);CHKERRQ(ierr);
524557cf195SMatthew G. Knepley   ierr = DMComputeL2Diff(fdm, 0.0, exactFuncs, exactCtxs, fu, &error);CHKERRQ(ierr);
525557cf195SMatthew G. Knepley   ierr = DMComputeL2GradientDiff(fdm, 0.0, exactFuncDers, exactCtxs, fu, n, &errorDer);CHKERRQ(ierr);
526557cf195SMatthew G. Knepley   /* Report result */
527557cf195SMatthew 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);}
528557cf195SMatthew G. Knepley   else                {ierr = PetscPrintf(comm, "%s tests pass for order %D at tolerance %g\n", testname, order, (double)tol);CHKERRQ(ierr);}
529557cf195SMatthew 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);}
530557cf195SMatthew G. Knepley   else                {ierr = PetscPrintf(comm, "%s tests pass for order %D derivatives at tolerance %g\n", testname, order, (double)tol);CHKERRQ(ierr);}
531557cf195SMatthew G. Knepley   PetscFunctionReturn(0);
532557cf195SMatthew G. Knepley }
533557cf195SMatthew G. Knepley 
534557cf195SMatthew G. Knepley static PetscErrorCode CheckTransfer(DM dm, InterpType inType, PetscInt order, AppCtx *user)
535557cf195SMatthew G. Knepley {
536557cf195SMatthew G. Knepley   PetscErrorCode (*exactFuncs[1]) (PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
537557cf195SMatthew G. Knepley   PetscErrorCode (*exactFuncDers[1]) (PetscInt, PetscReal, const PetscReal x[], const PetscReal n[], PetscInt, PetscScalar *u, void *ctx);
538*d6837840SMatthew G. Knepley   void           *exactCtxs[3];
539*d6837840SMatthew G. Knepley   DM              rdm = NULL, idm = NULL, fdm = NULL;
540557cf195SMatthew G. Knepley   Mat             Interp, InterpAdapt = NULL;
541*d6837840SMatthew G. Knepley   Vec             iu, fu, scaling = NULL;
542557cf195SMatthew G. Knepley   MPI_Comm        comm;
543*d6837840SMatthew G. Knepley   const char     *testname = "Unknown";
544557cf195SMatthew G. Knepley   char            checkname[PETSC_MAX_PATH_LEN];
545557cf195SMatthew G. Knepley   PetscErrorCode  ierr;
546557cf195SMatthew G. Knepley 
547557cf195SMatthew G. Knepley   PetscFunctionBeginUser;
548*d6837840SMatthew G. Knepley   exactCtxs[0] = exactCtxs[1] = exactCtxs[2] = user;
549557cf195SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
550557cf195SMatthew G. Knepley   ierr = DMRefine(dm, comm, &rdm);CHKERRQ(ierr);
551557cf195SMatthew G. Knepley   ierr = DMSetCoarseDM(rdm, dm);CHKERRQ(ierr);
552557cf195SMatthew G. Knepley   ierr = DMCopyDisc(dm, rdm);CHKERRQ(ierr);
553557cf195SMatthew G. Knepley   switch (inType) {
554557cf195SMatthew G. Knepley   case INTERPOLATION:
555557cf195SMatthew G. Knepley     testname = "Interpolation";
556557cf195SMatthew G. Knepley     idm = dm;
557557cf195SMatthew G. Knepley     fdm = rdm;
558557cf195SMatthew G. Knepley     break;
559557cf195SMatthew G. Knepley   case RESTRICTION:
560557cf195SMatthew G. Knepley     testname = "Restriction";
561557cf195SMatthew G. Knepley     idm = rdm;
562557cf195SMatthew G. Knepley     fdm = dm;
563557cf195SMatthew G. Knepley     break;
564557cf195SMatthew G. Knepley   case INJECTION:
565557cf195SMatthew G. Knepley     testname = "Injection";
566557cf195SMatthew G. Knepley     idm = rdm;
567557cf195SMatthew G. Knepley     fdm = dm;
568557cf195SMatthew G. Knepley     break;
569557cf195SMatthew G. Knepley   }
570557cf195SMatthew G. Knepley   ierr = DMGetGlobalVector(idm, &iu);CHKERRQ(ierr);
571557cf195SMatthew G. Knepley   ierr = DMGetGlobalVector(fdm, &fu);CHKERRQ(ierr);
572557cf195SMatthew G. Knepley   ierr = DMSetApplicationContext(dm, user);CHKERRQ(ierr);
573557cf195SMatthew G. Knepley   ierr = DMSetApplicationContext(rdm, user);CHKERRQ(ierr);
574557cf195SMatthew G. Knepley   /* Project function into initial FE function space */
575557cf195SMatthew G. Knepley   ierr = SetupFunctions(dm, PETSC_TRUE, order, 0, exactFuncs, exactFuncDers, user);CHKERRQ(ierr);
576557cf195SMatthew G. Knepley   ierr = DMProjectFunction(idm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, iu);CHKERRQ(ierr);
577557cf195SMatthew G. Knepley   /* Interpolate function into final FE function space */
578557cf195SMatthew G. Knepley   switch (inType) {
579557cf195SMatthew G. Knepley   case INTERPOLATION:
580557cf195SMatthew G. Knepley     ierr = DMCreateInterpolation(dm, rdm, &Interp, &scaling);CHKERRQ(ierr);
581557cf195SMatthew G. Knepley     ierr = MatInterpolate(Interp, iu, fu);CHKERRQ(ierr);
582557cf195SMatthew G. Knepley     break;
583557cf195SMatthew G. Knepley   case RESTRICTION:
584557cf195SMatthew G. Knepley     ierr = DMCreateInterpolation(dm, rdm, &Interp, &scaling);CHKERRQ(ierr);
585557cf195SMatthew G. Knepley     ierr = MatRestrict(Interp, iu, fu);CHKERRQ(ierr);
586557cf195SMatthew G. Knepley     ierr = VecPointwiseMult(fu, scaling, fu);CHKERRQ(ierr);
587557cf195SMatthew G. Knepley     break;
588557cf195SMatthew G. Knepley   case INJECTION:
589557cf195SMatthew G. Knepley     ierr = DMCreateInjection(dm, rdm, &Interp);CHKERRQ(ierr);
590557cf195SMatthew G. Knepley     ierr = MatRestrict(Interp, iu, fu);CHKERRQ(ierr);
591557cf195SMatthew G. Knepley     break;
592557cf195SMatthew G. Knepley   }
593557cf195SMatthew G. Knepley   ierr = CheckTransferError(fdm, PETSC_TRUE, order, 0, testname, fu, user);CHKERRQ(ierr);
594557cf195SMatthew G. Knepley   if (user->K && (inType == INTERPOLATION)) {
595557cf195SMatthew G. Knepley     KSP      smoother;
596557cf195SMatthew G. Knepley     Mat      A;
597557cf195SMatthew G. Knepley     Vec     *iV, *fV;
598557cf195SMatthew G. Knepley     PetscInt k, dim, d;
599557cf195SMatthew G. Knepley 
600557cf195SMatthew G. Knepley     ierr = PetscPrintf(comm, " Adapting interpolator using %s\n", user->usePoly ? "polynomials" : "harmonics");CHKERRQ(ierr);
601557cf195SMatthew G. Knepley     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
602557cf195SMatthew G. Knepley     ierr = PetscMalloc2(user->K*dim, &iV, user->K*dim, &fV);CHKERRQ(ierr);
603557cf195SMatthew G. Knepley     /* Project coarse modes into initial and final FE function space */
604557cf195SMatthew G. Knepley     for (k = 0; k < user->K; ++k) {
605557cf195SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
606557cf195SMatthew G. Knepley         ierr = DMGetGlobalVector(idm, &iV[k*dim+d]);CHKERRQ(ierr);
607557cf195SMatthew G. Knepley         ierr = DMGetGlobalVector(fdm, &fV[k*dim+d]);CHKERRQ(ierr);
608557cf195SMatthew G. Knepley         ierr = SetupFunctions(idm, user->usePoly, user->usePoly ? k : k+1, d, exactFuncs, exactFuncDers, user);CHKERRQ(ierr);
609557cf195SMatthew G. Knepley         ierr = DMProjectFunction(idm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, iV[k*dim+d]);CHKERRQ(ierr);
610557cf195SMatthew G. Knepley         ierr = DMProjectFunction(fdm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, fV[k*dim+d]);CHKERRQ(ierr);
611557cf195SMatthew G. Knepley       }
612557cf195SMatthew G. Knepley     }
613557cf195SMatthew G. Knepley     /* Adapt interpolator */
614557cf195SMatthew G. Knepley     ierr = DMCreateMatrix(rdm, &A);CHKERRQ(ierr);
615557cf195SMatthew G. Knepley     ierr = MatShift(A, 1.0);CHKERRQ(ierr);
616557cf195SMatthew G. Knepley     ierr = KSPCreate(comm, &smoother);CHKERRQ(ierr);
617557cf195SMatthew G. Knepley     ierr = KSPSetFromOptions(smoother);CHKERRQ(ierr);
618557cf195SMatthew G. Knepley     ierr = KSPSetOperators(smoother, A, A);CHKERRQ(ierr);
619557cf195SMatthew G. Knepley     ierr = DMAdaptInterpolator(dm, rdm, Interp, smoother, user->K*dim, fV, iV, &InterpAdapt, user);CHKERRQ(ierr);
620557cf195SMatthew G. Knepley     /* Interpolate function into final FE function space */
621557cf195SMatthew G. Knepley     ierr = PetscSNPrintf(checkname, PETSC_MAX_PATH_LEN, "  %s poly", testname);CHKERRQ(ierr);
622557cf195SMatthew G. Knepley     ierr = MatInterpolate(InterpAdapt, iu, fu);CHKERRQ(ierr);
623557cf195SMatthew G. Knepley     ierr = CheckTransferError(fdm, PETSC_TRUE, order, 0, checkname, fu, user);CHKERRQ(ierr);
624557cf195SMatthew G. Knepley     for (k = 0; k < user->K; ++k) {
625557cf195SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
626557cf195SMatthew G. Knepley         ierr = PetscSNPrintf(checkname, PETSC_MAX_PATH_LEN, "  %s trig (%D, %D)", testname, k, d);CHKERRQ(ierr);
627557cf195SMatthew G. Knepley         ierr = MatInterpolate(InterpAdapt, iV[k*dim+d], fV[k*dim+d]);CHKERRQ(ierr);
628557cf195SMatthew G. Knepley         ierr = CheckTransferError(fdm, PETSC_FALSE, k+1, d, checkname, fV[k*dim+d], user);CHKERRQ(ierr);
629557cf195SMatthew G. Knepley       }
630557cf195SMatthew G. Knepley     }
631557cf195SMatthew G. Knepley     /* Cleanup */
632557cf195SMatthew G. Knepley     ierr = KSPDestroy(&smoother);CHKERRQ(ierr);
633557cf195SMatthew G. Knepley     ierr = MatDestroy(&A);CHKERRQ(ierr);
634557cf195SMatthew G. Knepley     for (k = 0; k < user->K; ++k) {
635557cf195SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
636557cf195SMatthew G. Knepley         ierr = DMRestoreGlobalVector(idm, &iV[k*dim+d]);CHKERRQ(ierr);
637557cf195SMatthew G. Knepley         ierr = DMRestoreGlobalVector(fdm, &fV[k*dim+d]);CHKERRQ(ierr);
638557cf195SMatthew G. Knepley       }
639557cf195SMatthew G. Knepley     }
640557cf195SMatthew G. Knepley     ierr = PetscFree2(iV, fV);CHKERRQ(ierr);
641557cf195SMatthew G. Knepley     ierr = MatDestroy(&InterpAdapt);CHKERRQ(ierr);
642557cf195SMatthew G. Knepley   }
643557cf195SMatthew G. Knepley   ierr = DMRestoreGlobalVector(idm, &iu);CHKERRQ(ierr);
644557cf195SMatthew G. Knepley   ierr = DMRestoreGlobalVector(fdm, &fu);CHKERRQ(ierr);
645557cf195SMatthew G. Knepley   ierr = MatDestroy(&Interp);CHKERRQ(ierr);
646557cf195SMatthew G. Knepley   ierr = VecDestroy(&scaling);CHKERRQ(ierr);
647557cf195SMatthew G. Knepley   ierr = DMDestroy(&rdm);CHKERRQ(ierr);
648557cf195SMatthew G. Knepley   PetscFunctionReturn(0);
649557cf195SMatthew G. Knepley }
650557cf195SMatthew G. Knepley 
651557cf195SMatthew G. Knepley int main(int argc, char **argv)
652557cf195SMatthew G. Knepley {
653557cf195SMatthew G. Knepley   DM             dm;
654557cf195SMatthew G. Knepley   AppCtx         user;                 /* user-defined work context */
655557cf195SMatthew G. Knepley   PetscErrorCode ierr;
656557cf195SMatthew G. Knepley 
657557cf195SMatthew G. Knepley   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
658557cf195SMatthew G. Knepley   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
659557cf195SMatthew G. Knepley   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
660557cf195SMatthew G. Knepley   ierr = PetscFECreateDefault(PETSC_COMM_WORLD, user.dim, user.Nc, user.simplex, NULL, user.qorder, &user.fe);CHKERRQ(ierr);
661557cf195SMatthew G. Knepley   ierr = DMSetField(dm, 0, NULL, (PetscObject) user.fe);CHKERRQ(ierr);
662557cf195SMatthew G. Knepley   ierr = DMCreateDS(dm);CHKERRQ(ierr);
663557cf195SMatthew G. Knepley   ierr = CheckFunctions(dm, user.porder, &user);CHKERRQ(ierr);
664557cf195SMatthew G. Knepley   ierr = CheckTransfer(dm, INTERPOLATION, user.porder, &user);CHKERRQ(ierr);
665557cf195SMatthew G. Knepley   ierr = CheckTransfer(dm, INJECTION,  user.porder, &user);CHKERRQ(ierr);
666557cf195SMatthew G. Knepley   ierr = PetscFEDestroy(&user.fe);CHKERRQ(ierr);
667557cf195SMatthew G. Knepley   ierr = DMDestroy(&dm);CHKERRQ(ierr);
668557cf195SMatthew G. Knepley   ierr = PetscFinalize();
669557cf195SMatthew G. Knepley   return ierr;
670557cf195SMatthew G. Knepley }
671557cf195SMatthew G. Knepley 
672557cf195SMatthew G. Knepley /*TEST
673557cf195SMatthew G. Knepley 
674557cf195SMatthew G. Knepley   # TODO dim 3 will not work until I get composite elements in 3D (see plexrefine.c:34)
675557cf195SMatthew G. Knepley   # 2D/3D P_1 on a simplex
676557cf195SMatthew G. Knepley   test:
677557cf195SMatthew G. Knepley     suffix: p1
678557cf195SMatthew G. Knepley     requires: triangle ctetgen
679557cf195SMatthew G. Knepley     args: -dim {{2}separate output} -petscspace_degree 1 -num_comp 1 -qorder 1 -porder {{1}separate output}
680557cf195SMatthew G. Knepley   test:
681557cf195SMatthew G. Knepley     suffix: p1_pragmatic
682557cf195SMatthew G. Knepley     requires: triangle ctetgen pragmatic
683557cf195SMatthew G. Knepley     args: -dim {{2}separate output} -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder {{1 2}separate output}
684557cf195SMatthew G. Knepley   test:
685557cf195SMatthew G. Knepley     suffix: p1_adapt
686557cf195SMatthew G. Knepley     requires: triangle ctetgen
687557cf195SMatthew G. Knepley     args: -dim {{2}separate output} -dm_refine 3 -petscspace_degree 1 -qorder 1 -porder {{1 2}separate output}
688557cf195SMatthew G. Knepley 
689557cf195SMatthew G. Knepley   # TODO dim 3 will not work until I get composite elements in 3D (see plexrefine.c:34)
690557cf195SMatthew G. Knepley   # 2D/3D P_2 on a simplex
691557cf195SMatthew G. Knepley   test:
692557cf195SMatthew G. Knepley     suffix: p2
693557cf195SMatthew G. Knepley     requires: triangle ctetgen
694557cf195SMatthew G. Knepley     args: -dim {{2}separate output} -petscspace_degree 2 -qorder 2 -porder {{1 2 3}separate output}
695557cf195SMatthew G. Knepley   test:
696557cf195SMatthew G. Knepley     suffix: p2_pragmatic
697557cf195SMatthew G. Knepley     requires: triangle ctetgen pragmatic
698557cf195SMatthew G. Knepley     args: -dim {{2}separate output} -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder {{1 2 3}separate output}
699557cf195SMatthew G. Knepley 
700557cf195SMatthew G. Knepley   # TODO dim 3 will not work until I get composite elements in 3D (see plexrefine.c:34)
701557cf195SMatthew G. Knepley   # TODO This is broken. Check ex3 which worked
702557cf195SMatthew G. Knepley   # 2D/3D P_3 on a simplex
703557cf195SMatthew G. Knepley   test:
704*d6837840SMatthew G. Knepley     TODO: gll Lagrange nodes break this
705557cf195SMatthew G. Knepley     suffix: p3
706557cf195SMatthew G. Knepley     requires: triangle ctetgen !single
707557cf195SMatthew G. Knepley     args: -dim {{2}separate output} -petscspace_degree 3 -qorder 3 -porder {{1 2 3 4}separate output}
708557cf195SMatthew G. Knepley   test:
709*d6837840SMatthew G. Knepley     TODO: gll Lagrange nodes break this
710557cf195SMatthew G. Knepley     suffix: p3_pragmatic
711557cf195SMatthew G. Knepley     requires: triangle ctetgen pragmatic !single
712557cf195SMatthew G. Knepley     args: -dim {{2}separate output} -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -porder {{1 2 3 4}separate output}
713557cf195SMatthew G. Knepley 
714557cf195SMatthew G. Knepley   # 2D/3D Q_1 on a tensor cell
715557cf195SMatthew G. Knepley   test:
716557cf195SMatthew G. Knepley     suffix: q1
717557cf195SMatthew G. Knepley     requires: mpi_type_get_envelope
718557cf195SMatthew G. Knepley     args: -dim {{2 3}separate output} -simplex 0 -petscspace_degree 1 -qorder 1 -porder {{1 2}separate output}
719557cf195SMatthew G. Knepley 
720557cf195SMatthew G. Knepley   # 2D/3D Q_2 on a tensor cell
721557cf195SMatthew G. Knepley   test:
722557cf195SMatthew G. Knepley     suffix: q2
723*d6837840SMatthew G. Knepley     requires: mpi_type_get_envelope !single
724557cf195SMatthew G. Knepley     args: -dim {{2 3}separate output} -simplex 0 -petscspace_degree 2 -qorder 2 -porder {{1 2 3}separate output}
725557cf195SMatthew G. Knepley 
726557cf195SMatthew G. Knepley   # 2D/3D Q_3 on a tensor cell
727557cf195SMatthew G. Knepley   test:
728*d6837840SMatthew G. Knepley     TODO: gll Lagrange nodes break this
729557cf195SMatthew G. Knepley     suffix: q3
730557cf195SMatthew G. Knepley     requires: mpi_type_get_envelope !single
731557cf195SMatthew G. Knepley     args: -dim {{2 3}separate output} -simplex 0 -petscspace_degree 3 -qorder 3 -porder {{1 2 3 4}separate output}
732557cf195SMatthew G. Knepley 
733557cf195SMatthew G. Knepley   # 2D/3D P_1disc on a triangle/quadrilateral
734557cf195SMatthew G. Knepley   # TODO Missing injection functional for simplices
735557cf195SMatthew G. Knepley   test:
736557cf195SMatthew G. Knepley     suffix: p1d
737557cf195SMatthew G. Knepley     requires: triangle ctetgen
738557cf195SMatthew 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}
739557cf195SMatthew G. Knepley 
740557cf195SMatthew G. Knepley TEST*/
741