1c4762a1bSJed Brown static char help[] = "Checks the functionality of DMGetInterpolation() on deformed grids.\n\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown #include <petscdm.h>
4c4762a1bSJed Brown #include <petscdmda.h>
5c4762a1bSJed Brown
6c4762a1bSJed Brown typedef struct _n_CCmplx CCmplx;
7c4762a1bSJed Brown struct _n_CCmplx {
8c4762a1bSJed Brown PetscReal real;
9c4762a1bSJed Brown PetscReal imag;
10c4762a1bSJed Brown };
11c4762a1bSJed Brown
CCmplxPow(CCmplx a,PetscReal n)12d71ae5a4SJacob Faibussowitsch CCmplx CCmplxPow(CCmplx a, PetscReal n)
13d71ae5a4SJacob Faibussowitsch {
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 }
CCmplxExp(CCmplx a)22d71ae5a4SJacob Faibussowitsch CCmplx CCmplxExp(CCmplx a)
23d71ae5a4SJacob Faibussowitsch {
24c4762a1bSJed Brown CCmplx b;
25c4762a1bSJed Brown b.real = PetscExpReal(a.real) * PetscCosReal(a.imag);
26c4762a1bSJed Brown b.imag = PetscExpReal(a.real) * PetscSinReal(a.imag);
27c4762a1bSJed Brown return b;
28c4762a1bSJed Brown }
CCmplxSqrt(CCmplx a)29d71ae5a4SJacob Faibussowitsch CCmplx CCmplxSqrt(CCmplx a)
30d71ae5a4SJacob Faibussowitsch {
31c4762a1bSJed Brown CCmplx b;
32c4762a1bSJed Brown PetscReal r, theta;
33c4762a1bSJed Brown r = PetscSqrtReal(a.real * a.real + a.imag * a.imag);
34c4762a1bSJed Brown theta = PetscAtan2Real(a.imag, a.real);
35c4762a1bSJed Brown b.real = PetscSqrtReal(r) * PetscCosReal(0.5 * theta);
36c4762a1bSJed Brown b.imag = PetscSqrtReal(r) * PetscSinReal(0.5 * theta);
37c4762a1bSJed Brown return b;
38c4762a1bSJed Brown }
CCmplxAdd(CCmplx a,CCmplx c)39d71ae5a4SJacob Faibussowitsch CCmplx CCmplxAdd(CCmplx a, CCmplx c)
40d71ae5a4SJacob Faibussowitsch {
41c4762a1bSJed Brown CCmplx b;
42c4762a1bSJed Brown b.real = a.real + c.real;
43c4762a1bSJed Brown b.imag = a.imag + c.imag;
44c4762a1bSJed Brown return b;
45c4762a1bSJed Brown }
CCmplxRe(CCmplx a)46d71ae5a4SJacob Faibussowitsch PetscScalar CCmplxRe(CCmplx a)
47d71ae5a4SJacob Faibussowitsch {
48*300f1712SStefano Zampini return a.real;
49c4762a1bSJed Brown }
CCmplxIm(CCmplx a)50d71ae5a4SJacob Faibussowitsch PetscScalar CCmplxIm(CCmplx a)
51d71ae5a4SJacob Faibussowitsch {
52*300f1712SStefano Zampini return a.imag;
53c4762a1bSJed Brown }
54c4762a1bSJed Brown
DAApplyConformalMapping(DM da,PetscInt idx)55d71ae5a4SJacob Faibussowitsch PetscErrorCode DAApplyConformalMapping(DM da, PetscInt idx)
56d71ae5a4SJacob Faibussowitsch {
57c4762a1bSJed Brown PetscInt i, n;
58c4762a1bSJed Brown PetscInt sx, nx, sy, ny, sz, nz, dim;
59c4762a1bSJed Brown Vec Gcoords;
60c4762a1bSJed Brown PetscScalar *XX;
61c4762a1bSJed Brown PetscScalar xx, yy, zz;
62c4762a1bSJed Brown DM cda;
63c4762a1bSJed Brown
64c4762a1bSJed Brown PetscFunctionBeginUser;
65c4762a1bSJed Brown if (idx == 0) {
663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
67c4762a1bSJed Brown } else if (idx == 1) { /* dam break */
689566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0));
69c4762a1bSJed Brown } else if (idx == 2) { /* stagnation in a corner */
709566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, -1.0, 1.0, 0.0, 1.0, -1.0, 1.0));
71c4762a1bSJed Brown } else if (idx == 3) { /* nautilis */
729566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0));
731baa6e33SBarry Smith } else if (idx == 4) PetscCall(DMDASetUniformCoordinates(da, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0));
74c4762a1bSJed Brown
759566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &cda));
769566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &Gcoords));
77c4762a1bSJed Brown
789566063dSJacob Faibussowitsch PetscCall(VecGetArray(Gcoords, &XX));
799566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &sx, &sy, &sz, &nx, &ny, &nz));
809566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
819566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(Gcoords, &n));
82c4762a1bSJed Brown n = n / dim;
83c4762a1bSJed Brown
84c4762a1bSJed Brown for (i = 0; i < n; i++) {
85c4762a1bSJed Brown if ((dim == 3) && (idx != 2)) {
86c4762a1bSJed Brown PetscScalar Ni[8];
87c4762a1bSJed Brown PetscScalar xi = XX[dim * i];
88c4762a1bSJed Brown PetscScalar eta = XX[dim * i + 1];
89c4762a1bSJed Brown PetscScalar zeta = XX[dim * i + 2];
90c4762a1bSJed Brown PetscScalar xn[] = {-1.0, 1.0, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0};
91c4762a1bSJed Brown PetscScalar yn[] = {-1.0, -1.0, 1.0, 1.0, -1.0, -1.0, 1.0, 1.0};
92c4762a1bSJed Brown PetscScalar zn[] = {-0.1, -4.0, -0.2, -1.0, 0.1, 4.0, 0.2, 1.0};
93c4762a1bSJed Brown PetscInt p;
94c4762a1bSJed Brown
95c4762a1bSJed Brown Ni[0] = 0.125 * (1.0 - xi) * (1.0 - eta) * (1.0 - zeta);
96c4762a1bSJed Brown Ni[1] = 0.125 * (1.0 + xi) * (1.0 - eta) * (1.0 - zeta);
97c4762a1bSJed Brown Ni[2] = 0.125 * (1.0 - xi) * (1.0 + eta) * (1.0 - zeta);
98c4762a1bSJed Brown Ni[3] = 0.125 * (1.0 + xi) * (1.0 + eta) * (1.0 - zeta);
99c4762a1bSJed Brown
100c4762a1bSJed Brown Ni[4] = 0.125 * (1.0 - xi) * (1.0 - eta) * (1.0 + zeta);
101c4762a1bSJed Brown Ni[5] = 0.125 * (1.0 + xi) * (1.0 - eta) * (1.0 + zeta);
102c4762a1bSJed Brown Ni[6] = 0.125 * (1.0 - xi) * (1.0 + eta) * (1.0 + zeta);
103c4762a1bSJed Brown Ni[7] = 0.125 * (1.0 + xi) * (1.0 + eta) * (1.0 + zeta);
104c4762a1bSJed Brown
105c4762a1bSJed Brown xx = yy = zz = 0.0;
106c4762a1bSJed Brown for (p = 0; p < 8; p++) {
107c4762a1bSJed Brown xx += Ni[p] * xn[p];
108c4762a1bSJed Brown yy += Ni[p] * yn[p];
109c4762a1bSJed Brown zz += Ni[p] * zn[p];
110c4762a1bSJed Brown }
111c4762a1bSJed Brown XX[dim * i] = xx;
112c4762a1bSJed Brown XX[dim * i + 1] = yy;
113c4762a1bSJed Brown XX[dim * i + 2] = zz;
114c4762a1bSJed Brown }
115c4762a1bSJed Brown
116c4762a1bSJed Brown if (idx == 1) {
117c4762a1bSJed Brown CCmplx zeta, t1, t2;
118c4762a1bSJed Brown
119c4762a1bSJed Brown xx = XX[dim * i] - 0.8;
120c4762a1bSJed Brown yy = XX[dim * i + 1] + 1.5;
121c4762a1bSJed Brown
122c4762a1bSJed Brown zeta.real = PetscRealPart(xx);
123c4762a1bSJed Brown zeta.imag = PetscRealPart(yy);
124c4762a1bSJed Brown
125c4762a1bSJed Brown t1 = CCmplxPow(zeta, -1.0);
126c4762a1bSJed Brown t2 = CCmplxAdd(zeta, t1);
127c4762a1bSJed Brown
128c4762a1bSJed Brown XX[dim * i] = CCmplxRe(t2);
129c4762a1bSJed Brown XX[dim * i + 1] = CCmplxIm(t2);
130c4762a1bSJed Brown } else if (idx == 2) {
131c4762a1bSJed Brown CCmplx zeta, t1;
132c4762a1bSJed Brown
133c4762a1bSJed Brown xx = XX[dim * i];
134c4762a1bSJed Brown yy = XX[dim * i + 1];
135c4762a1bSJed Brown zeta.real = PetscRealPart(xx);
136c4762a1bSJed Brown zeta.imag = PetscRealPart(yy);
137c4762a1bSJed Brown
138c4762a1bSJed Brown t1 = CCmplxSqrt(zeta);
139c4762a1bSJed Brown XX[dim * i] = CCmplxRe(t1);
140c4762a1bSJed Brown XX[dim * i + 1] = CCmplxIm(t1);
141c4762a1bSJed Brown } else if (idx == 3) {
142c4762a1bSJed Brown CCmplx zeta, t1, t2;
143c4762a1bSJed Brown
144c4762a1bSJed Brown xx = XX[dim * i] - 0.8;
145c4762a1bSJed Brown yy = XX[dim * i + 1] + 1.5;
146c4762a1bSJed Brown
147c4762a1bSJed Brown zeta.real = PetscRealPart(xx);
148c4762a1bSJed Brown zeta.imag = PetscRealPart(yy);
149c4762a1bSJed Brown t1 = CCmplxPow(zeta, -1.0);
150c4762a1bSJed Brown t2 = CCmplxAdd(zeta, t1);
151c4762a1bSJed Brown XX[dim * i] = CCmplxRe(t2);
152c4762a1bSJed Brown XX[dim * i + 1] = CCmplxIm(t2);
153c4762a1bSJed Brown
154c4762a1bSJed Brown xx = XX[dim * i];
155c4762a1bSJed Brown yy = XX[dim * i + 1];
156c4762a1bSJed Brown zeta.real = PetscRealPart(xx);
157c4762a1bSJed Brown zeta.imag = PetscRealPart(yy);
158c4762a1bSJed Brown t1 = CCmplxExp(zeta);
159c4762a1bSJed Brown XX[dim * i] = CCmplxRe(t1);
160c4762a1bSJed Brown XX[dim * i + 1] = CCmplxIm(t1);
161c4762a1bSJed Brown
162c4762a1bSJed Brown xx = XX[dim * i] + 0.4;
163c4762a1bSJed Brown yy = XX[dim * i + 1];
164c4762a1bSJed Brown zeta.real = PetscRealPart(xx);
165c4762a1bSJed Brown zeta.imag = PetscRealPart(yy);
166c4762a1bSJed Brown t1 = CCmplxPow(zeta, 2.0);
167c4762a1bSJed Brown XX[dim * i] = CCmplxRe(t1);
168c4762a1bSJed Brown XX[dim * i + 1] = CCmplxIm(t1);
169c4762a1bSJed Brown } else if (idx == 4) {
170c4762a1bSJed Brown PetscScalar Ni[4];
171c4762a1bSJed Brown PetscScalar xi = XX[dim * i];
172c4762a1bSJed Brown PetscScalar eta = XX[dim * i + 1];
173c4762a1bSJed Brown PetscScalar xn[] = {0.0, 2.0, 0.2, 3.5};
174c4762a1bSJed Brown PetscScalar yn[] = {-1.3, 0.0, 2.0, 4.0};
175c4762a1bSJed Brown PetscInt p;
176c4762a1bSJed Brown
177c4762a1bSJed Brown Ni[0] = 0.25 * (1.0 - xi) * (1.0 - eta);
178c4762a1bSJed Brown Ni[1] = 0.25 * (1.0 + xi) * (1.0 - eta);
179c4762a1bSJed Brown Ni[2] = 0.25 * (1.0 - xi) * (1.0 + eta);
180c4762a1bSJed Brown Ni[3] = 0.25 * (1.0 + xi) * (1.0 + eta);
181c4762a1bSJed Brown
182c4762a1bSJed Brown xx = yy = 0.0;
183c4762a1bSJed Brown for (p = 0; p < 4; p++) {
184c4762a1bSJed Brown xx += Ni[p] * xn[p];
185c4762a1bSJed Brown yy += Ni[p] * yn[p];
186c4762a1bSJed Brown }
187c4762a1bSJed Brown XX[dim * i] = xx;
188c4762a1bSJed Brown XX[dim * i + 1] = yy;
189c4762a1bSJed Brown }
190c4762a1bSJed Brown }
1919566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Gcoords, &XX));
1923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
193c4762a1bSJed Brown }
194c4762a1bSJed Brown
DAApplyTrilinearMapping(DM da)195d71ae5a4SJacob Faibussowitsch PetscErrorCode DAApplyTrilinearMapping(DM da)
196d71ae5a4SJacob Faibussowitsch {
197c4762a1bSJed Brown PetscInt i, j, k;
198c4762a1bSJed Brown PetscInt sx, nx, sy, ny, sz, nz;
199c4762a1bSJed Brown Vec Gcoords;
200c4762a1bSJed Brown DMDACoor3d ***XX;
201c4762a1bSJed Brown PetscScalar xx, yy, zz;
202c4762a1bSJed Brown DM cda;
203c4762a1bSJed Brown
204c4762a1bSJed Brown PetscFunctionBeginUser;
2059566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0));
2069566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &cda));
2079566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &Gcoords));
208c4762a1bSJed Brown
2099566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda, Gcoords, &XX));
2109566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &sx, &sy, &sz, &nx, &ny, &nz));
211c4762a1bSJed Brown
212c4762a1bSJed Brown for (i = sx; i < sx + nx; i++) {
213c4762a1bSJed Brown for (j = sy; j < sy + ny; j++) {
214c4762a1bSJed Brown for (k = sz; k < sz + nz; k++) {
215c4762a1bSJed Brown PetscScalar Ni[8];
216c4762a1bSJed Brown PetscScalar xi = XX[k][j][i].x;
217c4762a1bSJed Brown PetscScalar eta = XX[k][j][i].y;
218c4762a1bSJed Brown PetscScalar zeta = XX[k][j][i].z;
219c4762a1bSJed Brown PetscScalar xn[] = {0.0, 2.0, 0.2, 3.5, 0.0, 2.1, 0.23, 3.125};
220c4762a1bSJed Brown PetscScalar yn[] = {-1.3, 0.0, 2.0, 4.0, -1.45, -0.1, 2.24, 3.79};
221c4762a1bSJed Brown PetscScalar zn[] = {0.0, 0.3, -0.1, 0.123, 0.956, 1.32, 1.12, 0.798};
222c4762a1bSJed Brown PetscInt p;
223c4762a1bSJed Brown
224c4762a1bSJed Brown Ni[0] = 0.125 * (1.0 - xi) * (1.0 - eta) * (1.0 - zeta);
225c4762a1bSJed Brown Ni[1] = 0.125 * (1.0 + xi) * (1.0 - eta) * (1.0 - zeta);
226c4762a1bSJed Brown Ni[2] = 0.125 * (1.0 - xi) * (1.0 + eta) * (1.0 - zeta);
227c4762a1bSJed Brown Ni[3] = 0.125 * (1.0 + xi) * (1.0 + eta) * (1.0 - zeta);
228c4762a1bSJed Brown
229c4762a1bSJed Brown Ni[4] = 0.125 * (1.0 - xi) * (1.0 - eta) * (1.0 + zeta);
230c4762a1bSJed Brown Ni[5] = 0.125 * (1.0 + xi) * (1.0 - eta) * (1.0 + zeta);
231c4762a1bSJed Brown Ni[6] = 0.125 * (1.0 - xi) * (1.0 + eta) * (1.0 + zeta);
232c4762a1bSJed Brown Ni[7] = 0.125 * (1.0 + xi) * (1.0 + eta) * (1.0 + zeta);
233c4762a1bSJed Brown
234c4762a1bSJed Brown xx = yy = zz = 0.0;
235c4762a1bSJed Brown for (p = 0; p < 8; p++) {
236c4762a1bSJed Brown xx += Ni[p] * xn[p];
237c4762a1bSJed Brown yy += Ni[p] * yn[p];
238c4762a1bSJed Brown zz += Ni[p] * zn[p];
239c4762a1bSJed Brown }
240c4762a1bSJed Brown XX[k][j][i].x = xx;
241c4762a1bSJed Brown XX[k][j][i].y = yy;
242c4762a1bSJed Brown XX[k][j][i].z = zz;
243c4762a1bSJed Brown }
244c4762a1bSJed Brown }
245c4762a1bSJed Brown }
2469566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda, Gcoords, &XX));
2473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
248c4762a1bSJed Brown }
249c4762a1bSJed Brown
DADefineXLinearField2D(DM da,Vec field)250d71ae5a4SJacob Faibussowitsch PetscErrorCode DADefineXLinearField2D(DM da, Vec field)
251d71ae5a4SJacob Faibussowitsch {
252c4762a1bSJed Brown PetscInt i, j;
253c4762a1bSJed Brown PetscInt sx, nx, sy, ny;
254c4762a1bSJed Brown Vec Gcoords;
255c4762a1bSJed Brown DMDACoor2d **XX;
256c4762a1bSJed Brown PetscScalar **FF;
257c4762a1bSJed Brown DM cda;
258c4762a1bSJed Brown
259c4762a1bSJed Brown PetscFunctionBeginUser;
2609566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &cda));
2619566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &Gcoords));
262c4762a1bSJed Brown
2639566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda, Gcoords, &XX));
2649566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, field, &FF));
265c4762a1bSJed Brown
2669566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &sx, &sy, 0, &nx, &ny, 0));
267c4762a1bSJed Brown
268c4762a1bSJed Brown for (i = sx; i < sx + nx; i++) {
269ad540459SPierre Jolivet 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;
270c4762a1bSJed Brown }
271c4762a1bSJed Brown
2729566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, field, &FF));
2739566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda, Gcoords, &XX));
2743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
275c4762a1bSJed Brown }
276c4762a1bSJed Brown
DADefineXLinearField3D(DM da,Vec field)277d71ae5a4SJacob Faibussowitsch PetscErrorCode DADefineXLinearField3D(DM da, Vec field)
278d71ae5a4SJacob Faibussowitsch {
279c4762a1bSJed Brown PetscInt i, j, k;
280c4762a1bSJed Brown PetscInt sx, nx, sy, ny, sz, nz;
281c4762a1bSJed Brown Vec Gcoords;
282c4762a1bSJed Brown DMDACoor3d ***XX;
283c4762a1bSJed Brown PetscScalar ***FF;
284c4762a1bSJed Brown DM cda;
285c4762a1bSJed Brown
286c4762a1bSJed Brown PetscFunctionBeginUser;
2879566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &cda));
2889566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &Gcoords));
289c4762a1bSJed Brown
2909566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda, Gcoords, &XX));
2919566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, field, &FF));
292c4762a1bSJed Brown
2939566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &sx, &sy, &sz, &nx, &ny, &nz));
294c4762a1bSJed Brown
295c4762a1bSJed Brown for (k = sz; k < sz + nz; k++) {
296c4762a1bSJed Brown for (j = sy; j < sy + ny; j++) {
297c4762a1bSJed Brown for (i = sx; i < sx + nx; i++) {
2989371c9d4SSatish 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 +
2999371c9d4SSatish Balay 3.79 * XX[k][j][i].x * XX[k][j][i].y * XX[k][j][i].z;
300c4762a1bSJed Brown }
301c4762a1bSJed Brown }
302c4762a1bSJed Brown }
303c4762a1bSJed Brown
3049566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, field, &FF));
3059566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda, Gcoords, &XX));
3063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
307c4762a1bSJed Brown }
308c4762a1bSJed Brown
da_test_RefineCoords1D(PetscInt mx)309d71ae5a4SJacob Faibussowitsch PetscErrorCode da_test_RefineCoords1D(PetscInt mx)
310d71ae5a4SJacob Faibussowitsch {
311c4762a1bSJed Brown DM dac, daf;
312c4762a1bSJed Brown PetscViewer vv;
313c4762a1bSJed Brown Vec ac, af;
314c4762a1bSJed Brown PetscInt Mx;
315c4762a1bSJed Brown Mat II, INTERP;
316c4762a1bSJed Brown Vec scale;
317c4762a1bSJed Brown PetscBool output = PETSC_FALSE;
318c4762a1bSJed Brown
319c4762a1bSJed Brown PetscFunctionBeginUser;
3209566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, mx + 1, 1, /* 1 dof */ 1, /* stencil = 1 */ NULL, &dac));
3219566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dac));
3229566063dSJacob Faibussowitsch PetscCall(DMSetUp(dac));
323c4762a1bSJed Brown
3249566063dSJacob Faibussowitsch PetscCall(DMRefine(dac, MPI_COMM_NULL, &daf));
3259566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, 0, &Mx, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
326c4762a1bSJed Brown Mx--;
327c4762a1bSJed Brown
3289566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(dac, -1.0, 1.0, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE));
3299566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(daf, -1.0, 1.0, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE));
330c4762a1bSJed Brown
331c4762a1bSJed Brown {
332c4762a1bSJed Brown DM cdaf, cdac;
333c4762a1bSJed Brown Vec coordsc, coordsf;
334c4762a1bSJed Brown
3359566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dac, &cdac));
3369566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(daf, &cdaf));
337c4762a1bSJed Brown
3389566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dac, &coordsc));
3399566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(daf, &coordsf));
340c4762a1bSJed Brown
3419566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(cdac, cdaf, &II, &scale));
3429566063dSJacob Faibussowitsch PetscCall(MatInterpolate(II, coordsc, coordsf));
3439566063dSJacob Faibussowitsch PetscCall(MatDestroy(&II));
3449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&scale));
345c4762a1bSJed Brown }
346c4762a1bSJed Brown
3479566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dac, daf, &INTERP, NULL));
348c4762a1bSJed Brown
3499566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dac, &ac));
3509566063dSJacob Faibussowitsch PetscCall(VecSet(ac, 66.99));
351c4762a1bSJed Brown
3529566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(daf, &af));
353c4762a1bSJed Brown
3549566063dSJacob Faibussowitsch PetscCall(MatMult(INTERP, ac, af));
355c4762a1bSJed Brown
356c4762a1bSJed Brown {
357c4762a1bSJed Brown Vec afexact;
358c4762a1bSJed Brown PetscReal nrm;
359c4762a1bSJed Brown PetscInt N;
360c4762a1bSJed Brown
3619566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(daf, &afexact));
3629566063dSJacob Faibussowitsch PetscCall(VecSet(afexact, 66.99));
3639566063dSJacob Faibussowitsch PetscCall(VecAXPY(afexact, -1.0, af)); /* af <= af - afinterp */
3649566063dSJacob Faibussowitsch PetscCall(VecNorm(afexact, NORM_2, &nrm));
3659566063dSJacob Faibussowitsch PetscCall(VecGetSize(afexact, &N));
36663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "=>%" PetscInt_FMT ", interp err = %1.4e\n", mx, Mx, (double)(nrm / PetscSqrtReal((PetscReal)N))));
3679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&afexact));
368c4762a1bSJed Brown }
369c4762a1bSJed Brown
3709566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-output", &output, NULL));
371c4762a1bSJed Brown if (output) {
3729566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "dac_1D.vtr", &vv));
3739566063dSJacob Faibussowitsch PetscCall(VecView(ac, vv));
3749566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vv));
375c4762a1bSJed Brown
3769566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "daf_1D.vtr", &vv));
3779566063dSJacob Faibussowitsch PetscCall(VecView(af, vv));
3789566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vv));
379c4762a1bSJed Brown }
380c4762a1bSJed Brown
3819566063dSJacob Faibussowitsch PetscCall(MatDestroy(&INTERP));
3829566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dac));
3839566063dSJacob Faibussowitsch PetscCall(DMDestroy(&daf));
3849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ac));
3859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&af));
3863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
387c4762a1bSJed Brown }
388c4762a1bSJed Brown
da_test_RefineCoords2D(PetscInt mx,PetscInt my)389d71ae5a4SJacob Faibussowitsch PetscErrorCode da_test_RefineCoords2D(PetscInt mx, PetscInt my)
390d71ae5a4SJacob Faibussowitsch {
391c4762a1bSJed Brown DM dac, daf;
392c4762a1bSJed Brown PetscViewer vv;
393c4762a1bSJed Brown Vec ac, af;
394c4762a1bSJed Brown PetscInt map_id, Mx, My;
395c4762a1bSJed Brown Mat II, INTERP;
396c4762a1bSJed Brown Vec scale;
397c4762a1bSJed Brown PetscBool output = PETSC_FALSE;
398c4762a1bSJed Brown
399c4762a1bSJed Brown PetscFunctionBeginUser;
4009566063dSJacob 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));
4019566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dac));
4029566063dSJacob Faibussowitsch PetscCall(DMSetUp(dac));
403c4762a1bSJed Brown
4049566063dSJacob Faibussowitsch PetscCall(DMRefine(dac, MPI_COMM_NULL, &daf));
4059566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, 0, &Mx, &My, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
4069371c9d4SSatish Balay Mx--;
4079371c9d4SSatish Balay My--;
408c4762a1bSJed Brown
4099566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(dac, -1.0, 1.0, -1.0, 1.0, PETSC_DECIDE, PETSC_DECIDE));
4109566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(daf, -1.0, 1.0, -1.0, 1.0, PETSC_DECIDE, PETSC_DECIDE));
411c4762a1bSJed Brown
412c4762a1bSJed Brown /* apply conformal mappings */
413c4762a1bSJed Brown map_id = 0;
4149566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-cmap", &map_id, NULL));
41548a46eb9SPierre Jolivet if (map_id >= 1) PetscCall(DAApplyConformalMapping(dac, map_id));
416c4762a1bSJed Brown
417c4762a1bSJed Brown {
418c4762a1bSJed Brown DM cdaf, cdac;
419c4762a1bSJed Brown Vec coordsc, coordsf;
420c4762a1bSJed Brown
4219566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dac, &cdac));
4229566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(daf, &cdaf));
423c4762a1bSJed Brown
4249566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dac, &coordsc));
4259566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(daf, &coordsf));
426c4762a1bSJed Brown
4279566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(cdac, cdaf, &II, &scale));
4289566063dSJacob Faibussowitsch PetscCall(MatInterpolate(II, coordsc, coordsf));
4299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&II));
4309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&scale));
431c4762a1bSJed Brown }
432c4762a1bSJed Brown
4339566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dac, daf, &INTERP, NULL));
434c4762a1bSJed Brown
4359566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dac, &ac));
4369566063dSJacob Faibussowitsch PetscCall(DADefineXLinearField2D(dac, ac));
437c4762a1bSJed Brown
4389566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(daf, &af));
4399566063dSJacob Faibussowitsch PetscCall(MatMult(INTERP, ac, af));
440c4762a1bSJed Brown
441c4762a1bSJed Brown {
442c4762a1bSJed Brown Vec afexact;
443c4762a1bSJed Brown PetscReal nrm;
444c4762a1bSJed Brown PetscInt N;
445c4762a1bSJed Brown
4469566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(daf, &afexact));
4479566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(afexact));
4489566063dSJacob Faibussowitsch PetscCall(DADefineXLinearField2D(daf, afexact));
4499566063dSJacob Faibussowitsch PetscCall(VecAXPY(afexact, -1.0, af)); /* af <= af - afinterp */
4509566063dSJacob Faibussowitsch PetscCall(VecNorm(afexact, NORM_2, &nrm));
4519566063dSJacob Faibussowitsch PetscCall(VecGetSize(afexact, &N));
45263a3b9bcSJacob 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))));
4539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&afexact));
454c4762a1bSJed Brown }
455c4762a1bSJed Brown
4569566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-output", &output, NULL));
457c4762a1bSJed Brown if (output) {
4589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "dac_2D.vtr", &vv));
4599566063dSJacob Faibussowitsch PetscCall(VecView(ac, vv));
4609566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vv));
461c4762a1bSJed Brown
4629566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "daf_2D.vtr", &vv));
4639566063dSJacob Faibussowitsch PetscCall(VecView(af, vv));
4649566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vv));
465c4762a1bSJed Brown }
466c4762a1bSJed Brown
4679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&INTERP));
4689566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dac));
4699566063dSJacob Faibussowitsch PetscCall(DMDestroy(&daf));
4709566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ac));
4719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&af));
4723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
473c4762a1bSJed Brown }
474c4762a1bSJed Brown
da_test_RefineCoords3D(PetscInt mx,PetscInt my,PetscInt mz)475d71ae5a4SJacob Faibussowitsch PetscErrorCode da_test_RefineCoords3D(PetscInt mx, PetscInt my, PetscInt mz)
476d71ae5a4SJacob Faibussowitsch {
477c4762a1bSJed Brown DM dac, daf;
478c4762a1bSJed Brown PetscViewer vv;
479c4762a1bSJed Brown Vec ac, af;
480c4762a1bSJed Brown PetscInt map_id, Mx, My, Mz;
481c4762a1bSJed Brown Mat II, INTERP;
482c4762a1bSJed Brown Vec scale;
483c4762a1bSJed Brown PetscBool output = PETSC_FALSE;
484c4762a1bSJed Brown
485c4762a1bSJed Brown PetscFunctionBeginUser;
486d0609cedSBarry 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 */
487d0609cedSBarry Smith 1, /* stencil = 1 */ NULL, NULL, NULL, &dac));
4889566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dac));
4899566063dSJacob Faibussowitsch PetscCall(DMSetUp(dac));
490c4762a1bSJed Brown
4919566063dSJacob Faibussowitsch PetscCall(DMRefine(dac, MPI_COMM_NULL, &daf));
4929566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, 0, &Mx, &My, &Mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
4939371c9d4SSatish Balay Mx--;
4949371c9d4SSatish Balay My--;
4959371c9d4SSatish Balay Mz--;
496c4762a1bSJed Brown
4979566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(dac, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0));
4989566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(daf, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0));
499c4762a1bSJed Brown
500c4762a1bSJed Brown /* apply trilinear mappings */
5019566063dSJacob Faibussowitsch /*PetscCall(DAApplyTrilinearMapping(dac));*/
502c4762a1bSJed Brown /* apply conformal mappings */
503c4762a1bSJed Brown map_id = 0;
5049566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-cmap", &map_id, NULL));
50548a46eb9SPierre Jolivet if (map_id >= 1) PetscCall(DAApplyConformalMapping(dac, map_id));
506c4762a1bSJed Brown
507c4762a1bSJed Brown {
508c4762a1bSJed Brown DM cdaf, cdac;
509c4762a1bSJed Brown Vec coordsc, coordsf;
510c4762a1bSJed Brown
5119566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dac, &cdac));
5129566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(daf, &cdaf));
513c4762a1bSJed Brown
5149566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dac, &coordsc));
5159566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(daf, &coordsf));
516c4762a1bSJed Brown
5179566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(cdac, cdaf, &II, &scale));
5189566063dSJacob Faibussowitsch PetscCall(MatInterpolate(II, coordsc, coordsf));
5199566063dSJacob Faibussowitsch PetscCall(MatDestroy(&II));
5209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&scale));
521c4762a1bSJed Brown }
522c4762a1bSJed Brown
5239566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dac, daf, &INTERP, NULL));
524c4762a1bSJed Brown
5259566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dac, &ac));
5269566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(ac));
5279566063dSJacob Faibussowitsch PetscCall(DADefineXLinearField3D(dac, ac));
528c4762a1bSJed Brown
5299566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(daf, &af));
5309566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(af));
531c4762a1bSJed Brown
5329566063dSJacob Faibussowitsch PetscCall(MatMult(INTERP, ac, af));
533c4762a1bSJed Brown
534c4762a1bSJed Brown {
535c4762a1bSJed Brown Vec afexact;
536c4762a1bSJed Brown PetscReal nrm;
537c4762a1bSJed Brown PetscInt N;
538c4762a1bSJed Brown
5399566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(daf, &afexact));
5409566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(afexact));
5419566063dSJacob Faibussowitsch PetscCall(DADefineXLinearField3D(daf, afexact));
5429566063dSJacob Faibussowitsch PetscCall(VecAXPY(afexact, -1.0, af)); /* af <= af - afinterp */
5439566063dSJacob Faibussowitsch PetscCall(VecNorm(afexact, NORM_2, &nrm));
5449566063dSJacob Faibussowitsch PetscCall(VecGetSize(afexact, &N));
54563a3b9bcSJacob 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))));
5469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&afexact));
547c4762a1bSJed Brown }
548c4762a1bSJed Brown
5499566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-output", &output, NULL));
550c4762a1bSJed Brown if (output) {
5519566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "dac_3D.vtr", &vv));
5529566063dSJacob Faibussowitsch PetscCall(VecView(ac, vv));
5539566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vv));
554c4762a1bSJed Brown
5559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "daf_3D.vtr", &vv));
5569566063dSJacob Faibussowitsch PetscCall(VecView(af, vv));
5579566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vv));
558c4762a1bSJed Brown }
559c4762a1bSJed Brown
5609566063dSJacob Faibussowitsch PetscCall(MatDestroy(&INTERP));
5619566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dac));
5629566063dSJacob Faibussowitsch PetscCall(DMDestroy(&daf));
5639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ac));
5649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&af));
5653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
566c4762a1bSJed Brown }
567c4762a1bSJed Brown
main(int argc,char ** argv)568d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
569d71ae5a4SJacob Faibussowitsch {
570c4762a1bSJed Brown PetscInt mx = 2, my = 2, mz = 2, l, nl, dim;
571c4762a1bSJed Brown
572327415f7SBarry Smith PetscFunctionBeginUser;
5739566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, 0, help));
5749566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &mx, 0));
5759566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &my, 0));
5769566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mz", &mz, 0));
577c4762a1bSJed Brown nl = 1;
5789566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nl", &nl, 0));
579c4762a1bSJed Brown dim = 2;
5809566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, 0));
581c4762a1bSJed Brown
582c4762a1bSJed Brown for (l = 0; l < nl; l++) {
5831baa6e33SBarry Smith if (dim == 1) PetscCall(da_test_RefineCoords1D(mx));
5841baa6e33SBarry Smith else if (dim == 2) PetscCall(da_test_RefineCoords2D(mx, my));
5851baa6e33SBarry Smith else if (dim == 3) PetscCall(da_test_RefineCoords3D(mx, my, mz));
586c4762a1bSJed Brown mx = mx * 2;
587c4762a1bSJed Brown my = my * 2;
588c4762a1bSJed Brown mz = mz * 2;
589c4762a1bSJed Brown }
5909566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
591b122ec5aSJacob Faibussowitsch return 0;
592c4762a1bSJed Brown }
593c4762a1bSJed Brown
594c4762a1bSJed Brown /*TEST
595c4762a1bSJed Brown
596c4762a1bSJed Brown test:
597c4762a1bSJed Brown suffix: 1d
598c4762a1bSJed Brown args: -mx 10 -nl 6 -dim 1
599c4762a1bSJed Brown
600c4762a1bSJed Brown test:
601c4762a1bSJed Brown suffix: 2d
60263037964SStefano Zampini output_file: output/ex36_2d.out
60363037964SStefano Zampini args: -mx 10 -my 10 -nl 6 -dim 2 -cmap {{0 1 2 3}}
604c4762a1bSJed Brown
605c4762a1bSJed Brown test:
606c4762a1bSJed Brown suffix: 2dp1
607c4762a1bSJed Brown nsize: 8
608c4762a1bSJed Brown args: -mx 10 -my 10 -nl 4 -dim 2 -cmap 3 -da_refine_x 3 -da_refine_y 4
609c4762a1bSJed Brown timeoutfactor: 2
610c4762a1bSJed Brown
611c4762a1bSJed Brown test:
612c4762a1bSJed Brown suffix: 2dp2
613c4762a1bSJed Brown nsize: 8
614c4762a1bSJed Brown args: -mx 10 -my 10 -nl 4 -dim 2 -cmap 3 -da_refine_x 3 -da_refine_y 1
615c4762a1bSJed Brown timeoutfactor: 2
616c4762a1bSJed Brown
617c4762a1bSJed Brown test:
618c4762a1bSJed Brown suffix: 3d
619c4762a1bSJed Brown args: -mx 5 -my 5 -mz 5 -nl 4 -dim 3 -cmap 3
620c4762a1bSJed Brown
621c4762a1bSJed Brown test:
622c4762a1bSJed Brown suffix: 3dp1
623c4762a1bSJed Brown nsize: 32
624c4762a1bSJed 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
625c4762a1bSJed Brown
626c4762a1bSJed Brown TEST*/
627