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 13d71ae5a4SJacob Faibussowitsch CCmplx CCmplxPow(CCmplx a, PetscReal n) 14d71ae5a4SJacob Faibussowitsch { 15c4762a1bSJed Brown CCmplx b; 16c4762a1bSJed Brown PetscReal r, theta; 17c4762a1bSJed Brown r = PetscSqrtReal(a.real * a.real + a.imag * a.imag); 18c4762a1bSJed Brown theta = PetscAtan2Real(a.imag, a.real); 19c4762a1bSJed Brown b.real = PetscPowReal(r, n) * PetscCosReal(n * theta); 20c4762a1bSJed Brown b.imag = PetscPowReal(r, n) * PetscSinReal(n * theta); 21c4762a1bSJed Brown return b; 22c4762a1bSJed Brown } 23d71ae5a4SJacob Faibussowitsch CCmplx CCmplxExp(CCmplx a) 24d71ae5a4SJacob Faibussowitsch { 25c4762a1bSJed Brown CCmplx b; 26c4762a1bSJed Brown b.real = PetscExpReal(a.real) * PetscCosReal(a.imag); 27c4762a1bSJed Brown b.imag = PetscExpReal(a.real) * PetscSinReal(a.imag); 28c4762a1bSJed Brown return b; 29c4762a1bSJed Brown } 30d71ae5a4SJacob Faibussowitsch CCmplx CCmplxSqrt(CCmplx a) 31d71ae5a4SJacob Faibussowitsch { 32c4762a1bSJed Brown CCmplx b; 33c4762a1bSJed Brown PetscReal r, theta; 34c4762a1bSJed Brown r = PetscSqrtReal(a.real * a.real + a.imag * a.imag); 35c4762a1bSJed Brown theta = PetscAtan2Real(a.imag, a.real); 36c4762a1bSJed Brown b.real = PetscSqrtReal(r) * PetscCosReal(0.5 * theta); 37c4762a1bSJed Brown b.imag = PetscSqrtReal(r) * PetscSinReal(0.5 * theta); 38c4762a1bSJed Brown return b; 39c4762a1bSJed Brown } 40d71ae5a4SJacob Faibussowitsch CCmplx CCmplxAdd(CCmplx a, CCmplx c) 41d71ae5a4SJacob Faibussowitsch { 42c4762a1bSJed Brown CCmplx b; 43c4762a1bSJed Brown b.real = a.real + c.real; 44c4762a1bSJed Brown b.imag = a.imag + c.imag; 45c4762a1bSJed Brown return b; 46c4762a1bSJed Brown } 47d71ae5a4SJacob Faibussowitsch PetscScalar CCmplxRe(CCmplx a) 48d71ae5a4SJacob Faibussowitsch { 49c4762a1bSJed Brown return (PetscScalar)a.real; 50c4762a1bSJed Brown } 51d71ae5a4SJacob Faibussowitsch PetscScalar CCmplxIm(CCmplx a) 52d71ae5a4SJacob Faibussowitsch { 53c4762a1bSJed Brown return (PetscScalar)a.imag; 54c4762a1bSJed Brown } 55c4762a1bSJed Brown 56d71ae5a4SJacob Faibussowitsch PetscErrorCode DAApplyConformalMapping(DM da, PetscInt idx) 57d71ae5a4SJacob Faibussowitsch { 58c4762a1bSJed Brown PetscInt i, n; 59c4762a1bSJed Brown PetscInt sx, nx, sy, ny, sz, nz, dim; 60c4762a1bSJed Brown Vec Gcoords; 61c4762a1bSJed Brown PetscScalar *XX; 62c4762a1bSJed Brown PetscScalar xx, yy, zz; 63c4762a1bSJed Brown DM cda; 64c4762a1bSJed Brown 65c4762a1bSJed Brown PetscFunctionBeginUser; 66c4762a1bSJed Brown if (idx == 0) { 67*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 68c4762a1bSJed Brown } else if (idx == 1) { /* dam break */ 699566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0)); 70c4762a1bSJed Brown } else if (idx == 2) { /* stagnation in a corner */ 719566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, -1.0, 1.0, 0.0, 1.0, -1.0, 1.0)); 72c4762a1bSJed Brown } else if (idx == 3) { /* nautilis */ 739566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0)); 741baa6e33SBarry Smith } else if (idx == 4) PetscCall(DMDASetUniformCoordinates(da, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0)); 75c4762a1bSJed Brown 769566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &cda)); 779566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &Gcoords)); 78c4762a1bSJed Brown 799566063dSJacob Faibussowitsch PetscCall(VecGetArray(Gcoords, &XX)); 809566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &sx, &sy, &sz, &nx, &ny, &nz)); 819566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 829566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(Gcoords, &n)); 83c4762a1bSJed Brown n = n / dim; 84c4762a1bSJed Brown 85c4762a1bSJed Brown for (i = 0; i < n; i++) { 86c4762a1bSJed Brown if ((dim == 3) && (idx != 2)) { 87c4762a1bSJed Brown PetscScalar Ni[8]; 88c4762a1bSJed Brown PetscScalar xi = XX[dim * i]; 89c4762a1bSJed Brown PetscScalar eta = XX[dim * i + 1]; 90c4762a1bSJed Brown PetscScalar zeta = XX[dim * i + 2]; 91c4762a1bSJed Brown PetscScalar xn[] = {-1.0, 1.0, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0}; 92c4762a1bSJed Brown PetscScalar yn[] = {-1.0, -1.0, 1.0, 1.0, -1.0, -1.0, 1.0, 1.0}; 93c4762a1bSJed Brown PetscScalar zn[] = {-0.1, -4.0, -0.2, -1.0, 0.1, 4.0, 0.2, 1.0}; 94c4762a1bSJed Brown PetscInt p; 95c4762a1bSJed Brown 96c4762a1bSJed Brown Ni[0] = 0.125 * (1.0 - xi) * (1.0 - eta) * (1.0 - zeta); 97c4762a1bSJed Brown Ni[1] = 0.125 * (1.0 + xi) * (1.0 - eta) * (1.0 - zeta); 98c4762a1bSJed Brown Ni[2] = 0.125 * (1.0 - xi) * (1.0 + eta) * (1.0 - zeta); 99c4762a1bSJed Brown Ni[3] = 0.125 * (1.0 + xi) * (1.0 + eta) * (1.0 - zeta); 100c4762a1bSJed Brown 101c4762a1bSJed Brown Ni[4] = 0.125 * (1.0 - xi) * (1.0 - eta) * (1.0 + zeta); 102c4762a1bSJed Brown Ni[5] = 0.125 * (1.0 + xi) * (1.0 - eta) * (1.0 + zeta); 103c4762a1bSJed Brown Ni[6] = 0.125 * (1.0 - xi) * (1.0 + eta) * (1.0 + zeta); 104c4762a1bSJed Brown Ni[7] = 0.125 * (1.0 + xi) * (1.0 + eta) * (1.0 + zeta); 105c4762a1bSJed Brown 106c4762a1bSJed Brown xx = yy = zz = 0.0; 107c4762a1bSJed Brown for (p = 0; p < 8; p++) { 108c4762a1bSJed Brown xx += Ni[p] * xn[p]; 109c4762a1bSJed Brown yy += Ni[p] * yn[p]; 110c4762a1bSJed Brown zz += Ni[p] * zn[p]; 111c4762a1bSJed Brown } 112c4762a1bSJed Brown XX[dim * i] = xx; 113c4762a1bSJed Brown XX[dim * i + 1] = yy; 114c4762a1bSJed Brown XX[dim * i + 2] = zz; 115c4762a1bSJed Brown } 116c4762a1bSJed Brown 117c4762a1bSJed Brown if (idx == 1) { 118c4762a1bSJed Brown CCmplx zeta, t1, t2; 119c4762a1bSJed Brown 120c4762a1bSJed Brown xx = XX[dim * i] - 0.8; 121c4762a1bSJed Brown yy = XX[dim * i + 1] + 1.5; 122c4762a1bSJed Brown 123c4762a1bSJed Brown zeta.real = PetscRealPart(xx); 124c4762a1bSJed Brown zeta.imag = PetscRealPart(yy); 125c4762a1bSJed Brown 126c4762a1bSJed Brown t1 = CCmplxPow(zeta, -1.0); 127c4762a1bSJed Brown t2 = CCmplxAdd(zeta, t1); 128c4762a1bSJed Brown 129c4762a1bSJed Brown XX[dim * i] = CCmplxRe(t2); 130c4762a1bSJed Brown XX[dim * i + 1] = CCmplxIm(t2); 131c4762a1bSJed Brown } else if (idx == 2) { 132c4762a1bSJed Brown CCmplx zeta, t1; 133c4762a1bSJed Brown 134c4762a1bSJed Brown xx = XX[dim * i]; 135c4762a1bSJed Brown yy = XX[dim * i + 1]; 136c4762a1bSJed Brown zeta.real = PetscRealPart(xx); 137c4762a1bSJed Brown zeta.imag = PetscRealPart(yy); 138c4762a1bSJed Brown 139c4762a1bSJed Brown t1 = CCmplxSqrt(zeta); 140c4762a1bSJed Brown XX[dim * i] = CCmplxRe(t1); 141c4762a1bSJed Brown XX[dim * i + 1] = CCmplxIm(t1); 142c4762a1bSJed Brown } else if (idx == 3) { 143c4762a1bSJed Brown CCmplx zeta, t1, t2; 144c4762a1bSJed Brown 145c4762a1bSJed Brown xx = XX[dim * i] - 0.8; 146c4762a1bSJed Brown yy = XX[dim * i + 1] + 1.5; 147c4762a1bSJed Brown 148c4762a1bSJed Brown zeta.real = PetscRealPart(xx); 149c4762a1bSJed Brown zeta.imag = PetscRealPart(yy); 150c4762a1bSJed Brown t1 = CCmplxPow(zeta, -1.0); 151c4762a1bSJed Brown t2 = CCmplxAdd(zeta, t1); 152c4762a1bSJed Brown XX[dim * i] = CCmplxRe(t2); 153c4762a1bSJed Brown XX[dim * i + 1] = CCmplxIm(t2); 154c4762a1bSJed Brown 155c4762a1bSJed Brown xx = XX[dim * i]; 156c4762a1bSJed Brown yy = XX[dim * i + 1]; 157c4762a1bSJed Brown zeta.real = PetscRealPart(xx); 158c4762a1bSJed Brown zeta.imag = PetscRealPart(yy); 159c4762a1bSJed Brown t1 = CCmplxExp(zeta); 160c4762a1bSJed Brown XX[dim * i] = CCmplxRe(t1); 161c4762a1bSJed Brown XX[dim * i + 1] = CCmplxIm(t1); 162c4762a1bSJed Brown 163c4762a1bSJed Brown xx = XX[dim * i] + 0.4; 164c4762a1bSJed Brown yy = XX[dim * i + 1]; 165c4762a1bSJed Brown zeta.real = PetscRealPart(xx); 166c4762a1bSJed Brown zeta.imag = PetscRealPart(yy); 167c4762a1bSJed Brown t1 = CCmplxPow(zeta, 2.0); 168c4762a1bSJed Brown XX[dim * i] = CCmplxRe(t1); 169c4762a1bSJed Brown XX[dim * i + 1] = CCmplxIm(t1); 170c4762a1bSJed Brown } else if (idx == 4) { 171c4762a1bSJed Brown PetscScalar Ni[4]; 172c4762a1bSJed Brown PetscScalar xi = XX[dim * i]; 173c4762a1bSJed Brown PetscScalar eta = XX[dim * i + 1]; 174c4762a1bSJed Brown PetscScalar xn[] = {0.0, 2.0, 0.2, 3.5}; 175c4762a1bSJed Brown PetscScalar yn[] = {-1.3, 0.0, 2.0, 4.0}; 176c4762a1bSJed Brown PetscInt p; 177c4762a1bSJed Brown 178c4762a1bSJed Brown Ni[0] = 0.25 * (1.0 - xi) * (1.0 - eta); 179c4762a1bSJed Brown Ni[1] = 0.25 * (1.0 + xi) * (1.0 - eta); 180c4762a1bSJed Brown Ni[2] = 0.25 * (1.0 - xi) * (1.0 + eta); 181c4762a1bSJed Brown Ni[3] = 0.25 * (1.0 + xi) * (1.0 + eta); 182c4762a1bSJed Brown 183c4762a1bSJed Brown xx = yy = 0.0; 184c4762a1bSJed Brown for (p = 0; p < 4; p++) { 185c4762a1bSJed Brown xx += Ni[p] * xn[p]; 186c4762a1bSJed Brown yy += Ni[p] * yn[p]; 187c4762a1bSJed Brown } 188c4762a1bSJed Brown XX[dim * i] = xx; 189c4762a1bSJed Brown XX[dim * i + 1] = yy; 190c4762a1bSJed Brown } 191c4762a1bSJed Brown } 1929566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Gcoords, &XX)); 193*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 194c4762a1bSJed Brown } 195c4762a1bSJed Brown 196d71ae5a4SJacob Faibussowitsch PetscErrorCode DAApplyTrilinearMapping(DM da) 197d71ae5a4SJacob Faibussowitsch { 198c4762a1bSJed Brown PetscInt i, j, k; 199c4762a1bSJed Brown PetscInt sx, nx, sy, ny, sz, nz; 200c4762a1bSJed Brown Vec Gcoords; 201c4762a1bSJed Brown DMDACoor3d ***XX; 202c4762a1bSJed Brown PetscScalar xx, yy, zz; 203c4762a1bSJed Brown DM cda; 204c4762a1bSJed Brown 205c4762a1bSJed Brown PetscFunctionBeginUser; 2069566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0)); 2079566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &cda)); 2089566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &Gcoords)); 209c4762a1bSJed Brown 2109566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda, Gcoords, &XX)); 2119566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &sx, &sy, &sz, &nx, &ny, &nz)); 212c4762a1bSJed Brown 213c4762a1bSJed Brown for (i = sx; i < sx + nx; i++) { 214c4762a1bSJed Brown for (j = sy; j < sy + ny; j++) { 215c4762a1bSJed Brown for (k = sz; k < sz + nz; k++) { 216c4762a1bSJed Brown PetscScalar Ni[8]; 217c4762a1bSJed Brown PetscScalar xi = XX[k][j][i].x; 218c4762a1bSJed Brown PetscScalar eta = XX[k][j][i].y; 219c4762a1bSJed Brown PetscScalar zeta = XX[k][j][i].z; 220c4762a1bSJed Brown PetscScalar xn[] = {0.0, 2.0, 0.2, 3.5, 0.0, 2.1, 0.23, 3.125}; 221c4762a1bSJed Brown PetscScalar yn[] = {-1.3, 0.0, 2.0, 4.0, -1.45, -0.1, 2.24, 3.79}; 222c4762a1bSJed Brown PetscScalar zn[] = {0.0, 0.3, -0.1, 0.123, 0.956, 1.32, 1.12, 0.798}; 223c4762a1bSJed Brown PetscInt p; 224c4762a1bSJed Brown 225c4762a1bSJed Brown Ni[0] = 0.125 * (1.0 - xi) * (1.0 - eta) * (1.0 - zeta); 226c4762a1bSJed Brown Ni[1] = 0.125 * (1.0 + xi) * (1.0 - eta) * (1.0 - zeta); 227c4762a1bSJed Brown Ni[2] = 0.125 * (1.0 - xi) * (1.0 + eta) * (1.0 - zeta); 228c4762a1bSJed Brown Ni[3] = 0.125 * (1.0 + xi) * (1.0 + eta) * (1.0 - zeta); 229c4762a1bSJed Brown 230c4762a1bSJed Brown Ni[4] = 0.125 * (1.0 - xi) * (1.0 - eta) * (1.0 + zeta); 231c4762a1bSJed Brown Ni[5] = 0.125 * (1.0 + xi) * (1.0 - eta) * (1.0 + zeta); 232c4762a1bSJed Brown Ni[6] = 0.125 * (1.0 - xi) * (1.0 + eta) * (1.0 + zeta); 233c4762a1bSJed Brown Ni[7] = 0.125 * (1.0 + xi) * (1.0 + eta) * (1.0 + zeta); 234c4762a1bSJed Brown 235c4762a1bSJed Brown xx = yy = zz = 0.0; 236c4762a1bSJed Brown for (p = 0; p < 8; p++) { 237c4762a1bSJed Brown xx += Ni[p] * xn[p]; 238c4762a1bSJed Brown yy += Ni[p] * yn[p]; 239c4762a1bSJed Brown zz += Ni[p] * zn[p]; 240c4762a1bSJed Brown } 241c4762a1bSJed Brown XX[k][j][i].x = xx; 242c4762a1bSJed Brown XX[k][j][i].y = yy; 243c4762a1bSJed Brown XX[k][j][i].z = zz; 244c4762a1bSJed Brown } 245c4762a1bSJed Brown } 246c4762a1bSJed Brown } 2479566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda, Gcoords, &XX)); 248*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 249c4762a1bSJed Brown } 250c4762a1bSJed Brown 251d71ae5a4SJacob Faibussowitsch PetscErrorCode DADefineXLinearField2D(DM da, Vec field) 252d71ae5a4SJacob Faibussowitsch { 253c4762a1bSJed Brown PetscInt i, j; 254c4762a1bSJed Brown PetscInt sx, nx, sy, ny; 255c4762a1bSJed Brown Vec Gcoords; 256c4762a1bSJed Brown DMDACoor2d **XX; 257c4762a1bSJed Brown PetscScalar **FF; 258c4762a1bSJed Brown DM cda; 259c4762a1bSJed Brown 260c4762a1bSJed Brown PetscFunctionBeginUser; 2619566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &cda)); 2629566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &Gcoords)); 263c4762a1bSJed Brown 2649566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda, Gcoords, &XX)); 2659566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, field, &FF)); 266c4762a1bSJed Brown 2679566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &sx, &sy, 0, &nx, &ny, 0)); 268c4762a1bSJed Brown 269c4762a1bSJed Brown for (i = sx; i < sx + nx; i++) { 270ad540459SPierre 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; 271c4762a1bSJed Brown } 272c4762a1bSJed Brown 2739566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, field, &FF)); 2749566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda, Gcoords, &XX)); 275*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 276c4762a1bSJed Brown } 277c4762a1bSJed Brown 278d71ae5a4SJacob Faibussowitsch PetscErrorCode DADefineXLinearField3D(DM da, Vec field) 279d71ae5a4SJacob Faibussowitsch { 280c4762a1bSJed Brown PetscInt i, j, k; 281c4762a1bSJed Brown PetscInt sx, nx, sy, ny, sz, nz; 282c4762a1bSJed Brown Vec Gcoords; 283c4762a1bSJed Brown DMDACoor3d ***XX; 284c4762a1bSJed Brown PetscScalar ***FF; 285c4762a1bSJed Brown DM cda; 286c4762a1bSJed Brown 287c4762a1bSJed Brown PetscFunctionBeginUser; 2889566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &cda)); 2899566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &Gcoords)); 290c4762a1bSJed Brown 2919566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda, Gcoords, &XX)); 2929566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, field, &FF)); 293c4762a1bSJed Brown 2949566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &sx, &sy, &sz, &nx, &ny, &nz)); 295c4762a1bSJed Brown 296c4762a1bSJed Brown for (k = sz; k < sz + nz; k++) { 297c4762a1bSJed Brown for (j = sy; j < sy + ny; j++) { 298c4762a1bSJed Brown for (i = sx; i < sx + nx; i++) { 2999371c9d4SSatish 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 + 3009371c9d4SSatish Balay 3.79 * XX[k][j][i].x * XX[k][j][i].y * XX[k][j][i].z; 301c4762a1bSJed Brown } 302c4762a1bSJed Brown } 303c4762a1bSJed Brown } 304c4762a1bSJed Brown 3059566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, field, &FF)); 3069566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda, Gcoords, &XX)); 307*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 308c4762a1bSJed Brown } 309c4762a1bSJed Brown 310d71ae5a4SJacob Faibussowitsch PetscErrorCode da_test_RefineCoords1D(PetscInt mx) 311d71ae5a4SJacob Faibussowitsch { 312c4762a1bSJed Brown DM dac, daf; 313c4762a1bSJed Brown PetscViewer vv; 314c4762a1bSJed Brown Vec ac, af; 315c4762a1bSJed Brown PetscInt Mx; 316c4762a1bSJed Brown Mat II, INTERP; 317c4762a1bSJed Brown Vec scale; 318c4762a1bSJed Brown PetscBool output = PETSC_FALSE; 319c4762a1bSJed Brown 320c4762a1bSJed Brown PetscFunctionBeginUser; 3219566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, mx + 1, 1, /* 1 dof */ 1, /* stencil = 1 */ NULL, &dac)); 3229566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dac)); 3239566063dSJacob Faibussowitsch PetscCall(DMSetUp(dac)); 324c4762a1bSJed Brown 3259566063dSJacob Faibussowitsch PetscCall(DMRefine(dac, MPI_COMM_NULL, &daf)); 3269566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, 0, &Mx, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 327c4762a1bSJed Brown Mx--; 328c4762a1bSJed Brown 3299566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(dac, -1.0, 1.0, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE)); 3309566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(daf, -1.0, 1.0, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE)); 331c4762a1bSJed Brown 332c4762a1bSJed Brown { 333c4762a1bSJed Brown DM cdaf, cdac; 334c4762a1bSJed Brown Vec coordsc, coordsf; 335c4762a1bSJed Brown 3369566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dac, &cdac)); 3379566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(daf, &cdaf)); 338c4762a1bSJed Brown 3399566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dac, &coordsc)); 3409566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(daf, &coordsf)); 341c4762a1bSJed Brown 3429566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(cdac, cdaf, &II, &scale)); 3439566063dSJacob Faibussowitsch PetscCall(MatInterpolate(II, coordsc, coordsf)); 3449566063dSJacob Faibussowitsch PetscCall(MatDestroy(&II)); 3459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&scale)); 346c4762a1bSJed Brown } 347c4762a1bSJed Brown 3489566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dac, daf, &INTERP, NULL)); 349c4762a1bSJed Brown 3509566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dac, &ac)); 3519566063dSJacob Faibussowitsch PetscCall(VecSet(ac, 66.99)); 352c4762a1bSJed Brown 3539566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(daf, &af)); 354c4762a1bSJed Brown 3559566063dSJacob Faibussowitsch PetscCall(MatMult(INTERP, ac, af)); 356c4762a1bSJed Brown 357c4762a1bSJed Brown { 358c4762a1bSJed Brown Vec afexact; 359c4762a1bSJed Brown PetscReal nrm; 360c4762a1bSJed Brown PetscInt N; 361c4762a1bSJed Brown 3629566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(daf, &afexact)); 3639566063dSJacob Faibussowitsch PetscCall(VecSet(afexact, 66.99)); 3649566063dSJacob Faibussowitsch PetscCall(VecAXPY(afexact, -1.0, af)); /* af <= af - afinterp */ 3659566063dSJacob Faibussowitsch PetscCall(VecNorm(afexact, NORM_2, &nrm)); 3669566063dSJacob Faibussowitsch PetscCall(VecGetSize(afexact, &N)); 36763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "=>%" PetscInt_FMT ", interp err = %1.4e\n", mx, Mx, (double)(nrm / PetscSqrtReal((PetscReal)N)))); 3689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&afexact)); 369c4762a1bSJed Brown } 370c4762a1bSJed Brown 3719566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-output", &output, NULL)); 372c4762a1bSJed Brown if (output) { 3739566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "dac_1D.vtr", &vv)); 3749566063dSJacob Faibussowitsch PetscCall(VecView(ac, vv)); 3759566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vv)); 376c4762a1bSJed Brown 3779566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "daf_1D.vtr", &vv)); 3789566063dSJacob Faibussowitsch PetscCall(VecView(af, vv)); 3799566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vv)); 380c4762a1bSJed Brown } 381c4762a1bSJed Brown 3829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&INTERP)); 3839566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dac)); 3849566063dSJacob Faibussowitsch PetscCall(DMDestroy(&daf)); 3859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ac)); 3869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&af)); 387*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 388c4762a1bSJed Brown } 389c4762a1bSJed Brown 390d71ae5a4SJacob Faibussowitsch PetscErrorCode da_test_RefineCoords2D(PetscInt mx, PetscInt my) 391d71ae5a4SJacob Faibussowitsch { 392c4762a1bSJed Brown DM dac, daf; 393c4762a1bSJed Brown PetscViewer vv; 394c4762a1bSJed Brown Vec ac, af; 395c4762a1bSJed Brown PetscInt map_id, Mx, My; 396c4762a1bSJed Brown Mat II, INTERP; 397c4762a1bSJed Brown Vec scale; 398c4762a1bSJed Brown PetscBool output = PETSC_FALSE; 399c4762a1bSJed Brown 400c4762a1bSJed Brown PetscFunctionBeginUser; 4019566063dSJacob 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)); 4029566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dac)); 4039566063dSJacob Faibussowitsch PetscCall(DMSetUp(dac)); 404c4762a1bSJed Brown 4059566063dSJacob Faibussowitsch PetscCall(DMRefine(dac, MPI_COMM_NULL, &daf)); 4069566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, 0, &Mx, &My, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 4079371c9d4SSatish Balay Mx--; 4089371c9d4SSatish Balay My--; 409c4762a1bSJed Brown 4109566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(dac, -1.0, 1.0, -1.0, 1.0, PETSC_DECIDE, PETSC_DECIDE)); 4119566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(daf, -1.0, 1.0, -1.0, 1.0, PETSC_DECIDE, PETSC_DECIDE)); 412c4762a1bSJed Brown 413c4762a1bSJed Brown /* apply conformal mappings */ 414c4762a1bSJed Brown map_id = 0; 4159566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-cmap", &map_id, NULL)); 41648a46eb9SPierre Jolivet if (map_id >= 1) PetscCall(DAApplyConformalMapping(dac, map_id)); 417c4762a1bSJed Brown 418c4762a1bSJed Brown { 419c4762a1bSJed Brown DM cdaf, cdac; 420c4762a1bSJed Brown Vec coordsc, coordsf; 421c4762a1bSJed Brown 4229566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dac, &cdac)); 4239566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(daf, &cdaf)); 424c4762a1bSJed Brown 4259566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dac, &coordsc)); 4269566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(daf, &coordsf)); 427c4762a1bSJed Brown 4289566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(cdac, cdaf, &II, &scale)); 4299566063dSJacob Faibussowitsch PetscCall(MatInterpolate(II, coordsc, coordsf)); 4309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&II)); 4319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&scale)); 432c4762a1bSJed Brown } 433c4762a1bSJed Brown 4349566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dac, daf, &INTERP, NULL)); 435c4762a1bSJed Brown 4369566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dac, &ac)); 4379566063dSJacob Faibussowitsch PetscCall(DADefineXLinearField2D(dac, ac)); 438c4762a1bSJed Brown 4399566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(daf, &af)); 4409566063dSJacob Faibussowitsch PetscCall(MatMult(INTERP, ac, af)); 441c4762a1bSJed Brown 442c4762a1bSJed Brown { 443c4762a1bSJed Brown Vec afexact; 444c4762a1bSJed Brown PetscReal nrm; 445c4762a1bSJed Brown PetscInt N; 446c4762a1bSJed Brown 4479566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(daf, &afexact)); 4489566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(afexact)); 4499566063dSJacob Faibussowitsch PetscCall(DADefineXLinearField2D(daf, afexact)); 4509566063dSJacob Faibussowitsch PetscCall(VecAXPY(afexact, -1.0, af)); /* af <= af - afinterp */ 4519566063dSJacob Faibussowitsch PetscCall(VecNorm(afexact, NORM_2, &nrm)); 4529566063dSJacob Faibussowitsch PetscCall(VecGetSize(afexact, &N)); 45363a3b9bcSJacob 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)))); 4549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&afexact)); 455c4762a1bSJed Brown } 456c4762a1bSJed Brown 4579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-output", &output, NULL)); 458c4762a1bSJed Brown if (output) { 4599566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "dac_2D.vtr", &vv)); 4609566063dSJacob Faibussowitsch PetscCall(VecView(ac, vv)); 4619566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vv)); 462c4762a1bSJed Brown 4639566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "daf_2D.vtr", &vv)); 4649566063dSJacob Faibussowitsch PetscCall(VecView(af, vv)); 4659566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vv)); 466c4762a1bSJed Brown } 467c4762a1bSJed Brown 4689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&INTERP)); 4699566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dac)); 4709566063dSJacob Faibussowitsch PetscCall(DMDestroy(&daf)); 4719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ac)); 4729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&af)); 473*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 474c4762a1bSJed Brown } 475c4762a1bSJed Brown 476d71ae5a4SJacob Faibussowitsch PetscErrorCode da_test_RefineCoords3D(PetscInt mx, PetscInt my, PetscInt mz) 477d71ae5a4SJacob Faibussowitsch { 478c4762a1bSJed Brown DM dac, daf; 479c4762a1bSJed Brown PetscViewer vv; 480c4762a1bSJed Brown Vec ac, af; 481c4762a1bSJed Brown PetscInt map_id, Mx, My, Mz; 482c4762a1bSJed Brown Mat II, INTERP; 483c4762a1bSJed Brown Vec scale; 484c4762a1bSJed Brown PetscBool output = PETSC_FALSE; 485c4762a1bSJed Brown 486c4762a1bSJed Brown PetscFunctionBeginUser; 487d0609cedSBarry 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 */ 488d0609cedSBarry Smith 1, /* stencil = 1 */ NULL, NULL, NULL, &dac)); 4899566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dac)); 4909566063dSJacob Faibussowitsch PetscCall(DMSetUp(dac)); 491c4762a1bSJed Brown 4929566063dSJacob Faibussowitsch PetscCall(DMRefine(dac, MPI_COMM_NULL, &daf)); 4939566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf, 0, &Mx, &My, &Mz, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 4949371c9d4SSatish Balay Mx--; 4959371c9d4SSatish Balay My--; 4969371c9d4SSatish Balay Mz--; 497c4762a1bSJed Brown 4989566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(dac, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0)); 4999566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(daf, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0)); 500c4762a1bSJed Brown 501c4762a1bSJed Brown /* apply trilinear mappings */ 5029566063dSJacob Faibussowitsch /*PetscCall(DAApplyTrilinearMapping(dac));*/ 503c4762a1bSJed Brown /* apply conformal mappings */ 504c4762a1bSJed Brown map_id = 0; 5059566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-cmap", &map_id, NULL)); 50648a46eb9SPierre Jolivet if (map_id >= 1) PetscCall(DAApplyConformalMapping(dac, map_id)); 507c4762a1bSJed Brown 508c4762a1bSJed Brown { 509c4762a1bSJed Brown DM cdaf, cdac; 510c4762a1bSJed Brown Vec coordsc, coordsf; 511c4762a1bSJed Brown 5129566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dac, &cdac)); 5139566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(daf, &cdaf)); 514c4762a1bSJed Brown 5159566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dac, &coordsc)); 5169566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(daf, &coordsf)); 517c4762a1bSJed Brown 5189566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(cdac, cdaf, &II, &scale)); 5199566063dSJacob Faibussowitsch PetscCall(MatInterpolate(II, coordsc, coordsf)); 5209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&II)); 5219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&scale)); 522c4762a1bSJed Brown } 523c4762a1bSJed Brown 5249566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dac, daf, &INTERP, NULL)); 525c4762a1bSJed Brown 5269566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dac, &ac)); 5279566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(ac)); 5289566063dSJacob Faibussowitsch PetscCall(DADefineXLinearField3D(dac, ac)); 529c4762a1bSJed Brown 5309566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(daf, &af)); 5319566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(af)); 532c4762a1bSJed Brown 5339566063dSJacob Faibussowitsch PetscCall(MatMult(INTERP, ac, af)); 534c4762a1bSJed Brown 535c4762a1bSJed Brown { 536c4762a1bSJed Brown Vec afexact; 537c4762a1bSJed Brown PetscReal nrm; 538c4762a1bSJed Brown PetscInt N; 539c4762a1bSJed Brown 5409566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(daf, &afexact)); 5419566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(afexact)); 5429566063dSJacob Faibussowitsch PetscCall(DADefineXLinearField3D(daf, afexact)); 5439566063dSJacob Faibussowitsch PetscCall(VecAXPY(afexact, -1.0, af)); /* af <= af - afinterp */ 5449566063dSJacob Faibussowitsch PetscCall(VecNorm(afexact, NORM_2, &nrm)); 5459566063dSJacob Faibussowitsch PetscCall(VecGetSize(afexact, &N)); 54663a3b9bcSJacob 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)))); 5479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&afexact)); 548c4762a1bSJed Brown } 549c4762a1bSJed Brown 5509566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-output", &output, NULL)); 551c4762a1bSJed Brown if (output) { 5529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "dac_3D.vtr", &vv)); 5539566063dSJacob Faibussowitsch PetscCall(VecView(ac, vv)); 5549566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vv)); 555c4762a1bSJed Brown 5569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "daf_3D.vtr", &vv)); 5579566063dSJacob Faibussowitsch PetscCall(VecView(af, vv)); 5589566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vv)); 559c4762a1bSJed Brown } 560c4762a1bSJed Brown 5619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&INTERP)); 5629566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dac)); 5639566063dSJacob Faibussowitsch PetscCall(DMDestroy(&daf)); 5649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ac)); 5659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&af)); 566*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 567c4762a1bSJed Brown } 568c4762a1bSJed Brown 569d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 570d71ae5a4SJacob Faibussowitsch { 571c4762a1bSJed Brown PetscInt mx = 2, my = 2, mz = 2, l, nl, dim; 572c4762a1bSJed Brown 573327415f7SBarry Smith PetscFunctionBeginUser; 5749566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, 0, help)); 5759566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &mx, 0)); 5769566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &my, 0)); 5779566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mz", &mz, 0)); 578c4762a1bSJed Brown nl = 1; 5799566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nl", &nl, 0)); 580c4762a1bSJed Brown dim = 2; 5819566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, 0)); 582c4762a1bSJed Brown 583c4762a1bSJed Brown for (l = 0; l < nl; l++) { 5841baa6e33SBarry Smith if (dim == 1) PetscCall(da_test_RefineCoords1D(mx)); 5851baa6e33SBarry Smith else if (dim == 2) PetscCall(da_test_RefineCoords2D(mx, my)); 5861baa6e33SBarry Smith else if (dim == 3) PetscCall(da_test_RefineCoords3D(mx, my, mz)); 587c4762a1bSJed Brown mx = mx * 2; 588c4762a1bSJed Brown my = my * 2; 589c4762a1bSJed Brown mz = mz * 2; 590c4762a1bSJed Brown } 5919566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 592b122ec5aSJacob Faibussowitsch return 0; 593c4762a1bSJed Brown } 594c4762a1bSJed Brown 595c4762a1bSJed Brown /*TEST 596c4762a1bSJed Brown 597c4762a1bSJed Brown test: 598c4762a1bSJed Brown suffix: 1d 599c4762a1bSJed Brown args: -mx 10 -nl 6 -dim 1 600c4762a1bSJed Brown 601c4762a1bSJed Brown test: 602c4762a1bSJed Brown suffix: 2d 60363037964SStefano Zampini output_file: output/ex36_2d.out 60463037964SStefano Zampini args: -mx 10 -my 10 -nl 6 -dim 2 -cmap {{0 1 2 3}} 605c4762a1bSJed Brown 606c4762a1bSJed Brown test: 607c4762a1bSJed Brown suffix: 2dp1 608c4762a1bSJed Brown nsize: 8 609c4762a1bSJed Brown args: -mx 10 -my 10 -nl 4 -dim 2 -cmap 3 -da_refine_x 3 -da_refine_y 4 610c4762a1bSJed Brown timeoutfactor: 2 611c4762a1bSJed Brown 612c4762a1bSJed Brown test: 613c4762a1bSJed Brown suffix: 2dp2 614c4762a1bSJed Brown nsize: 8 615c4762a1bSJed Brown args: -mx 10 -my 10 -nl 4 -dim 2 -cmap 3 -da_refine_x 3 -da_refine_y 1 616c4762a1bSJed Brown timeoutfactor: 2 617c4762a1bSJed Brown 618c4762a1bSJed Brown test: 619c4762a1bSJed Brown suffix: 3d 620c4762a1bSJed Brown args: -mx 5 -my 5 -mz 5 -nl 4 -dim 3 -cmap 3 621c4762a1bSJed Brown 622c4762a1bSJed Brown test: 623c4762a1bSJed Brown suffix: 3dp1 624c4762a1bSJed Brown nsize: 32 625c4762a1bSJed 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 626c4762a1bSJed Brown 627c4762a1bSJed Brown TEST*/ 628