xref: /petsc/src/dm/tests/ex36.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Checks the functionality of DMGetInterpolation() on deformed grids.\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscdm.h>
5c4762a1bSJed Brown #include <petscdmda.h>
6c4762a1bSJed Brown 
7c4762a1bSJed Brown typedef struct _n_CCmplx CCmplx;
8c4762a1bSJed Brown struct _n_CCmplx {
9c4762a1bSJed Brown   PetscReal real;
10c4762a1bSJed Brown   PetscReal imag;
11c4762a1bSJed Brown };
12c4762a1bSJed Brown 
139371c9d4SSatish Balay CCmplx CCmplxPow(CCmplx a, PetscReal n) {
14c4762a1bSJed Brown   CCmplx    b;
15c4762a1bSJed Brown   PetscReal r, theta;
16c4762a1bSJed Brown   r      = PetscSqrtReal(a.real * a.real + a.imag * a.imag);
17c4762a1bSJed Brown   theta  = PetscAtan2Real(a.imag, a.real);
18c4762a1bSJed Brown   b.real = PetscPowReal(r, n) * PetscCosReal(n * theta);
19c4762a1bSJed Brown   b.imag = PetscPowReal(r, n) * PetscSinReal(n * theta);
20c4762a1bSJed Brown   return b;
21c4762a1bSJed Brown }
229371c9d4SSatish Balay CCmplx CCmplxExp(CCmplx a) {
23c4762a1bSJed Brown   CCmplx b;
24c4762a1bSJed Brown   b.real = PetscExpReal(a.real) * PetscCosReal(a.imag);
25c4762a1bSJed Brown   b.imag = PetscExpReal(a.real) * PetscSinReal(a.imag);
26c4762a1bSJed Brown   return b;
27c4762a1bSJed Brown }
289371c9d4SSatish Balay CCmplx CCmplxSqrt(CCmplx a) {
29c4762a1bSJed Brown   CCmplx    b;
30c4762a1bSJed Brown   PetscReal r, theta;
31c4762a1bSJed Brown   r      = PetscSqrtReal(a.real * a.real + a.imag * a.imag);
32c4762a1bSJed Brown   theta  = PetscAtan2Real(a.imag, a.real);
33c4762a1bSJed Brown   b.real = PetscSqrtReal(r) * PetscCosReal(0.5 * theta);
34c4762a1bSJed Brown   b.imag = PetscSqrtReal(r) * PetscSinReal(0.5 * theta);
35c4762a1bSJed Brown   return b;
36c4762a1bSJed Brown }
379371c9d4SSatish Balay CCmplx CCmplxAdd(CCmplx a, CCmplx c) {
38c4762a1bSJed Brown   CCmplx b;
39c4762a1bSJed Brown   b.real = a.real + c.real;
40c4762a1bSJed Brown   b.imag = a.imag + c.imag;
41c4762a1bSJed Brown   return b;
42c4762a1bSJed Brown }
439371c9d4SSatish Balay PetscScalar CCmplxRe(CCmplx a) {
44c4762a1bSJed Brown   return (PetscScalar)a.real;
45c4762a1bSJed Brown }
469371c9d4SSatish Balay PetscScalar CCmplxIm(CCmplx a) {
47c4762a1bSJed Brown   return (PetscScalar)a.imag;
48c4762a1bSJed Brown }
49c4762a1bSJed Brown 
509371c9d4SSatish Balay PetscErrorCode DAApplyConformalMapping(DM da, PetscInt idx) {
51c4762a1bSJed Brown   PetscInt     i, n;
52c4762a1bSJed Brown   PetscInt     sx, nx, sy, ny, sz, nz, dim;
53c4762a1bSJed Brown   Vec          Gcoords;
54c4762a1bSJed Brown   PetscScalar *XX;
55c4762a1bSJed Brown   PetscScalar  xx, yy, zz;
56c4762a1bSJed Brown   DM           cda;
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   PetscFunctionBeginUser;
59c4762a1bSJed Brown   if (idx == 0) {
60c4762a1bSJed Brown     PetscFunctionReturn(0);
61c4762a1bSJed Brown   } else if (idx == 1) { /* dam break */
629566063dSJacob Faibussowitsch     PetscCall(DMDASetUniformCoordinates(da, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0));
63c4762a1bSJed Brown   } else if (idx == 2) { /* stagnation in a corner */
649566063dSJacob Faibussowitsch     PetscCall(DMDASetUniformCoordinates(da, -1.0, 1.0, 0.0, 1.0, -1.0, 1.0));
65c4762a1bSJed Brown   } else if (idx == 3) { /* nautilis */
669566063dSJacob Faibussowitsch     PetscCall(DMDASetUniformCoordinates(da, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0));
671baa6e33SBarry Smith   } else if (idx == 4) PetscCall(DMDASetUniformCoordinates(da, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0));
68c4762a1bSJed Brown 
699566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(da, &cda));
709566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(da, &Gcoords));
71c4762a1bSJed Brown 
729566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Gcoords, &XX));
739566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &sx, &sy, &sz, &nx, &ny, &nz));
749566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
759566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(Gcoords, &n));
76c4762a1bSJed Brown   n = n / dim;
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   for (i = 0; i < n; i++) {
79c4762a1bSJed Brown     if ((dim == 3) && (idx != 2)) {
80c4762a1bSJed Brown       PetscScalar Ni[8];
81c4762a1bSJed Brown       PetscScalar xi   = XX[dim * i];
82c4762a1bSJed Brown       PetscScalar eta  = XX[dim * i + 1];
83c4762a1bSJed Brown       PetscScalar zeta = XX[dim * i + 2];
84c4762a1bSJed Brown       PetscScalar xn[] = {-1.0, 1.0, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0};
85c4762a1bSJed Brown       PetscScalar yn[] = {-1.0, -1.0, 1.0, 1.0, -1.0, -1.0, 1.0, 1.0};
86c4762a1bSJed Brown       PetscScalar zn[] = {-0.1, -4.0, -0.2, -1.0, 0.1, 4.0, 0.2, 1.0};
87c4762a1bSJed Brown       PetscInt    p;
88c4762a1bSJed Brown 
89c4762a1bSJed Brown       Ni[0] = 0.125 * (1.0 - xi) * (1.0 - eta) * (1.0 - zeta);
90c4762a1bSJed Brown       Ni[1] = 0.125 * (1.0 + xi) * (1.0 - eta) * (1.0 - zeta);
91c4762a1bSJed Brown       Ni[2] = 0.125 * (1.0 - xi) * (1.0 + eta) * (1.0 - zeta);
92c4762a1bSJed Brown       Ni[3] = 0.125 * (1.0 + xi) * (1.0 + eta) * (1.0 - zeta);
93c4762a1bSJed Brown 
94c4762a1bSJed Brown       Ni[4] = 0.125 * (1.0 - xi) * (1.0 - eta) * (1.0 + zeta);
95c4762a1bSJed Brown       Ni[5] = 0.125 * (1.0 + xi) * (1.0 - eta) * (1.0 + zeta);
96c4762a1bSJed Brown       Ni[6] = 0.125 * (1.0 - xi) * (1.0 + eta) * (1.0 + zeta);
97c4762a1bSJed Brown       Ni[7] = 0.125 * (1.0 + xi) * (1.0 + eta) * (1.0 + zeta);
98c4762a1bSJed Brown 
99c4762a1bSJed Brown       xx = yy = zz = 0.0;
100c4762a1bSJed Brown       for (p = 0; p < 8; p++) {
101c4762a1bSJed Brown         xx += Ni[p] * xn[p];
102c4762a1bSJed Brown         yy += Ni[p] * yn[p];
103c4762a1bSJed Brown         zz += Ni[p] * zn[p];
104c4762a1bSJed Brown       }
105c4762a1bSJed Brown       XX[dim * i]     = xx;
106c4762a1bSJed Brown       XX[dim * i + 1] = yy;
107c4762a1bSJed Brown       XX[dim * i + 2] = zz;
108c4762a1bSJed Brown     }
109c4762a1bSJed Brown 
110c4762a1bSJed Brown     if (idx == 1) {
111c4762a1bSJed Brown       CCmplx zeta, t1, t2;
112c4762a1bSJed Brown 
113c4762a1bSJed Brown       xx = XX[dim * i] - 0.8;
114c4762a1bSJed Brown       yy = XX[dim * i + 1] + 1.5;
115c4762a1bSJed Brown 
116c4762a1bSJed Brown       zeta.real = PetscRealPart(xx);
117c4762a1bSJed Brown       zeta.imag = PetscRealPart(yy);
118c4762a1bSJed Brown 
119c4762a1bSJed Brown       t1 = CCmplxPow(zeta, -1.0);
120c4762a1bSJed Brown       t2 = CCmplxAdd(zeta, t1);
121c4762a1bSJed Brown 
122c4762a1bSJed Brown       XX[dim * i]     = CCmplxRe(t2);
123c4762a1bSJed Brown       XX[dim * i + 1] = CCmplxIm(t2);
124c4762a1bSJed Brown     } else if (idx == 2) {
125c4762a1bSJed Brown       CCmplx zeta, t1;
126c4762a1bSJed Brown 
127c4762a1bSJed Brown       xx        = XX[dim * i];
128c4762a1bSJed Brown       yy        = XX[dim * i + 1];
129c4762a1bSJed Brown       zeta.real = PetscRealPart(xx);
130c4762a1bSJed Brown       zeta.imag = PetscRealPart(yy);
131c4762a1bSJed Brown 
132c4762a1bSJed Brown       t1              = CCmplxSqrt(zeta);
133c4762a1bSJed Brown       XX[dim * i]     = CCmplxRe(t1);
134c4762a1bSJed Brown       XX[dim * i + 1] = CCmplxIm(t1);
135c4762a1bSJed Brown     } else if (idx == 3) {
136c4762a1bSJed Brown       CCmplx zeta, t1, t2;
137c4762a1bSJed Brown 
138c4762a1bSJed Brown       xx = XX[dim * i] - 0.8;
139c4762a1bSJed Brown       yy = XX[dim * i + 1] + 1.5;
140c4762a1bSJed Brown 
141c4762a1bSJed Brown       zeta.real       = PetscRealPart(xx);
142c4762a1bSJed Brown       zeta.imag       = PetscRealPart(yy);
143c4762a1bSJed Brown       t1              = CCmplxPow(zeta, -1.0);
144c4762a1bSJed Brown       t2              = CCmplxAdd(zeta, t1);
145c4762a1bSJed Brown       XX[dim * i]     = CCmplxRe(t2);
146c4762a1bSJed Brown       XX[dim * i + 1] = CCmplxIm(t2);
147c4762a1bSJed Brown 
148c4762a1bSJed Brown       xx              = XX[dim * i];
149c4762a1bSJed Brown       yy              = XX[dim * i + 1];
150c4762a1bSJed Brown       zeta.real       = PetscRealPart(xx);
151c4762a1bSJed Brown       zeta.imag       = PetscRealPart(yy);
152c4762a1bSJed Brown       t1              = CCmplxExp(zeta);
153c4762a1bSJed Brown       XX[dim * i]     = CCmplxRe(t1);
154c4762a1bSJed Brown       XX[dim * i + 1] = CCmplxIm(t1);
155c4762a1bSJed Brown 
156c4762a1bSJed Brown       xx              = XX[dim * i] + 0.4;
157c4762a1bSJed Brown       yy              = XX[dim * i + 1];
158c4762a1bSJed Brown       zeta.real       = PetscRealPart(xx);
159c4762a1bSJed Brown       zeta.imag       = PetscRealPart(yy);
160c4762a1bSJed Brown       t1              = CCmplxPow(zeta, 2.0);
161c4762a1bSJed Brown       XX[dim * i]     = CCmplxRe(t1);
162c4762a1bSJed Brown       XX[dim * i + 1] = CCmplxIm(t1);
163c4762a1bSJed Brown     } else if (idx == 4) {
164c4762a1bSJed Brown       PetscScalar Ni[4];
165c4762a1bSJed Brown       PetscScalar xi   = XX[dim * i];
166c4762a1bSJed Brown       PetscScalar eta  = XX[dim * i + 1];
167c4762a1bSJed Brown       PetscScalar xn[] = {0.0, 2.0, 0.2, 3.5};
168c4762a1bSJed Brown       PetscScalar yn[] = {-1.3, 0.0, 2.0, 4.0};
169c4762a1bSJed Brown       PetscInt    p;
170c4762a1bSJed Brown 
171c4762a1bSJed Brown       Ni[0] = 0.25 * (1.0 - xi) * (1.0 - eta);
172c4762a1bSJed Brown       Ni[1] = 0.25 * (1.0 + xi) * (1.0 - eta);
173c4762a1bSJed Brown       Ni[2] = 0.25 * (1.0 - xi) * (1.0 + eta);
174c4762a1bSJed Brown       Ni[3] = 0.25 * (1.0 + xi) * (1.0 + eta);
175c4762a1bSJed Brown 
176c4762a1bSJed Brown       xx = yy = 0.0;
177c4762a1bSJed Brown       for (p = 0; p < 4; p++) {
178c4762a1bSJed Brown         xx += Ni[p] * xn[p];
179c4762a1bSJed Brown         yy += Ni[p] * yn[p];
180c4762a1bSJed Brown       }
181c4762a1bSJed Brown       XX[dim * i]     = xx;
182c4762a1bSJed Brown       XX[dim * i + 1] = yy;
183c4762a1bSJed Brown     }
184c4762a1bSJed Brown   }
1859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Gcoords, &XX));
186c4762a1bSJed Brown   PetscFunctionReturn(0);
187c4762a1bSJed Brown }
188c4762a1bSJed Brown 
1899371c9d4SSatish Balay PetscErrorCode DAApplyTrilinearMapping(DM da) {
190c4762a1bSJed Brown   PetscInt      i, j, k;
191c4762a1bSJed Brown   PetscInt      sx, nx, sy, ny, sz, nz;
192c4762a1bSJed Brown   Vec           Gcoords;
193c4762a1bSJed Brown   DMDACoor3d ***XX;
194c4762a1bSJed Brown   PetscScalar   xx, yy, zz;
195c4762a1bSJed Brown   DM            cda;
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   PetscFunctionBeginUser;
1989566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0));
1999566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(da, &cda));
2009566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(da, &Gcoords));
201c4762a1bSJed Brown 
2029566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(cda, Gcoords, &XX));
2039566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &sx, &sy, &sz, &nx, &ny, &nz));
204c4762a1bSJed Brown 
205c4762a1bSJed Brown   for (i = sx; i < sx + nx; i++) {
206c4762a1bSJed Brown     for (j = sy; j < sy + ny; j++) {
207c4762a1bSJed Brown       for (k = sz; k < sz + nz; k++) {
208c4762a1bSJed Brown         PetscScalar Ni[8];
209c4762a1bSJed Brown         PetscScalar xi   = XX[k][j][i].x;
210c4762a1bSJed Brown         PetscScalar eta  = XX[k][j][i].y;
211c4762a1bSJed Brown         PetscScalar zeta = XX[k][j][i].z;
212c4762a1bSJed Brown         PetscScalar xn[] = {0.0, 2.0, 0.2, 3.5, 0.0, 2.1, 0.23, 3.125};
213c4762a1bSJed Brown         PetscScalar yn[] = {-1.3, 0.0, 2.0, 4.0, -1.45, -0.1, 2.24, 3.79};
214c4762a1bSJed Brown         PetscScalar zn[] = {0.0, 0.3, -0.1, 0.123, 0.956, 1.32, 1.12, 0.798};
215c4762a1bSJed Brown         PetscInt    p;
216c4762a1bSJed Brown 
217c4762a1bSJed Brown         Ni[0] = 0.125 * (1.0 - xi) * (1.0 - eta) * (1.0 - zeta);
218c4762a1bSJed Brown         Ni[1] = 0.125 * (1.0 + xi) * (1.0 - eta) * (1.0 - zeta);
219c4762a1bSJed Brown         Ni[2] = 0.125 * (1.0 - xi) * (1.0 + eta) * (1.0 - zeta);
220c4762a1bSJed Brown         Ni[3] = 0.125 * (1.0 + xi) * (1.0 + eta) * (1.0 - zeta);
221c4762a1bSJed Brown 
222c4762a1bSJed Brown         Ni[4] = 0.125 * (1.0 - xi) * (1.0 - eta) * (1.0 + zeta);
223c4762a1bSJed Brown         Ni[5] = 0.125 * (1.0 + xi) * (1.0 - eta) * (1.0 + zeta);
224c4762a1bSJed Brown         Ni[6] = 0.125 * (1.0 - xi) * (1.0 + eta) * (1.0 + zeta);
225c4762a1bSJed Brown         Ni[7] = 0.125 * (1.0 + xi) * (1.0 + eta) * (1.0 + zeta);
226c4762a1bSJed Brown 
227c4762a1bSJed Brown         xx = yy = zz = 0.0;
228c4762a1bSJed Brown         for (p = 0; p < 8; p++) {
229c4762a1bSJed Brown           xx += Ni[p] * xn[p];
230c4762a1bSJed Brown           yy += Ni[p] * yn[p];
231c4762a1bSJed Brown           zz += Ni[p] * zn[p];
232c4762a1bSJed Brown         }
233c4762a1bSJed Brown         XX[k][j][i].x = xx;
234c4762a1bSJed Brown         XX[k][j][i].y = yy;
235c4762a1bSJed Brown         XX[k][j][i].z = zz;
236c4762a1bSJed Brown       }
237c4762a1bSJed Brown     }
238c4762a1bSJed Brown   }
2399566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(cda, Gcoords, &XX));
240c4762a1bSJed Brown   PetscFunctionReturn(0);
241c4762a1bSJed Brown }
242c4762a1bSJed Brown 
2439371c9d4SSatish Balay PetscErrorCode DADefineXLinearField2D(DM da, Vec field) {
244c4762a1bSJed Brown   PetscInt      i, j;
245c4762a1bSJed Brown   PetscInt      sx, nx, sy, ny;
246c4762a1bSJed Brown   Vec           Gcoords;
247c4762a1bSJed Brown   DMDACoor2d  **XX;
248c4762a1bSJed Brown   PetscScalar **FF;
249c4762a1bSJed Brown   DM            cda;
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   PetscFunctionBeginUser;
2529566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(da, &cda));
2539566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(da, &Gcoords));
254c4762a1bSJed Brown 
2559566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(cda, Gcoords, &XX));
2569566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, field, &FF));
257c4762a1bSJed Brown 
2589566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &sx, &sy, 0, &nx, &ny, 0));
259c4762a1bSJed Brown 
260c4762a1bSJed Brown   for (i = sx; i < sx + nx; i++) {
2619371c9d4SSatish Balay     for (j = sy; j < sy + ny; j++) { FF[j][i] = 10.0 + 3.0 * XX[j][i].x + 5.5 * XX[j][i].y + 8.003 * XX[j][i].x * XX[j][i].y; }
262c4762a1bSJed Brown   }
263c4762a1bSJed Brown 
2649566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, field, &FF));
2659566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(cda, Gcoords, &XX));
266c4762a1bSJed Brown   PetscFunctionReturn(0);
267c4762a1bSJed Brown }
268c4762a1bSJed Brown 
2699371c9d4SSatish Balay PetscErrorCode DADefineXLinearField3D(DM da, Vec field) {
270c4762a1bSJed Brown   PetscInt       i, j, k;
271c4762a1bSJed Brown   PetscInt       sx, nx, sy, ny, sz, nz;
272c4762a1bSJed Brown   Vec            Gcoords;
273c4762a1bSJed Brown   DMDACoor3d  ***XX;
274c4762a1bSJed Brown   PetscScalar ***FF;
275c4762a1bSJed Brown   DM             cda;
276c4762a1bSJed Brown 
277c4762a1bSJed Brown   PetscFunctionBeginUser;
2789566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(da, &cda));
2799566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(da, &Gcoords));
280c4762a1bSJed Brown 
2819566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(cda, Gcoords, &XX));
2829566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, field, &FF));
283c4762a1bSJed Brown 
2849566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &sx, &sy, &sz, &nx, &ny, &nz));
285c4762a1bSJed Brown 
286c4762a1bSJed Brown   for (k = sz; k < sz + nz; k++) {
287c4762a1bSJed Brown     for (j = sy; j < sy + ny; j++) {
288c4762a1bSJed Brown       for (i = sx; i < sx + nx; i++) {
2899371c9d4SSatish Balay         FF[k][j][i] = 10.0 + 4.05 * XX[k][j][i].x + 5.50 * XX[k][j][i].y + 1.33 * XX[k][j][i].z + 2.03 * XX[k][j][i].x * XX[k][j][i].y + 0.03 * XX[k][j][i].x * XX[k][j][i].z + 0.83 * XX[k][j][i].y * XX[k][j][i].z +
2909371c9d4SSatish Balay                       3.79 * XX[k][j][i].x * XX[k][j][i].y * XX[k][j][i].z;
291c4762a1bSJed Brown       }
292c4762a1bSJed Brown     }
293c4762a1bSJed Brown   }
294c4762a1bSJed Brown 
2959566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, field, &FF));
2969566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(cda, Gcoords, &XX));
297c4762a1bSJed Brown   PetscFunctionReturn(0);
298c4762a1bSJed Brown }
299c4762a1bSJed Brown 
3009371c9d4SSatish Balay PetscErrorCode da_test_RefineCoords1D(PetscInt mx) {
301c4762a1bSJed Brown   DM          dac, daf;
302c4762a1bSJed Brown   PetscViewer vv;
303c4762a1bSJed Brown   Vec         ac, af;
304c4762a1bSJed Brown   PetscInt    Mx;
305c4762a1bSJed Brown   Mat         II, INTERP;
306c4762a1bSJed Brown   Vec         scale;
307c4762a1bSJed Brown   PetscBool   output = PETSC_FALSE;
308c4762a1bSJed Brown 
309c4762a1bSJed Brown   PetscFunctionBeginUser;
3109566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, mx + 1, 1, /* 1 dof */ 1, /* stencil = 1 */ NULL, &dac));
3119566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dac));
3129566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dac));
313c4762a1bSJed Brown 
3149566063dSJacob Faibussowitsch   PetscCall(DMRefine(dac, MPI_COMM_NULL, &daf));
3159566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf, 0, &Mx, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
316c4762a1bSJed Brown   Mx--;
317c4762a1bSJed Brown 
3189566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(dac, -1.0, 1.0, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE));
3199566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(daf, -1.0, 1.0, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE));
320c4762a1bSJed Brown 
321c4762a1bSJed Brown   {
322c4762a1bSJed Brown     DM  cdaf, cdac;
323c4762a1bSJed Brown     Vec coordsc, coordsf;
324c4762a1bSJed Brown 
3259566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dac, &cdac));
3269566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(daf, &cdaf));
327c4762a1bSJed Brown 
3289566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(dac, &coordsc));
3299566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(daf, &coordsf));
330c4762a1bSJed Brown 
3319566063dSJacob Faibussowitsch     PetscCall(DMCreateInterpolation(cdac, cdaf, &II, &scale));
3329566063dSJacob Faibussowitsch     PetscCall(MatInterpolate(II, coordsc, coordsf));
3339566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&II));
3349566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&scale));
335c4762a1bSJed Brown   }
336c4762a1bSJed Brown 
3379566063dSJacob Faibussowitsch   PetscCall(DMCreateInterpolation(dac, daf, &INTERP, NULL));
338c4762a1bSJed Brown 
3399566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dac, &ac));
3409566063dSJacob Faibussowitsch   PetscCall(VecSet(ac, 66.99));
341c4762a1bSJed Brown 
3429566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(daf, &af));
343c4762a1bSJed Brown 
3449566063dSJacob Faibussowitsch   PetscCall(MatMult(INTERP, ac, af));
345c4762a1bSJed Brown 
346c4762a1bSJed Brown   {
347c4762a1bSJed Brown     Vec       afexact;
348c4762a1bSJed Brown     PetscReal nrm;
349c4762a1bSJed Brown     PetscInt  N;
350c4762a1bSJed Brown 
3519566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(daf, &afexact));
3529566063dSJacob Faibussowitsch     PetscCall(VecSet(afexact, 66.99));
3539566063dSJacob Faibussowitsch     PetscCall(VecAXPY(afexact, -1.0, af)); /* af <= af - afinterp */
3549566063dSJacob Faibussowitsch     PetscCall(VecNorm(afexact, NORM_2, &nrm));
3559566063dSJacob Faibussowitsch     PetscCall(VecGetSize(afexact, &N));
35663a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "=>%" PetscInt_FMT ", interp err = %1.4e\n", mx, Mx, (double)(nrm / PetscSqrtReal((PetscReal)N))));
3579566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&afexact));
358c4762a1bSJed Brown   }
359c4762a1bSJed Brown 
3609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-output", &output, NULL));
361c4762a1bSJed Brown   if (output) {
3629566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "dac_1D.vtr", &vv));
3639566063dSJacob Faibussowitsch     PetscCall(VecView(ac, vv));
3649566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vv));
365c4762a1bSJed Brown 
3669566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "daf_1D.vtr", &vv));
3679566063dSJacob Faibussowitsch     PetscCall(VecView(af, vv));
3689566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vv));
369c4762a1bSJed Brown   }
370c4762a1bSJed Brown 
3719566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&INTERP));
3729566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dac));
3739566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&daf));
3749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ac));
3759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&af));
376c4762a1bSJed Brown   PetscFunctionReturn(0);
377c4762a1bSJed Brown }
378c4762a1bSJed Brown 
3799371c9d4SSatish Balay PetscErrorCode da_test_RefineCoords2D(PetscInt mx, PetscInt my) {
380c4762a1bSJed Brown   DM          dac, daf;
381c4762a1bSJed Brown   PetscViewer vv;
382c4762a1bSJed Brown   Vec         ac, af;
383c4762a1bSJed Brown   PetscInt    map_id, Mx, My;
384c4762a1bSJed Brown   Mat         II, INTERP;
385c4762a1bSJed Brown   Vec         scale;
386c4762a1bSJed Brown   PetscBool   output = PETSC_FALSE;
387c4762a1bSJed Brown 
388c4762a1bSJed Brown   PetscFunctionBeginUser;
3899566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, mx + 1, my + 1, PETSC_DECIDE, PETSC_DECIDE, 1, /* 1 dof */ 1, /* stencil = 1 */ NULL, NULL, &dac));
3909566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dac));
3919566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dac));
392c4762a1bSJed Brown 
3939566063dSJacob Faibussowitsch   PetscCall(DMRefine(dac, MPI_COMM_NULL, &daf));
3949566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf, 0, &Mx, &My, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
3959371c9d4SSatish Balay   Mx--;
3969371c9d4SSatish Balay   My--;
397c4762a1bSJed Brown 
3989566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(dac, -1.0, 1.0, -1.0, 1.0, PETSC_DECIDE, PETSC_DECIDE));
3999566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(daf, -1.0, 1.0, -1.0, 1.0, PETSC_DECIDE, PETSC_DECIDE));
400c4762a1bSJed Brown 
401c4762a1bSJed Brown   /* apply conformal mappings */
402c4762a1bSJed Brown   map_id = 0;
4039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-cmap", &map_id, NULL));
404*48a46eb9SPierre Jolivet   if (map_id >= 1) PetscCall(DAApplyConformalMapping(dac, map_id));
405c4762a1bSJed Brown 
406c4762a1bSJed Brown   {
407c4762a1bSJed Brown     DM  cdaf, cdac;
408c4762a1bSJed Brown     Vec coordsc, coordsf;
409c4762a1bSJed Brown 
4109566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dac, &cdac));
4119566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(daf, &cdaf));
412c4762a1bSJed Brown 
4139566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(dac, &coordsc));
4149566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(daf, &coordsf));
415c4762a1bSJed Brown 
4169566063dSJacob Faibussowitsch     PetscCall(DMCreateInterpolation(cdac, cdaf, &II, &scale));
4179566063dSJacob Faibussowitsch     PetscCall(MatInterpolate(II, coordsc, coordsf));
4189566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&II));
4199566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&scale));
420c4762a1bSJed Brown   }
421c4762a1bSJed Brown 
4229566063dSJacob Faibussowitsch   PetscCall(DMCreateInterpolation(dac, daf, &INTERP, NULL));
423c4762a1bSJed Brown 
4249566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dac, &ac));
4259566063dSJacob Faibussowitsch   PetscCall(DADefineXLinearField2D(dac, ac));
426c4762a1bSJed Brown 
4279566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(daf, &af));
4289566063dSJacob Faibussowitsch   PetscCall(MatMult(INTERP, ac, af));
429c4762a1bSJed Brown 
430c4762a1bSJed Brown   {
431c4762a1bSJed Brown     Vec       afexact;
432c4762a1bSJed Brown     PetscReal nrm;
433c4762a1bSJed Brown     PetscInt  N;
434c4762a1bSJed Brown 
4359566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(daf, &afexact));
4369566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(afexact));
4379566063dSJacob Faibussowitsch     PetscCall(DADefineXLinearField2D(daf, afexact));
4389566063dSJacob Faibussowitsch     PetscCall(VecAXPY(afexact, -1.0, af)); /* af <= af - afinterp */
4399566063dSJacob Faibussowitsch     PetscCall(VecNorm(afexact, NORM_2, &nrm));
4409566063dSJacob Faibussowitsch     PetscCall(VecGetSize(afexact, &N));
44163a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[%" PetscInt_FMT " x %" PetscInt_FMT "]=>[%" PetscInt_FMT " x %" PetscInt_FMT "], interp err = %1.4e\n", mx, my, Mx, My, (double)(nrm / PetscSqrtReal((PetscReal)N))));
4429566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&afexact));
443c4762a1bSJed Brown   }
444c4762a1bSJed Brown 
4459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-output", &output, NULL));
446c4762a1bSJed Brown   if (output) {
4479566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "dac_2D.vtr", &vv));
4489566063dSJacob Faibussowitsch     PetscCall(VecView(ac, vv));
4499566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vv));
450c4762a1bSJed Brown 
4519566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "daf_2D.vtr", &vv));
4529566063dSJacob Faibussowitsch     PetscCall(VecView(af, vv));
4539566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vv));
454c4762a1bSJed Brown   }
455c4762a1bSJed Brown 
4569566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&INTERP));
4579566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dac));
4589566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&daf));
4599566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ac));
4609566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&af));
461c4762a1bSJed Brown   PetscFunctionReturn(0);
462c4762a1bSJed Brown }
463c4762a1bSJed Brown 
4649371c9d4SSatish Balay PetscErrorCode da_test_RefineCoords3D(PetscInt mx, PetscInt my, PetscInt mz) {
465c4762a1bSJed Brown   DM          dac, daf;
466c4762a1bSJed Brown   PetscViewer vv;
467c4762a1bSJed Brown   Vec         ac, af;
468c4762a1bSJed Brown   PetscInt    map_id, Mx, My, Mz;
469c4762a1bSJed Brown   Mat         II, INTERP;
470c4762a1bSJed Brown   Vec         scale;
471c4762a1bSJed Brown   PetscBool   output = PETSC_FALSE;
472c4762a1bSJed Brown 
473c4762a1bSJed Brown   PetscFunctionBeginUser;
474d0609cedSBarry Smith   PetscCall(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, mx + 1, my + 1, mz + 1, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, 1, /* 1 dof */
475d0609cedSBarry Smith                          1, /* stencil = 1 */ NULL, NULL, NULL, &dac));
4769566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dac));
4779566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dac));
478c4762a1bSJed Brown 
4799566063dSJacob Faibussowitsch   PetscCall(DMRefine(dac, MPI_COMM_NULL, &daf));
4809566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(daf, 0, &Mx, &My, &Mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
4819371c9d4SSatish Balay   Mx--;
4829371c9d4SSatish Balay   My--;
4839371c9d4SSatish Balay   Mz--;
484c4762a1bSJed Brown 
4859566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(dac, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0));
4869566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(daf, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0));
487c4762a1bSJed Brown 
488c4762a1bSJed Brown   /* apply trilinear mappings */
4899566063dSJacob Faibussowitsch   /*PetscCall(DAApplyTrilinearMapping(dac));*/
490c4762a1bSJed Brown   /* apply conformal mappings */
491c4762a1bSJed Brown   map_id = 0;
4929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-cmap", &map_id, NULL));
493*48a46eb9SPierre Jolivet   if (map_id >= 1) PetscCall(DAApplyConformalMapping(dac, map_id));
494c4762a1bSJed Brown 
495c4762a1bSJed Brown   {
496c4762a1bSJed Brown     DM  cdaf, cdac;
497c4762a1bSJed Brown     Vec coordsc, coordsf;
498c4762a1bSJed Brown 
4999566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dac, &cdac));
5009566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(daf, &cdaf));
501c4762a1bSJed Brown 
5029566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(dac, &coordsc));
5039566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(daf, &coordsf));
504c4762a1bSJed Brown 
5059566063dSJacob Faibussowitsch     PetscCall(DMCreateInterpolation(cdac, cdaf, &II, &scale));
5069566063dSJacob Faibussowitsch     PetscCall(MatInterpolate(II, coordsc, coordsf));
5079566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&II));
5089566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&scale));
509c4762a1bSJed Brown   }
510c4762a1bSJed Brown 
5119566063dSJacob Faibussowitsch   PetscCall(DMCreateInterpolation(dac, daf, &INTERP, NULL));
512c4762a1bSJed Brown 
5139566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dac, &ac));
5149566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(ac));
5159566063dSJacob Faibussowitsch   PetscCall(DADefineXLinearField3D(dac, ac));
516c4762a1bSJed Brown 
5179566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(daf, &af));
5189566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(af));
519c4762a1bSJed Brown 
5209566063dSJacob Faibussowitsch   PetscCall(MatMult(INTERP, ac, af));
521c4762a1bSJed Brown 
522c4762a1bSJed Brown   {
523c4762a1bSJed Brown     Vec       afexact;
524c4762a1bSJed Brown     PetscReal nrm;
525c4762a1bSJed Brown     PetscInt  N;
526c4762a1bSJed Brown 
5279566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(daf, &afexact));
5289566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(afexact));
5299566063dSJacob Faibussowitsch     PetscCall(DADefineXLinearField3D(daf, afexact));
5309566063dSJacob Faibussowitsch     PetscCall(VecAXPY(afexact, -1.0, af)); /* af <= af - afinterp */
5319566063dSJacob Faibussowitsch     PetscCall(VecNorm(afexact, NORM_2, &nrm));
5329566063dSJacob Faibussowitsch     PetscCall(VecGetSize(afexact, &N));
53363a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[%" PetscInt_FMT " x %" PetscInt_FMT " x %" PetscInt_FMT "]=>[%" PetscInt_FMT " x %" PetscInt_FMT " x %" PetscInt_FMT "], interp err = %1.4e\n", mx, my, mz, Mx, My, Mz, (double)(nrm / PetscSqrtReal((PetscReal)N))));
5349566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&afexact));
535c4762a1bSJed Brown   }
536c4762a1bSJed Brown 
5379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-output", &output, NULL));
538c4762a1bSJed Brown   if (output) {
5399566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "dac_3D.vtr", &vv));
5409566063dSJacob Faibussowitsch     PetscCall(VecView(ac, vv));
5419566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vv));
542c4762a1bSJed Brown 
5439566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "daf_3D.vtr", &vv));
5449566063dSJacob Faibussowitsch     PetscCall(VecView(af, vv));
5459566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&vv));
546c4762a1bSJed Brown   }
547c4762a1bSJed Brown 
5489566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&INTERP));
5499566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dac));
5509566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&daf));
5519566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ac));
5529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&af));
553c4762a1bSJed Brown   PetscFunctionReturn(0);
554c4762a1bSJed Brown }
555c4762a1bSJed Brown 
5569371c9d4SSatish Balay int main(int argc, char **argv) {
557c4762a1bSJed Brown   PetscInt mx = 2, my = 2, mz = 2, l, nl, dim;
558c4762a1bSJed Brown 
559327415f7SBarry Smith   PetscFunctionBeginUser;
5609566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, 0, help));
5619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &mx, 0));
5629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &my, 0));
5639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mz", &mz, 0));
564c4762a1bSJed Brown   nl = 1;
5659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nl", &nl, 0));
566c4762a1bSJed Brown   dim = 2;
5679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, 0));
568c4762a1bSJed Brown 
569c4762a1bSJed Brown   for (l = 0; l < nl; l++) {
5701baa6e33SBarry Smith     if (dim == 1) PetscCall(da_test_RefineCoords1D(mx));
5711baa6e33SBarry Smith     else if (dim == 2) PetscCall(da_test_RefineCoords2D(mx, my));
5721baa6e33SBarry Smith     else if (dim == 3) PetscCall(da_test_RefineCoords3D(mx, my, mz));
573c4762a1bSJed Brown     mx = mx * 2;
574c4762a1bSJed Brown     my = my * 2;
575c4762a1bSJed Brown     mz = mz * 2;
576c4762a1bSJed Brown   }
5779566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
578b122ec5aSJacob Faibussowitsch   return 0;
579c4762a1bSJed Brown }
580c4762a1bSJed Brown 
581c4762a1bSJed Brown /*TEST
582c4762a1bSJed Brown 
583c4762a1bSJed Brown    test:
584c4762a1bSJed Brown       suffix: 1d
585c4762a1bSJed Brown       args: -mx 10 -nl 6 -dim 1
586c4762a1bSJed Brown 
587c4762a1bSJed Brown    test:
588c4762a1bSJed Brown       suffix: 2d
58963037964SStefano Zampini       output_file: output/ex36_2d.out
59063037964SStefano Zampini       args: -mx 10 -my 10 -nl 6 -dim 2 -cmap {{0 1 2 3}}
591c4762a1bSJed Brown 
592c4762a1bSJed Brown    test:
593c4762a1bSJed Brown       suffix: 2dp1
594c4762a1bSJed Brown       nsize: 8
595c4762a1bSJed Brown       args: -mx 10 -my 10 -nl 4 -dim 2 -cmap 3 -da_refine_x 3 -da_refine_y 4
596c4762a1bSJed Brown       timeoutfactor: 2
597c4762a1bSJed Brown 
598c4762a1bSJed Brown    test:
599c4762a1bSJed Brown       suffix: 2dp2
600c4762a1bSJed Brown       nsize: 8
601c4762a1bSJed Brown       args: -mx 10 -my 10 -nl 4 -dim 2 -cmap 3 -da_refine_x 3 -da_refine_y 1
602c4762a1bSJed Brown       timeoutfactor: 2
603c4762a1bSJed Brown 
604c4762a1bSJed Brown    test:
605c4762a1bSJed Brown       suffix: 3d
606c4762a1bSJed Brown       args: -mx 5 -my 5 -mz 5 -nl 4 -dim 3 -cmap 3
607c4762a1bSJed Brown 
608c4762a1bSJed Brown    test:
609c4762a1bSJed Brown       suffix: 3dp1
610c4762a1bSJed Brown       nsize: 32
611c4762a1bSJed Brown       args: -mx 5 -my 5 -mz 5 -nl 3 -dim 3 -cmap 1 -da_refine_x 1 -da_refine_y 3 -da_refine_z 4
612c4762a1bSJed Brown 
613c4762a1bSJed Brown TEST*/
614