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 13c4762a1bSJed Brown CCmplx CCmplxPow(CCmplx a,PetscReal n) 14c4762a1bSJed Brown { 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 } 23c4762a1bSJed Brown CCmplx CCmplxExp(CCmplx a) 24c4762a1bSJed Brown { 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 } 30c4762a1bSJed Brown CCmplx CCmplxSqrt(CCmplx a) 31c4762a1bSJed Brown { 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 } 40c4762a1bSJed Brown CCmplx CCmplxAdd(CCmplx a,CCmplx c) 41c4762a1bSJed Brown { 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 } 47c4762a1bSJed Brown PetscScalar CCmplxRe(CCmplx a) 48c4762a1bSJed Brown { 49c4762a1bSJed Brown return (PetscScalar)a.real; 50c4762a1bSJed Brown } 51c4762a1bSJed Brown PetscScalar CCmplxIm(CCmplx a) 52c4762a1bSJed Brown { 53c4762a1bSJed Brown return (PetscScalar)a.imag; 54c4762a1bSJed Brown } 55c4762a1bSJed Brown 56c4762a1bSJed Brown PetscErrorCode DAApplyConformalMapping(DM da,PetscInt idx) 57c4762a1bSJed Brown { 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) { 67c4762a1bSJed Brown PetscFunctionReturn(0); 68c4762a1bSJed Brown } else if (idx==1) { /* dam break */ 695f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetUniformCoordinates(da, -1.0,1.0, -1.0,1.0, -1.0,1.0)); 70c4762a1bSJed Brown } else if (idx==2) { /* stagnation in a corner */ 715f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetUniformCoordinates(da, -1.0,1.0, 0.0,1.0, -1.0,1.0)); 72c4762a1bSJed Brown } else if (idx==3) { /* nautilis */ 735f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetUniformCoordinates(da, -1.0,1.0, -1.0,1.0, -1.0,1.0)); 74c4762a1bSJed Brown } else if (idx==4) { 755f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetUniformCoordinates(da, -1.0,1.0, -1.0,1.0, -1.0,1.0)); 76c4762a1bSJed Brown } 77c4762a1bSJed Brown 785f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(da,&cda)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinates(da,&Gcoords)); 80c4762a1bSJed Brown 815f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(Gcoords,&XX)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&sx,&sy,&sz,&nx,&ny,&nz)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da, &dim, 0,0,0, 0,0,0, 0,0,0,0,0,0)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(Gcoords,&n)); 85c4762a1bSJed Brown n = n / dim; 86c4762a1bSJed Brown 87c4762a1bSJed Brown for (i=0; i<n; i++) { 88c4762a1bSJed Brown if ((dim==3) && (idx!=2)) { 89c4762a1bSJed Brown PetscScalar Ni[8]; 90c4762a1bSJed Brown PetscScalar xi = XX[dim*i]; 91c4762a1bSJed Brown PetscScalar eta = XX[dim*i+1]; 92c4762a1bSJed Brown PetscScalar zeta = XX[dim*i+2]; 93c4762a1bSJed Brown PetscScalar xn[] = {-1.0,1.0,-1.0,1.0, -1.0,1.0,-1.0,1.0 }; 94c4762a1bSJed Brown PetscScalar yn[] = {-1.0,-1.0,1.0,1.0, -1.0,-1.0,1.0,1.0 }; 95c4762a1bSJed Brown PetscScalar zn[] = {-0.1,-4.0,-0.2,-1.0, 0.1,4.0,0.2,1.0 }; 96c4762a1bSJed Brown PetscInt p; 97c4762a1bSJed Brown 98c4762a1bSJed Brown Ni[0] = 0.125*(1.0-xi)*(1.0-eta)*(1.0-zeta); 99c4762a1bSJed Brown Ni[1] = 0.125*(1.0+xi)*(1.0-eta)*(1.0-zeta); 100c4762a1bSJed Brown Ni[2] = 0.125*(1.0-xi)*(1.0+eta)*(1.0-zeta); 101c4762a1bSJed Brown Ni[3] = 0.125*(1.0+xi)*(1.0+eta)*(1.0-zeta); 102c4762a1bSJed Brown 103c4762a1bSJed Brown Ni[4] = 0.125*(1.0-xi)*(1.0-eta)*(1.0+zeta); 104c4762a1bSJed Brown Ni[5] = 0.125*(1.0+xi)*(1.0-eta)*(1.0+zeta); 105c4762a1bSJed Brown Ni[6] = 0.125*(1.0-xi)*(1.0+eta)*(1.0+zeta); 106c4762a1bSJed Brown Ni[7] = 0.125*(1.0+xi)*(1.0+eta)*(1.0+zeta); 107c4762a1bSJed Brown 108c4762a1bSJed Brown xx = yy = zz = 0.0; 109c4762a1bSJed Brown for (p=0; p<8; p++) { 110c4762a1bSJed Brown xx += Ni[p]*xn[p]; 111c4762a1bSJed Brown yy += Ni[p]*yn[p]; 112c4762a1bSJed Brown zz += Ni[p]*zn[p]; 113c4762a1bSJed Brown } 114c4762a1bSJed Brown XX[dim*i] = xx; 115c4762a1bSJed Brown XX[dim*i+1] = yy; 116c4762a1bSJed Brown XX[dim*i+2] = zz; 117c4762a1bSJed Brown } 118c4762a1bSJed Brown 119c4762a1bSJed Brown if (idx==1) { 120c4762a1bSJed Brown CCmplx zeta,t1,t2; 121c4762a1bSJed Brown 122c4762a1bSJed Brown xx = XX[dim*i] - 0.8; 123c4762a1bSJed Brown yy = XX[dim*i+1] + 1.5; 124c4762a1bSJed Brown 125c4762a1bSJed Brown zeta.real = PetscRealPart(xx); 126c4762a1bSJed Brown zeta.imag = PetscRealPart(yy); 127c4762a1bSJed Brown 128c4762a1bSJed Brown t1 = CCmplxPow(zeta,-1.0); 129c4762a1bSJed Brown t2 = CCmplxAdd(zeta,t1); 130c4762a1bSJed Brown 131c4762a1bSJed Brown XX[dim*i] = CCmplxRe(t2); 132c4762a1bSJed Brown XX[dim*i+1] = CCmplxIm(t2); 133c4762a1bSJed Brown } else if (idx==2) { 134c4762a1bSJed Brown CCmplx zeta,t1; 135c4762a1bSJed Brown 136c4762a1bSJed Brown xx = XX[dim*i]; 137c4762a1bSJed Brown yy = XX[dim*i+1]; 138c4762a1bSJed Brown zeta.real = PetscRealPart(xx); 139c4762a1bSJed Brown zeta.imag = PetscRealPart(yy); 140c4762a1bSJed Brown 141c4762a1bSJed Brown t1 = CCmplxSqrt(zeta); 142c4762a1bSJed Brown XX[dim*i] = CCmplxRe(t1); 143c4762a1bSJed Brown XX[dim*i+1] = CCmplxIm(t1); 144c4762a1bSJed Brown } else if (idx==3) { 145c4762a1bSJed Brown CCmplx zeta,t1,t2; 146c4762a1bSJed Brown 147c4762a1bSJed Brown xx = XX[dim*i] - 0.8; 148c4762a1bSJed Brown yy = XX[dim*i+1] + 1.5; 149c4762a1bSJed Brown 150c4762a1bSJed Brown zeta.real = PetscRealPart(xx); 151c4762a1bSJed Brown zeta.imag = PetscRealPart(yy); 152c4762a1bSJed Brown t1 = CCmplxPow(zeta,-1.0); 153c4762a1bSJed Brown t2 = CCmplxAdd(zeta,t1); 154c4762a1bSJed Brown XX[dim*i] = CCmplxRe(t2); 155c4762a1bSJed Brown XX[dim*i+1] = CCmplxIm(t2); 156c4762a1bSJed Brown 157c4762a1bSJed Brown xx = XX[dim*i]; 158c4762a1bSJed Brown yy = XX[dim*i+1]; 159c4762a1bSJed Brown zeta.real = PetscRealPart(xx); 160c4762a1bSJed Brown zeta.imag = PetscRealPart(yy); 161c4762a1bSJed Brown t1 = CCmplxExp(zeta); 162c4762a1bSJed Brown XX[dim*i] = CCmplxRe(t1); 163c4762a1bSJed Brown XX[dim*i+1] = CCmplxIm(t1); 164c4762a1bSJed Brown 165c4762a1bSJed Brown xx = XX[dim*i] + 0.4; 166c4762a1bSJed Brown yy = XX[dim*i+1]; 167c4762a1bSJed Brown zeta.real = PetscRealPart(xx); 168c4762a1bSJed Brown zeta.imag = PetscRealPart(yy); 169c4762a1bSJed Brown t1 = CCmplxPow(zeta,2.0); 170c4762a1bSJed Brown XX[dim*i] = CCmplxRe(t1); 171c4762a1bSJed Brown XX[dim*i+1] = CCmplxIm(t1); 172c4762a1bSJed Brown } else if (idx==4) { 173c4762a1bSJed Brown PetscScalar Ni[4]; 174c4762a1bSJed Brown PetscScalar xi = XX[dim*i]; 175c4762a1bSJed Brown PetscScalar eta = XX[dim*i+1]; 176c4762a1bSJed Brown PetscScalar xn[] = {0.0,2.0,0.2,3.5}; 177c4762a1bSJed Brown PetscScalar yn[] = {-1.3,0.0,2.0,4.0}; 178c4762a1bSJed Brown PetscInt p; 179c4762a1bSJed Brown 180c4762a1bSJed Brown Ni[0] = 0.25*(1.0-xi)*(1.0-eta); 181c4762a1bSJed Brown Ni[1] = 0.25*(1.0+xi)*(1.0-eta); 182c4762a1bSJed Brown Ni[2] = 0.25*(1.0-xi)*(1.0+eta); 183c4762a1bSJed Brown Ni[3] = 0.25*(1.0+xi)*(1.0+eta); 184c4762a1bSJed Brown 185c4762a1bSJed Brown xx = yy = 0.0; 186c4762a1bSJed Brown for (p=0; p<4; p++) { 187c4762a1bSJed Brown xx += Ni[p]*xn[p]; 188c4762a1bSJed Brown yy += Ni[p]*yn[p]; 189c4762a1bSJed Brown } 190c4762a1bSJed Brown XX[dim*i] = xx; 191c4762a1bSJed Brown XX[dim*i+1] = yy; 192c4762a1bSJed Brown } 193c4762a1bSJed Brown } 1945f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(Gcoords,&XX)); 195c4762a1bSJed Brown PetscFunctionReturn(0); 196c4762a1bSJed Brown } 197c4762a1bSJed Brown 198c4762a1bSJed Brown PetscErrorCode DAApplyTrilinearMapping(DM da) 199c4762a1bSJed Brown { 200c4762a1bSJed Brown PetscInt i,j,k; 201c4762a1bSJed Brown PetscInt sx,nx,sy,ny,sz,nz; 202c4762a1bSJed Brown Vec Gcoords; 203c4762a1bSJed Brown DMDACoor3d ***XX; 204c4762a1bSJed Brown PetscScalar xx,yy,zz; 205c4762a1bSJed Brown DM cda; 206c4762a1bSJed Brown 207c4762a1bSJed Brown PetscFunctionBeginUser; 2085f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetUniformCoordinates(da, -1.0,1.0, -1.0,1.0, -1.0,1.0)); 2095f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(da,&cda)); 2105f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinates(da,&Gcoords)); 211c4762a1bSJed Brown 2125f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(cda,Gcoords,&XX)); 2135f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&sx,&sy,&sz,&nx,&ny,&nz)); 214c4762a1bSJed Brown 215c4762a1bSJed Brown for (i=sx; i<sx+nx; i++) { 216c4762a1bSJed Brown for (j=sy; j<sy+ny; j++) { 217c4762a1bSJed Brown for (k=sz; k<sz+nz; k++) { 218c4762a1bSJed Brown PetscScalar Ni[8]; 219c4762a1bSJed Brown PetscScalar xi = XX[k][j][i].x; 220c4762a1bSJed Brown PetscScalar eta = XX[k][j][i].y; 221c4762a1bSJed Brown PetscScalar zeta = XX[k][j][i].z; 222c4762a1bSJed Brown PetscScalar xn[] = {0.0,2.0,0.2,3.5, 0.0,2.1,0.23,3.125 }; 223c4762a1bSJed Brown PetscScalar yn[] = {-1.3,0.0,2.0,4.0, -1.45,-0.1,2.24,3.79 }; 224c4762a1bSJed Brown PetscScalar zn[] = {0.0,0.3,-0.1,0.123, 0.956,1.32,1.12,0.798 }; 225c4762a1bSJed Brown PetscInt p; 226c4762a1bSJed Brown 227c4762a1bSJed Brown Ni[0] = 0.125*(1.0-xi)*(1.0-eta)*(1.0-zeta); 228c4762a1bSJed Brown Ni[1] = 0.125*(1.0+xi)*(1.0-eta)*(1.0-zeta); 229c4762a1bSJed Brown Ni[2] = 0.125*(1.0-xi)*(1.0+eta)*(1.0-zeta); 230c4762a1bSJed Brown Ni[3] = 0.125*(1.0+xi)*(1.0+eta)*(1.0-zeta); 231c4762a1bSJed Brown 232c4762a1bSJed Brown Ni[4] = 0.125*(1.0-xi)*(1.0-eta)*(1.0+zeta); 233c4762a1bSJed Brown Ni[5] = 0.125*(1.0+xi)*(1.0-eta)*(1.0+zeta); 234c4762a1bSJed Brown Ni[6] = 0.125*(1.0-xi)*(1.0+eta)*(1.0+zeta); 235c4762a1bSJed Brown Ni[7] = 0.125*(1.0+xi)*(1.0+eta)*(1.0+zeta); 236c4762a1bSJed Brown 237c4762a1bSJed Brown xx = yy = zz = 0.0; 238c4762a1bSJed Brown for (p=0; p<8; p++) { 239c4762a1bSJed Brown xx += Ni[p]*xn[p]; 240c4762a1bSJed Brown yy += Ni[p]*yn[p]; 241c4762a1bSJed Brown zz += Ni[p]*zn[p]; 242c4762a1bSJed Brown } 243c4762a1bSJed Brown XX[k][j][i].x = xx; 244c4762a1bSJed Brown XX[k][j][i].y = yy; 245c4762a1bSJed Brown XX[k][j][i].z = zz; 246c4762a1bSJed Brown } 247c4762a1bSJed Brown } 248c4762a1bSJed Brown } 2495f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(cda,Gcoords,&XX)); 250c4762a1bSJed Brown PetscFunctionReturn(0); 251c4762a1bSJed Brown } 252c4762a1bSJed Brown 253c4762a1bSJed Brown PetscErrorCode DADefineXLinearField2D(DM da,Vec field) 254c4762a1bSJed Brown { 255c4762a1bSJed Brown PetscInt i,j; 256c4762a1bSJed Brown PetscInt sx,nx,sy,ny; 257c4762a1bSJed Brown Vec Gcoords; 258c4762a1bSJed Brown DMDACoor2d **XX; 259c4762a1bSJed Brown PetscScalar **FF; 260c4762a1bSJed Brown DM cda; 261c4762a1bSJed Brown 262c4762a1bSJed Brown PetscFunctionBeginUser; 2635f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(da,&cda)); 2645f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinates(da,&Gcoords)); 265c4762a1bSJed Brown 2665f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(cda,Gcoords,&XX)); 2675f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,field,&FF)); 268c4762a1bSJed Brown 2695f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&sx,&sy,0,&nx,&ny,0)); 270c4762a1bSJed Brown 271c4762a1bSJed Brown for (i=sx; i<sx+nx; i++) { 272c4762a1bSJed Brown for (j=sy; j<sy+ny; j++) { 273c4762a1bSJed Brown 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; 274c4762a1bSJed Brown } 275c4762a1bSJed Brown } 276c4762a1bSJed Brown 2775f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,field,&FF)); 2785f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(cda,Gcoords,&XX)); 279c4762a1bSJed Brown PetscFunctionReturn(0); 280c4762a1bSJed Brown } 281c4762a1bSJed Brown 282c4762a1bSJed Brown PetscErrorCode DADefineXLinearField3D(DM da,Vec field) 283c4762a1bSJed Brown { 284c4762a1bSJed Brown PetscInt i,j,k; 285c4762a1bSJed Brown PetscInt sx,nx,sy,ny,sz,nz; 286c4762a1bSJed Brown Vec Gcoords; 287c4762a1bSJed Brown DMDACoor3d ***XX; 288c4762a1bSJed Brown PetscScalar ***FF; 289c4762a1bSJed Brown DM cda; 290c4762a1bSJed Brown 291c4762a1bSJed Brown PetscFunctionBeginUser; 2925f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(da,&cda)); 2935f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinates(da,&Gcoords)); 294c4762a1bSJed Brown 2955f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(cda,Gcoords,&XX)); 2965f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,field,&FF)); 297c4762a1bSJed Brown 2985f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&sx,&sy,&sz,&nx,&ny,&nz)); 299c4762a1bSJed Brown 300c4762a1bSJed Brown for (k=sz; k<sz+nz; k++) { 301c4762a1bSJed Brown for (j=sy; j<sy+ny; j++) { 302c4762a1bSJed Brown for (i=sx; i<sx+nx; i++) { 303c4762a1bSJed Brown FF[k][j][i] = 10.0 304c4762a1bSJed Brown + 4.05 * XX[k][j][i].x 305c4762a1bSJed Brown + 5.50 * XX[k][j][i].y 306c4762a1bSJed Brown + 1.33 * XX[k][j][i].z 307c4762a1bSJed Brown + 2.03 * XX[k][j][i].x * XX[k][j][i].y 308c4762a1bSJed Brown + 0.03 * XX[k][j][i].x * XX[k][j][i].z 309c4762a1bSJed Brown + 0.83 * XX[k][j][i].y * XX[k][j][i].z 310c4762a1bSJed Brown + 3.79 * XX[k][j][i].x * XX[k][j][i].y * XX[k][j][i].z; 311c4762a1bSJed Brown } 312c4762a1bSJed Brown } 313c4762a1bSJed Brown } 314c4762a1bSJed Brown 3155f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,field,&FF)); 3165f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(cda,Gcoords,&XX)); 317c4762a1bSJed Brown PetscFunctionReturn(0); 318c4762a1bSJed Brown } 319c4762a1bSJed Brown 320c4762a1bSJed Brown PetscErrorCode da_test_RefineCoords1D(PetscInt mx) 321c4762a1bSJed Brown { 322c4762a1bSJed Brown DM dac,daf; 323c4762a1bSJed Brown PetscViewer vv; 324c4762a1bSJed Brown Vec ac,af; 325c4762a1bSJed Brown PetscInt Mx; 326c4762a1bSJed Brown Mat II,INTERP; 327c4762a1bSJed Brown Vec scale; 328c4762a1bSJed Brown PetscBool output = PETSC_FALSE; 329c4762a1bSJed Brown 330c4762a1bSJed Brown PetscFunctionBeginUser; 3315f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE,mx+1,1, /* 1 dof */1, /* stencil = 1 */NULL,&dac)); 3325f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(dac)); 3335f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(dac)); 334c4762a1bSJed Brown 3355f80ce2aSJacob Faibussowitsch CHKERRQ(DMRefine(dac,MPI_COMM_NULL,&daf)); 3365f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(daf,0,&Mx,0,0,0,0,0,0,0,0,0,0,0)); 337c4762a1bSJed Brown Mx--; 338c4762a1bSJed Brown 3395f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetUniformCoordinates(dac, -1.0,1.0, PETSC_DECIDE,PETSC_DECIDE, PETSC_DECIDE,PETSC_DECIDE)); 3405f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetUniformCoordinates(daf, -1.0,1.0, PETSC_DECIDE,PETSC_DECIDE, PETSC_DECIDE,PETSC_DECIDE)); 341c4762a1bSJed Brown 342c4762a1bSJed Brown { 343c4762a1bSJed Brown DM cdaf,cdac; 344c4762a1bSJed Brown Vec coordsc,coordsf; 345c4762a1bSJed Brown 3465f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(dac,&cdac)); 3475f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(daf,&cdaf)); 348c4762a1bSJed Brown 3495f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinates(dac,&coordsc)); 3505f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinates(daf,&coordsf)); 351c4762a1bSJed Brown 3525f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateInterpolation(cdac,cdaf,&II,&scale)); 3535f80ce2aSJacob Faibussowitsch CHKERRQ(MatInterpolate(II,coordsc,coordsf)); 3545f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&II)); 3555f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&scale)); 356c4762a1bSJed Brown } 357c4762a1bSJed Brown 3585f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateInterpolation(dac,daf,&INTERP,NULL)); 359c4762a1bSJed Brown 3605f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dac,&ac)); 3615f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(ac,66.99)); 362c4762a1bSJed Brown 3635f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(daf,&af)); 364c4762a1bSJed Brown 3655f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(INTERP,ac, af)); 366c4762a1bSJed Brown 367c4762a1bSJed Brown { 368c4762a1bSJed Brown Vec afexact; 369c4762a1bSJed Brown PetscReal nrm; 370c4762a1bSJed Brown PetscInt N; 371c4762a1bSJed Brown 3725f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(daf,&afexact)); 3735f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(afexact,66.99)); 3745f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(afexact,-1.0,af)); /* af <= af - afinterp */ 3755f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(afexact,NORM_2,&nrm)); 3765f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(afexact,&N)); 3775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%D=>%D, interp err = %1.4e\n",mx,Mx,(double)(nrm/PetscSqrtReal((PetscReal)N)))); 3785f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&afexact)); 379c4762a1bSJed Brown } 380c4762a1bSJed Brown 3815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-output",&output,NULL)); 382c4762a1bSJed Brown if (output) { 3835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "dac_1D.vtr", &vv)); 3845f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(ac, vv)); 3855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&vv)); 386c4762a1bSJed Brown 3875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "daf_1D.vtr", &vv)); 3885f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(af, vv)); 3895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&vv)); 390c4762a1bSJed Brown } 391c4762a1bSJed Brown 3925f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&INTERP)); 3935f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dac)); 3945f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&daf)); 3955f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ac)); 3965f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&af)); 397c4762a1bSJed Brown PetscFunctionReturn(0); 398c4762a1bSJed Brown } 399c4762a1bSJed Brown 400c4762a1bSJed Brown PetscErrorCode da_test_RefineCoords2D(PetscInt mx,PetscInt my) 401c4762a1bSJed Brown { 402c4762a1bSJed Brown DM dac,daf; 403c4762a1bSJed Brown PetscViewer vv; 404c4762a1bSJed Brown Vec ac,af; 405c4762a1bSJed Brown PetscInt map_id,Mx,My; 406c4762a1bSJed Brown Mat II,INTERP; 407c4762a1bSJed Brown Vec scale; 408c4762a1bSJed Brown PetscBool output = PETSC_FALSE; 409c4762a1bSJed Brown 410c4762a1bSJed Brown PetscFunctionBeginUser; 4115f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 4125f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(dac)); 4135f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(dac)); 414c4762a1bSJed Brown 4155f80ce2aSJacob Faibussowitsch CHKERRQ(DMRefine(dac,MPI_COMM_NULL,&daf)); 4165f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(daf,0,&Mx,&My,0,0,0,0,0,0,0,0,0,0)); 417c4762a1bSJed Brown Mx--; My--; 418c4762a1bSJed Brown 4195f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetUniformCoordinates(dac, -1.0,1.0, -1.0,1.0, PETSC_DECIDE,PETSC_DECIDE)); 4205f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetUniformCoordinates(daf, -1.0,1.0, -1.0,1.0, PETSC_DECIDE,PETSC_DECIDE)); 421c4762a1bSJed Brown 422c4762a1bSJed Brown /* apply conformal mappings */ 423c4762a1bSJed Brown map_id = 0; 4245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-cmap", &map_id,NULL)); 425c4762a1bSJed Brown if (map_id >= 1) { 4265f80ce2aSJacob Faibussowitsch CHKERRQ(DAApplyConformalMapping(dac,map_id)); 427c4762a1bSJed Brown } 428c4762a1bSJed Brown 429c4762a1bSJed Brown { 430c4762a1bSJed Brown DM cdaf,cdac; 431c4762a1bSJed Brown Vec coordsc,coordsf; 432c4762a1bSJed Brown 4335f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(dac,&cdac)); 4345f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(daf,&cdaf)); 435c4762a1bSJed Brown 4365f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinates(dac,&coordsc)); 4375f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinates(daf,&coordsf)); 438c4762a1bSJed Brown 4395f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateInterpolation(cdac,cdaf,&II,&scale)); 4405f80ce2aSJacob Faibussowitsch CHKERRQ(MatInterpolate(II,coordsc,coordsf)); 4415f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&II)); 4425f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&scale)); 443c4762a1bSJed Brown } 444c4762a1bSJed Brown 4455f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateInterpolation(dac,daf,&INTERP,NULL)); 446c4762a1bSJed Brown 4475f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dac,&ac)); 4485f80ce2aSJacob Faibussowitsch CHKERRQ(DADefineXLinearField2D(dac,ac)); 449c4762a1bSJed Brown 4505f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(daf,&af)); 4515f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(INTERP,ac, af)); 452c4762a1bSJed Brown 453c4762a1bSJed Brown { 454c4762a1bSJed Brown Vec afexact; 455c4762a1bSJed Brown PetscReal nrm; 456c4762a1bSJed Brown PetscInt N; 457c4762a1bSJed Brown 4585f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(daf,&afexact)); 4595f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(afexact)); 4605f80ce2aSJacob Faibussowitsch CHKERRQ(DADefineXLinearField2D(daf,afexact)); 4615f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(afexact,-1.0,af)); /* af <= af - afinterp */ 4625f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(afexact,NORM_2,&nrm)); 4635f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(afexact,&N)); 4645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"[%D x %D]=>[%D x %D], interp err = %1.4e\n",mx,my,Mx,My,(double)(nrm/PetscSqrtReal((PetscReal)N)))); 4655f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&afexact)); 466c4762a1bSJed Brown } 467c4762a1bSJed Brown 4685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-output",&output,NULL)); 469c4762a1bSJed Brown if (output) { 4705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "dac_2D.vtr", &vv)); 4715f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(ac, vv)); 4725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&vv)); 473c4762a1bSJed Brown 4745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "daf_2D.vtr", &vv)); 4755f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(af, vv)); 4765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&vv)); 477c4762a1bSJed Brown } 478c4762a1bSJed Brown 4795f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&INTERP)); 4805f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dac)); 4815f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&daf)); 4825f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ac)); 4835f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&af)); 484c4762a1bSJed Brown PetscFunctionReturn(0); 485c4762a1bSJed Brown } 486c4762a1bSJed Brown 487c4762a1bSJed Brown PetscErrorCode da_test_RefineCoords3D(PetscInt mx,PetscInt my,PetscInt mz) 488c4762a1bSJed Brown { 489c4762a1bSJed Brown PetscErrorCode ierr; 490c4762a1bSJed Brown DM dac,daf; 491c4762a1bSJed Brown PetscViewer vv; 492c4762a1bSJed Brown Vec ac,af; 493c4762a1bSJed Brown PetscInt map_id,Mx,My,Mz; 494c4762a1bSJed Brown Mat II,INTERP; 495c4762a1bSJed Brown Vec scale; 496c4762a1bSJed Brown PetscBool output = PETSC_FALSE; 497c4762a1bSJed Brown 498c4762a1bSJed Brown PetscFunctionBeginUser; 499c4762a1bSJed Brown ierr = 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 */ 500c4762a1bSJed Brown 1, /* stencil = 1 */NULL,NULL,NULL,&dac);CHKERRQ(ierr); 5015f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(dac)); 5025f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(dac)); 503c4762a1bSJed Brown 5045f80ce2aSJacob Faibussowitsch CHKERRQ(DMRefine(dac,MPI_COMM_NULL,&daf)); 5055f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(daf,0,&Mx,&My,&Mz,0,0,0,0,0,0,0,0,0)); 506c4762a1bSJed Brown Mx--; My--; Mz--; 507c4762a1bSJed Brown 5085f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetUniformCoordinates(dac, -1.0,1.0, -1.0,1.0, -1.0,1.0)); 5095f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetUniformCoordinates(daf, -1.0,1.0, -1.0,1.0, -1.0,1.0)); 510c4762a1bSJed Brown 511c4762a1bSJed Brown /* apply trilinear mappings */ 5125f80ce2aSJacob Faibussowitsch /*CHKERRQ(DAApplyTrilinearMapping(dac));*/ 513c4762a1bSJed Brown /* apply conformal mappings */ 514c4762a1bSJed Brown map_id = 0; 5155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-cmap", &map_id,NULL)); 516c4762a1bSJed Brown if (map_id >= 1) { 5175f80ce2aSJacob Faibussowitsch CHKERRQ(DAApplyConformalMapping(dac,map_id)); 518c4762a1bSJed Brown } 519c4762a1bSJed Brown 520c4762a1bSJed Brown { 521c4762a1bSJed Brown DM cdaf,cdac; 522c4762a1bSJed Brown Vec coordsc,coordsf; 523c4762a1bSJed Brown 5245f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(dac,&cdac)); 5255f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(daf,&cdaf)); 526c4762a1bSJed Brown 5275f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinates(dac,&coordsc)); 5285f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinates(daf,&coordsf)); 529c4762a1bSJed Brown 5305f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateInterpolation(cdac,cdaf,&II,&scale)); 5315f80ce2aSJacob Faibussowitsch CHKERRQ(MatInterpolate(II,coordsc,coordsf)); 5325f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&II)); 5335f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&scale)); 534c4762a1bSJed Brown } 535c4762a1bSJed Brown 5365f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateInterpolation(dac,daf,&INTERP,NULL)); 537c4762a1bSJed Brown 5385f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dac,&ac)); 5395f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(ac)); 5405f80ce2aSJacob Faibussowitsch CHKERRQ(DADefineXLinearField3D(dac,ac)); 541c4762a1bSJed Brown 5425f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(daf,&af)); 5435f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(af)); 544c4762a1bSJed Brown 5455f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(INTERP,ac, af)); 546c4762a1bSJed Brown 547c4762a1bSJed Brown { 548c4762a1bSJed Brown Vec afexact; 549c4762a1bSJed Brown PetscReal nrm; 550c4762a1bSJed Brown PetscInt N; 551c4762a1bSJed Brown 5525f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(daf,&afexact)); 5535f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(afexact)); 5545f80ce2aSJacob Faibussowitsch CHKERRQ(DADefineXLinearField3D(daf,afexact)); 5555f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(afexact,-1.0,af)); /* af <= af - afinterp */ 5565f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(afexact,NORM_2,&nrm)); 5575f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(afexact,&N)); 5585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"[%D x %D x %D]=>[%D x %D x %D], interp err = %1.4e\n",mx,my,mz,Mx,My,Mz,(double)(nrm/PetscSqrtReal((PetscReal)N)))); 5595f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&afexact)); 560c4762a1bSJed Brown } 561c4762a1bSJed Brown 5625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-output",&output,NULL)); 563c4762a1bSJed Brown if (output) { 5645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "dac_3D.vtr", &vv)); 5655f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(ac, vv)); 5665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&vv)); 567c4762a1bSJed Brown 5685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "daf_3D.vtr", &vv)); 5695f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(af, vv)); 5705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&vv)); 571c4762a1bSJed Brown } 572c4762a1bSJed Brown 5735f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&INTERP)); 5745f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dac)); 5755f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&daf)); 5765f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ac)); 5775f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&af)); 578c4762a1bSJed Brown PetscFunctionReturn(0); 579c4762a1bSJed Brown } 580c4762a1bSJed Brown 581c4762a1bSJed Brown int main(int argc,char **argv) 582c4762a1bSJed Brown { 583c4762a1bSJed Brown PetscInt mx = 2,my = 2,mz = 2,l,nl,dim; 584c4762a1bSJed Brown 585*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,0,help)); 5865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mx", &mx, 0)); 5875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-my", &my, 0)); 5885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mz", &mz, 0)); 589c4762a1bSJed Brown nl = 1; 5905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nl", &nl, 0)); 591c4762a1bSJed Brown dim = 2; 5925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-dim", &dim, 0)); 593c4762a1bSJed Brown 594c4762a1bSJed Brown for (l=0; l<nl; l++) { 595b5675b0fSBarry Smith if (dim==1) { 5965f80ce2aSJacob Faibussowitsch CHKERRQ(da_test_RefineCoords1D(mx)); 597b5675b0fSBarry Smith } else if (dim==2) { 5985f80ce2aSJacob Faibussowitsch CHKERRQ(da_test_RefineCoords2D(mx,my)); 599b5675b0fSBarry Smith } else if (dim==3) { 6005f80ce2aSJacob Faibussowitsch CHKERRQ(da_test_RefineCoords3D(mx,my,mz)); 601b5675b0fSBarry Smith } 602c4762a1bSJed Brown mx = mx * 2; 603c4762a1bSJed Brown my = my * 2; 604c4762a1bSJed Brown mz = mz * 2; 605c4762a1bSJed Brown } 606*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 607*b122ec5aSJacob Faibussowitsch return 0; 608c4762a1bSJed Brown } 609c4762a1bSJed Brown 610c4762a1bSJed Brown /*TEST 611c4762a1bSJed Brown 612c4762a1bSJed Brown test: 613c4762a1bSJed Brown suffix: 1d 614c4762a1bSJed Brown args: -mx 10 -nl 6 -dim 1 615c4762a1bSJed Brown 616c4762a1bSJed Brown test: 617c4762a1bSJed Brown suffix: 2d 61863037964SStefano Zampini output_file: output/ex36_2d.out 61963037964SStefano Zampini args: -mx 10 -my 10 -nl 6 -dim 2 -cmap {{0 1 2 3}} 620c4762a1bSJed Brown 621c4762a1bSJed Brown test: 622c4762a1bSJed Brown suffix: 2dp1 623c4762a1bSJed Brown nsize: 8 624c4762a1bSJed Brown args: -mx 10 -my 10 -nl 4 -dim 2 -cmap 3 -da_refine_x 3 -da_refine_y 4 625c4762a1bSJed Brown timeoutfactor: 2 626c4762a1bSJed Brown 627c4762a1bSJed Brown test: 628c4762a1bSJed Brown suffix: 2dp2 629c4762a1bSJed Brown nsize: 8 630c4762a1bSJed Brown args: -mx 10 -my 10 -nl 4 -dim 2 -cmap 3 -da_refine_x 3 -da_refine_y 1 631c4762a1bSJed Brown timeoutfactor: 2 632c4762a1bSJed Brown 633c4762a1bSJed Brown test: 634c4762a1bSJed Brown suffix: 3d 635c4762a1bSJed Brown args: -mx 5 -my 5 -mz 5 -nl 4 -dim 3 -cmap 3 636c4762a1bSJed Brown 637c4762a1bSJed Brown test: 638c4762a1bSJed Brown suffix: 3dp1 639c4762a1bSJed Brown nsize: 32 640c4762a1bSJed 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 641c4762a1bSJed Brown 642c4762a1bSJed Brown TEST*/ 643