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 PetscErrorCode ierr; 59c4762a1bSJed Brown PetscInt i,n; 60c4762a1bSJed Brown PetscInt sx,nx,sy,ny,sz,nz,dim; 61c4762a1bSJed Brown Vec Gcoords; 62c4762a1bSJed Brown PetscScalar *XX; 63c4762a1bSJed Brown PetscScalar xx,yy,zz; 64c4762a1bSJed Brown DM cda; 65c4762a1bSJed Brown 66c4762a1bSJed Brown PetscFunctionBeginUser; 67c4762a1bSJed Brown if (idx==0) { 68c4762a1bSJed Brown PetscFunctionReturn(0); 69c4762a1bSJed Brown } else if (idx==1) { /* dam break */ 70c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(da, -1.0,1.0, -1.0,1.0, -1.0,1.0);CHKERRQ(ierr); 71c4762a1bSJed Brown } else if (idx==2) { /* stagnation in a corner */ 72c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(da, -1.0,1.0, 0.0,1.0, -1.0,1.0);CHKERRQ(ierr); 73c4762a1bSJed Brown } else if (idx==3) { /* nautilis */ 74c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(da, -1.0,1.0, -1.0,1.0, -1.0,1.0);CHKERRQ(ierr); 75c4762a1bSJed Brown } else if (idx==4) { 76c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(da, -1.0,1.0, -1.0,1.0, -1.0,1.0);CHKERRQ(ierr); 77c4762a1bSJed Brown } 78c4762a1bSJed Brown 79c4762a1bSJed Brown ierr = DMGetCoordinateDM(da,&cda);CHKERRQ(ierr); 80c4762a1bSJed Brown ierr = DMGetCoordinates(da,&Gcoords);CHKERRQ(ierr); 81c4762a1bSJed Brown 82c4762a1bSJed Brown ierr = VecGetArray(Gcoords,&XX);CHKERRQ(ierr); 83c4762a1bSJed Brown ierr = DMDAGetCorners(da,&sx,&sy,&sz,&nx,&ny,&nz);CHKERRQ(ierr); 84c4762a1bSJed Brown ierr = DMDAGetInfo(da, &dim, 0,0,0, 0,0,0, 0,0,0,0,0,0);CHKERRQ(ierr); 85c4762a1bSJed Brown ierr = VecGetLocalSize(Gcoords,&n);CHKERRQ(ierr); 86c4762a1bSJed Brown n = n / dim; 87c4762a1bSJed Brown 88c4762a1bSJed Brown for (i=0; i<n; i++) { 89c4762a1bSJed Brown if ((dim==3) && (idx!=2)) { 90c4762a1bSJed Brown PetscScalar Ni[8]; 91c4762a1bSJed Brown PetscScalar xi = XX[dim*i]; 92c4762a1bSJed Brown PetscScalar eta = XX[dim*i+1]; 93c4762a1bSJed Brown PetscScalar zeta = XX[dim*i+2]; 94c4762a1bSJed Brown PetscScalar xn[] = {-1.0,1.0,-1.0,1.0, -1.0,1.0,-1.0,1.0 }; 95c4762a1bSJed Brown PetscScalar yn[] = {-1.0,-1.0,1.0,1.0, -1.0,-1.0,1.0,1.0 }; 96c4762a1bSJed Brown PetscScalar zn[] = {-0.1,-4.0,-0.2,-1.0, 0.1,4.0,0.2,1.0 }; 97c4762a1bSJed Brown PetscInt p; 98c4762a1bSJed Brown 99c4762a1bSJed Brown Ni[0] = 0.125*(1.0-xi)*(1.0-eta)*(1.0-zeta); 100c4762a1bSJed Brown Ni[1] = 0.125*(1.0+xi)*(1.0-eta)*(1.0-zeta); 101c4762a1bSJed Brown Ni[2] = 0.125*(1.0-xi)*(1.0+eta)*(1.0-zeta); 102c4762a1bSJed Brown Ni[3] = 0.125*(1.0+xi)*(1.0+eta)*(1.0-zeta); 103c4762a1bSJed Brown 104c4762a1bSJed Brown Ni[4] = 0.125*(1.0-xi)*(1.0-eta)*(1.0+zeta); 105c4762a1bSJed Brown Ni[5] = 0.125*(1.0+xi)*(1.0-eta)*(1.0+zeta); 106c4762a1bSJed Brown Ni[6] = 0.125*(1.0-xi)*(1.0+eta)*(1.0+zeta); 107c4762a1bSJed Brown Ni[7] = 0.125*(1.0+xi)*(1.0+eta)*(1.0+zeta); 108c4762a1bSJed Brown 109c4762a1bSJed Brown xx = yy = zz = 0.0; 110c4762a1bSJed Brown for (p=0; p<8; p++) { 111c4762a1bSJed Brown xx += Ni[p]*xn[p]; 112c4762a1bSJed Brown yy += Ni[p]*yn[p]; 113c4762a1bSJed Brown zz += Ni[p]*zn[p]; 114c4762a1bSJed Brown } 115c4762a1bSJed Brown XX[dim*i] = xx; 116c4762a1bSJed Brown XX[dim*i+1] = yy; 117c4762a1bSJed Brown XX[dim*i+2] = zz; 118c4762a1bSJed Brown } 119c4762a1bSJed Brown 120c4762a1bSJed Brown if (idx==1) { 121c4762a1bSJed Brown CCmplx zeta,t1,t2; 122c4762a1bSJed Brown 123c4762a1bSJed Brown xx = XX[dim*i] - 0.8; 124c4762a1bSJed Brown yy = XX[dim*i+1] + 1.5; 125c4762a1bSJed Brown 126c4762a1bSJed Brown zeta.real = PetscRealPart(xx); 127c4762a1bSJed Brown zeta.imag = PetscRealPart(yy); 128c4762a1bSJed Brown 129c4762a1bSJed Brown t1 = CCmplxPow(zeta,-1.0); 130c4762a1bSJed Brown t2 = CCmplxAdd(zeta,t1); 131c4762a1bSJed Brown 132c4762a1bSJed Brown XX[dim*i] = CCmplxRe(t2); 133c4762a1bSJed Brown XX[dim*i+1] = CCmplxIm(t2); 134c4762a1bSJed Brown } else if (idx==2) { 135c4762a1bSJed Brown CCmplx zeta,t1; 136c4762a1bSJed Brown 137c4762a1bSJed Brown xx = XX[dim*i]; 138c4762a1bSJed Brown yy = XX[dim*i+1]; 139c4762a1bSJed Brown zeta.real = PetscRealPart(xx); 140c4762a1bSJed Brown zeta.imag = PetscRealPart(yy); 141c4762a1bSJed Brown 142c4762a1bSJed Brown t1 = CCmplxSqrt(zeta); 143c4762a1bSJed Brown XX[dim*i] = CCmplxRe(t1); 144c4762a1bSJed Brown XX[dim*i+1] = CCmplxIm(t1); 145c4762a1bSJed Brown } else if (idx==3) { 146c4762a1bSJed Brown CCmplx zeta,t1,t2; 147c4762a1bSJed Brown 148c4762a1bSJed Brown xx = XX[dim*i] - 0.8; 149c4762a1bSJed Brown yy = XX[dim*i+1] + 1.5; 150c4762a1bSJed Brown 151c4762a1bSJed Brown zeta.real = PetscRealPart(xx); 152c4762a1bSJed Brown zeta.imag = PetscRealPart(yy); 153c4762a1bSJed Brown t1 = CCmplxPow(zeta,-1.0); 154c4762a1bSJed Brown t2 = CCmplxAdd(zeta,t1); 155c4762a1bSJed Brown XX[dim*i] = CCmplxRe(t2); 156c4762a1bSJed Brown XX[dim*i+1] = CCmplxIm(t2); 157c4762a1bSJed Brown 158c4762a1bSJed Brown xx = XX[dim*i]; 159c4762a1bSJed Brown yy = XX[dim*i+1]; 160c4762a1bSJed Brown zeta.real = PetscRealPart(xx); 161c4762a1bSJed Brown zeta.imag = PetscRealPart(yy); 162c4762a1bSJed Brown t1 = CCmplxExp(zeta); 163c4762a1bSJed Brown XX[dim*i] = CCmplxRe(t1); 164c4762a1bSJed Brown XX[dim*i+1] = CCmplxIm(t1); 165c4762a1bSJed Brown 166c4762a1bSJed Brown xx = XX[dim*i] + 0.4; 167c4762a1bSJed Brown yy = XX[dim*i+1]; 168c4762a1bSJed Brown zeta.real = PetscRealPart(xx); 169c4762a1bSJed Brown zeta.imag = PetscRealPart(yy); 170c4762a1bSJed Brown t1 = CCmplxPow(zeta,2.0); 171c4762a1bSJed Brown XX[dim*i] = CCmplxRe(t1); 172c4762a1bSJed Brown XX[dim*i+1] = CCmplxIm(t1); 173c4762a1bSJed Brown } else if (idx==4) { 174c4762a1bSJed Brown PetscScalar Ni[4]; 175c4762a1bSJed Brown PetscScalar xi = XX[dim*i]; 176c4762a1bSJed Brown PetscScalar eta = XX[dim*i+1]; 177c4762a1bSJed Brown PetscScalar xn[] = {0.0,2.0,0.2,3.5}; 178c4762a1bSJed Brown PetscScalar yn[] = {-1.3,0.0,2.0,4.0}; 179c4762a1bSJed Brown PetscInt p; 180c4762a1bSJed Brown 181c4762a1bSJed Brown Ni[0] = 0.25*(1.0-xi)*(1.0-eta); 182c4762a1bSJed Brown Ni[1] = 0.25*(1.0+xi)*(1.0-eta); 183c4762a1bSJed Brown Ni[2] = 0.25*(1.0-xi)*(1.0+eta); 184c4762a1bSJed Brown Ni[3] = 0.25*(1.0+xi)*(1.0+eta); 185c4762a1bSJed Brown 186c4762a1bSJed Brown xx = yy = 0.0; 187c4762a1bSJed Brown for (p=0; p<4; p++) { 188c4762a1bSJed Brown xx += Ni[p]*xn[p]; 189c4762a1bSJed Brown yy += Ni[p]*yn[p]; 190c4762a1bSJed Brown } 191c4762a1bSJed Brown XX[dim*i] = xx; 192c4762a1bSJed Brown XX[dim*i+1] = yy; 193c4762a1bSJed Brown } 194c4762a1bSJed Brown } 195c4762a1bSJed Brown ierr = VecRestoreArray(Gcoords,&XX);CHKERRQ(ierr); 196c4762a1bSJed Brown PetscFunctionReturn(0); 197c4762a1bSJed Brown } 198c4762a1bSJed Brown 199c4762a1bSJed Brown PetscErrorCode DAApplyTrilinearMapping(DM da) 200c4762a1bSJed Brown { 201c4762a1bSJed Brown PetscErrorCode ierr; 202c4762a1bSJed Brown PetscInt i,j,k; 203c4762a1bSJed Brown PetscInt sx,nx,sy,ny,sz,nz; 204c4762a1bSJed Brown Vec Gcoords; 205c4762a1bSJed Brown DMDACoor3d ***XX; 206c4762a1bSJed Brown PetscScalar xx,yy,zz; 207c4762a1bSJed Brown DM cda; 208c4762a1bSJed Brown 209c4762a1bSJed Brown PetscFunctionBeginUser; 210c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(da, -1.0,1.0, -1.0,1.0, -1.0,1.0);CHKERRQ(ierr); 211c4762a1bSJed Brown ierr = DMGetCoordinateDM(da,&cda);CHKERRQ(ierr); 212c4762a1bSJed Brown ierr = DMGetCoordinates(da,&Gcoords);CHKERRQ(ierr); 213c4762a1bSJed Brown 214c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(cda,Gcoords,&XX);CHKERRQ(ierr); 215c4762a1bSJed Brown ierr = DMDAGetCorners(da,&sx,&sy,&sz,&nx,&ny,&nz);CHKERRQ(ierr); 216c4762a1bSJed Brown 217c4762a1bSJed Brown for (i=sx; i<sx+nx; i++) { 218c4762a1bSJed Brown for (j=sy; j<sy+ny; j++) { 219c4762a1bSJed Brown for (k=sz; k<sz+nz; k++) { 220c4762a1bSJed Brown PetscScalar Ni[8]; 221c4762a1bSJed Brown PetscScalar xi = XX[k][j][i].x; 222c4762a1bSJed Brown PetscScalar eta = XX[k][j][i].y; 223c4762a1bSJed Brown PetscScalar zeta = XX[k][j][i].z; 224c4762a1bSJed Brown PetscScalar xn[] = {0.0,2.0,0.2,3.5, 0.0,2.1,0.23,3.125 }; 225c4762a1bSJed Brown PetscScalar yn[] = {-1.3,0.0,2.0,4.0, -1.45,-0.1,2.24,3.79 }; 226c4762a1bSJed Brown PetscScalar zn[] = {0.0,0.3,-0.1,0.123, 0.956,1.32,1.12,0.798 }; 227c4762a1bSJed Brown PetscInt p; 228c4762a1bSJed Brown 229c4762a1bSJed Brown Ni[0] = 0.125*(1.0-xi)*(1.0-eta)*(1.0-zeta); 230c4762a1bSJed Brown Ni[1] = 0.125*(1.0+xi)*(1.0-eta)*(1.0-zeta); 231c4762a1bSJed Brown Ni[2] = 0.125*(1.0-xi)*(1.0+eta)*(1.0-zeta); 232c4762a1bSJed Brown Ni[3] = 0.125*(1.0+xi)*(1.0+eta)*(1.0-zeta); 233c4762a1bSJed Brown 234c4762a1bSJed Brown Ni[4] = 0.125*(1.0-xi)*(1.0-eta)*(1.0+zeta); 235c4762a1bSJed Brown Ni[5] = 0.125*(1.0+xi)*(1.0-eta)*(1.0+zeta); 236c4762a1bSJed Brown Ni[6] = 0.125*(1.0-xi)*(1.0+eta)*(1.0+zeta); 237c4762a1bSJed Brown Ni[7] = 0.125*(1.0+xi)*(1.0+eta)*(1.0+zeta); 238c4762a1bSJed Brown 239c4762a1bSJed Brown xx = yy = zz = 0.0; 240c4762a1bSJed Brown for (p=0; p<8; p++) { 241c4762a1bSJed Brown xx += Ni[p]*xn[p]; 242c4762a1bSJed Brown yy += Ni[p]*yn[p]; 243c4762a1bSJed Brown zz += Ni[p]*zn[p]; 244c4762a1bSJed Brown } 245c4762a1bSJed Brown XX[k][j][i].x = xx; 246c4762a1bSJed Brown XX[k][j][i].y = yy; 247c4762a1bSJed Brown XX[k][j][i].z = zz; 248c4762a1bSJed Brown } 249c4762a1bSJed Brown } 250c4762a1bSJed Brown } 251c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(cda,Gcoords,&XX);CHKERRQ(ierr); 252c4762a1bSJed Brown PetscFunctionReturn(0); 253c4762a1bSJed Brown } 254c4762a1bSJed Brown 255c4762a1bSJed Brown PetscErrorCode DADefineXLinearField2D(DM da,Vec field) 256c4762a1bSJed Brown { 257c4762a1bSJed Brown PetscErrorCode ierr; 258c4762a1bSJed Brown PetscInt i,j; 259c4762a1bSJed Brown PetscInt sx,nx,sy,ny; 260c4762a1bSJed Brown Vec Gcoords; 261c4762a1bSJed Brown DMDACoor2d **XX; 262c4762a1bSJed Brown PetscScalar **FF; 263c4762a1bSJed Brown DM cda; 264c4762a1bSJed Brown 265c4762a1bSJed Brown PetscFunctionBeginUser; 266c4762a1bSJed Brown ierr = DMGetCoordinateDM(da,&cda);CHKERRQ(ierr); 267c4762a1bSJed Brown ierr = DMGetCoordinates(da,&Gcoords);CHKERRQ(ierr); 268c4762a1bSJed Brown 269c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(cda,Gcoords,&XX);CHKERRQ(ierr); 270c4762a1bSJed Brown ierr = DMDAVecGetArray(da,field,&FF);CHKERRQ(ierr); 271c4762a1bSJed Brown 272c4762a1bSJed Brown ierr = DMDAGetCorners(da,&sx,&sy,0,&nx,&ny,0);CHKERRQ(ierr); 273c4762a1bSJed Brown 274c4762a1bSJed Brown for (i=sx; i<sx+nx; i++) { 275c4762a1bSJed Brown for (j=sy; j<sy+ny; j++) { 276c4762a1bSJed 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; 277c4762a1bSJed Brown } 278c4762a1bSJed Brown } 279c4762a1bSJed Brown 280c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,field,&FF);CHKERRQ(ierr); 281c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(cda,Gcoords,&XX);CHKERRQ(ierr); 282c4762a1bSJed Brown PetscFunctionReturn(0); 283c4762a1bSJed Brown } 284c4762a1bSJed Brown 285c4762a1bSJed Brown PetscErrorCode DADefineXLinearField3D(DM da,Vec field) 286c4762a1bSJed Brown { 287c4762a1bSJed Brown PetscErrorCode ierr; 288c4762a1bSJed Brown PetscInt i,j,k; 289c4762a1bSJed Brown PetscInt sx,nx,sy,ny,sz,nz; 290c4762a1bSJed Brown Vec Gcoords; 291c4762a1bSJed Brown DMDACoor3d ***XX; 292c4762a1bSJed Brown PetscScalar ***FF; 293c4762a1bSJed Brown DM cda; 294c4762a1bSJed Brown 295c4762a1bSJed Brown PetscFunctionBeginUser; 296c4762a1bSJed Brown ierr = DMGetCoordinateDM(da,&cda);CHKERRQ(ierr); 297c4762a1bSJed Brown ierr = DMGetCoordinates(da,&Gcoords);CHKERRQ(ierr); 298c4762a1bSJed Brown 299c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(cda,Gcoords,&XX);CHKERRQ(ierr); 300c4762a1bSJed Brown ierr = DMDAVecGetArray(da,field,&FF);CHKERRQ(ierr); 301c4762a1bSJed Brown 302c4762a1bSJed Brown ierr = DMDAGetCorners(da,&sx,&sy,&sz,&nx,&ny,&nz);CHKERRQ(ierr); 303c4762a1bSJed Brown 304c4762a1bSJed Brown for (k=sz; k<sz+nz; k++) { 305c4762a1bSJed Brown for (j=sy; j<sy+ny; j++) { 306c4762a1bSJed Brown for (i=sx; i<sx+nx; i++) { 307c4762a1bSJed Brown FF[k][j][i] = 10.0 308c4762a1bSJed Brown + 4.05 * XX[k][j][i].x 309c4762a1bSJed Brown + 5.50 * XX[k][j][i].y 310c4762a1bSJed Brown + 1.33 * XX[k][j][i].z 311c4762a1bSJed Brown + 2.03 * XX[k][j][i].x * XX[k][j][i].y 312c4762a1bSJed Brown + 0.03 * XX[k][j][i].x * XX[k][j][i].z 313c4762a1bSJed Brown + 0.83 * XX[k][j][i].y * XX[k][j][i].z 314c4762a1bSJed Brown + 3.79 * XX[k][j][i].x * XX[k][j][i].y * XX[k][j][i].z; 315c4762a1bSJed Brown } 316c4762a1bSJed Brown } 317c4762a1bSJed Brown } 318c4762a1bSJed Brown 319c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,field,&FF);CHKERRQ(ierr); 320c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(cda,Gcoords,&XX);CHKERRQ(ierr); 321c4762a1bSJed Brown PetscFunctionReturn(0); 322c4762a1bSJed Brown } 323c4762a1bSJed Brown 324c4762a1bSJed Brown PetscErrorCode da_test_RefineCoords1D(PetscInt mx) 325c4762a1bSJed Brown { 326c4762a1bSJed Brown PetscErrorCode ierr; 327c4762a1bSJed Brown DM dac,daf; 328c4762a1bSJed Brown PetscViewer vv; 329c4762a1bSJed Brown Vec ac,af; 330c4762a1bSJed Brown PetscInt Mx; 331c4762a1bSJed Brown Mat II,INTERP; 332c4762a1bSJed Brown Vec scale; 333c4762a1bSJed Brown PetscBool output = PETSC_FALSE; 334c4762a1bSJed Brown 335c4762a1bSJed Brown PetscFunctionBeginUser; 336c4762a1bSJed Brown ierr = DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE,mx+1,1, /* 1 dof */1, /* stencil = 1 */NULL,&dac);CHKERRQ(ierr); 337c4762a1bSJed Brown ierr = DMSetFromOptions(dac);CHKERRQ(ierr); 338c4762a1bSJed Brown ierr = DMSetUp(dac);CHKERRQ(ierr); 339c4762a1bSJed Brown 340c4762a1bSJed Brown ierr = DMRefine(dac,MPI_COMM_NULL,&daf);CHKERRQ(ierr); 341c4762a1bSJed Brown ierr = DMDAGetInfo(daf,0,&Mx,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 342c4762a1bSJed Brown Mx--; 343c4762a1bSJed Brown 344c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(dac, -1.0,1.0, PETSC_DECIDE,PETSC_DECIDE, PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 345c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(daf, -1.0,1.0, PETSC_DECIDE,PETSC_DECIDE, PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 346c4762a1bSJed Brown 347c4762a1bSJed Brown { 348c4762a1bSJed Brown DM cdaf,cdac; 349c4762a1bSJed Brown Vec coordsc,coordsf; 350c4762a1bSJed Brown 351c4762a1bSJed Brown ierr = DMGetCoordinateDM(dac,&cdac);CHKERRQ(ierr); 352c4762a1bSJed Brown ierr = DMGetCoordinateDM(daf,&cdaf);CHKERRQ(ierr); 353c4762a1bSJed Brown 354c4762a1bSJed Brown ierr = DMGetCoordinates(dac,&coordsc);CHKERRQ(ierr); 355c4762a1bSJed Brown ierr = DMGetCoordinates(daf,&coordsf);CHKERRQ(ierr); 356c4762a1bSJed Brown 357c4762a1bSJed Brown ierr = DMCreateInterpolation(cdac,cdaf,&II,&scale);CHKERRQ(ierr); 358c4762a1bSJed Brown ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr); 359c4762a1bSJed Brown ierr = MatDestroy(&II);CHKERRQ(ierr); 360c4762a1bSJed Brown ierr = VecDestroy(&scale);CHKERRQ(ierr); 361c4762a1bSJed Brown } 362c4762a1bSJed Brown 363c4762a1bSJed Brown ierr = DMCreateInterpolation(dac,daf,&INTERP,NULL);CHKERRQ(ierr); 364c4762a1bSJed Brown 365c4762a1bSJed Brown ierr = DMCreateGlobalVector(dac,&ac);CHKERRQ(ierr); 366c4762a1bSJed Brown ierr = VecSet(ac,66.99);CHKERRQ(ierr); 367c4762a1bSJed Brown 368c4762a1bSJed Brown ierr = DMCreateGlobalVector(daf,&af);CHKERRQ(ierr); 369c4762a1bSJed Brown 370c4762a1bSJed Brown ierr = MatMult(INTERP,ac, af);CHKERRQ(ierr); 371c4762a1bSJed Brown 372c4762a1bSJed Brown { 373c4762a1bSJed Brown Vec afexact; 374c4762a1bSJed Brown PetscReal nrm; 375c4762a1bSJed Brown PetscInt N; 376c4762a1bSJed Brown 377c4762a1bSJed Brown ierr = DMCreateGlobalVector(daf,&afexact);CHKERRQ(ierr); 378c4762a1bSJed Brown ierr = VecSet(afexact,66.99);CHKERRQ(ierr); 379c4762a1bSJed Brown ierr = VecAXPY(afexact,-1.0,af);CHKERRQ(ierr); /* af <= af - afinterp */ 380c4762a1bSJed Brown ierr = VecNorm(afexact,NORM_2,&nrm);CHKERRQ(ierr); 381c4762a1bSJed Brown ierr = VecGetSize(afexact,&N);CHKERRQ(ierr); 382*b5675b0fSBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"%D=>%D, interp err = %1.4e\n",mx,Mx,(double)(nrm/PetscSqrtReal((PetscReal)N)));CHKERRQ(ierr); 383c4762a1bSJed Brown ierr = VecDestroy(&afexact);CHKERRQ(ierr); 384c4762a1bSJed Brown } 385c4762a1bSJed Brown 386*b5675b0fSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-output",&output,NULL);CHKERRQ(ierr); 387c4762a1bSJed Brown if (output) { 388c4762a1bSJed Brown ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD, "dac_1D.vtk", &vv);CHKERRQ(ierr); 389c4762a1bSJed Brown ierr = PetscViewerPushFormat(vv, PETSC_VIEWER_ASCII_VTK);CHKERRQ(ierr); 390c4762a1bSJed Brown ierr = DMView(dac, vv);CHKERRQ(ierr); 391c4762a1bSJed Brown ierr = VecView(ac, vv);CHKERRQ(ierr); 392c4762a1bSJed Brown ierr = PetscViewerPopFormat(vv);CHKERRQ(ierr); 393c4762a1bSJed Brown ierr = PetscViewerDestroy(&vv);CHKERRQ(ierr); 394c4762a1bSJed Brown 395c4762a1bSJed Brown ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD, "daf_1D.vtk", &vv);CHKERRQ(ierr); 396c4762a1bSJed Brown ierr = PetscViewerPushFormat(vv, PETSC_VIEWER_ASCII_VTK);CHKERRQ(ierr); 397c4762a1bSJed Brown ierr = DMView(daf, vv);CHKERRQ(ierr); 398c4762a1bSJed Brown ierr = VecView(af, vv);CHKERRQ(ierr); 399c4762a1bSJed Brown ierr = PetscViewerPopFormat(vv);CHKERRQ(ierr); 400c4762a1bSJed Brown ierr = PetscViewerDestroy(&vv);CHKERRQ(ierr); 401c4762a1bSJed Brown } 402c4762a1bSJed Brown 403c4762a1bSJed Brown ierr = MatDestroy(&INTERP);CHKERRQ(ierr); 404c4762a1bSJed Brown ierr = DMDestroy(&dac);CHKERRQ(ierr); 405c4762a1bSJed Brown ierr = DMDestroy(&daf);CHKERRQ(ierr); 406c4762a1bSJed Brown ierr = VecDestroy(&ac);CHKERRQ(ierr); 407c4762a1bSJed Brown ierr = VecDestroy(&af);CHKERRQ(ierr); 408c4762a1bSJed Brown PetscFunctionReturn(0); 409c4762a1bSJed Brown } 410c4762a1bSJed Brown 411c4762a1bSJed Brown PetscErrorCode da_test_RefineCoords2D(PetscInt mx,PetscInt my) 412c4762a1bSJed Brown { 413c4762a1bSJed Brown PetscErrorCode ierr; 414c4762a1bSJed Brown DM dac,daf; 415c4762a1bSJed Brown PetscViewer vv; 416c4762a1bSJed Brown Vec ac,af; 417c4762a1bSJed Brown PetscInt map_id,Mx,My; 418c4762a1bSJed Brown Mat II,INTERP; 419c4762a1bSJed Brown Vec scale; 420c4762a1bSJed Brown PetscBool output = PETSC_FALSE; 421c4762a1bSJed Brown 422c4762a1bSJed Brown PetscFunctionBeginUser; 423c4762a1bSJed Brown ierr = 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);CHKERRQ(ierr); 424c4762a1bSJed Brown ierr = DMSetFromOptions(dac);CHKERRQ(ierr); 425c4762a1bSJed Brown ierr = DMSetUp(dac);CHKERRQ(ierr); 426c4762a1bSJed Brown 427c4762a1bSJed Brown ierr = DMRefine(dac,MPI_COMM_NULL,&daf);CHKERRQ(ierr); 428c4762a1bSJed Brown ierr = DMDAGetInfo(daf,0,&Mx,&My,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 429c4762a1bSJed Brown Mx--; My--; 430c4762a1bSJed Brown 431c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(dac, -1.0,1.0, -1.0,1.0, PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 432c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(daf, -1.0,1.0, -1.0,1.0, PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 433c4762a1bSJed Brown 434c4762a1bSJed Brown /* apply conformal mappings */ 435c4762a1bSJed Brown map_id = 0; 436c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-cmap", &map_id,NULL);CHKERRQ(ierr); 437c4762a1bSJed Brown if (map_id >= 1) { 438c4762a1bSJed Brown ierr = DAApplyConformalMapping(dac,map_id);CHKERRQ(ierr); 439c4762a1bSJed Brown } 440c4762a1bSJed Brown 441c4762a1bSJed Brown { 442c4762a1bSJed Brown DM cdaf,cdac; 443c4762a1bSJed Brown Vec coordsc,coordsf; 444c4762a1bSJed Brown 445c4762a1bSJed Brown ierr = DMGetCoordinateDM(dac,&cdac);CHKERRQ(ierr); 446c4762a1bSJed Brown ierr = DMGetCoordinateDM(daf,&cdaf);CHKERRQ(ierr); 447c4762a1bSJed Brown 448c4762a1bSJed Brown ierr = DMGetCoordinates(dac,&coordsc);CHKERRQ(ierr); 449c4762a1bSJed Brown ierr = DMGetCoordinates(daf,&coordsf);CHKERRQ(ierr); 450c4762a1bSJed Brown 451c4762a1bSJed Brown ierr = DMCreateInterpolation(cdac,cdaf,&II,&scale);CHKERRQ(ierr); 452c4762a1bSJed Brown ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr); 453c4762a1bSJed Brown ierr = MatDestroy(&II);CHKERRQ(ierr); 454c4762a1bSJed Brown ierr = VecDestroy(&scale);CHKERRQ(ierr); 455c4762a1bSJed Brown } 456c4762a1bSJed Brown 457c4762a1bSJed Brown 458c4762a1bSJed Brown ierr = DMCreateInterpolation(dac,daf,&INTERP,NULL);CHKERRQ(ierr); 459c4762a1bSJed Brown 460c4762a1bSJed Brown ierr = DMCreateGlobalVector(dac,&ac);CHKERRQ(ierr); 461c4762a1bSJed Brown ierr = DADefineXLinearField2D(dac,ac);CHKERRQ(ierr); 462c4762a1bSJed Brown 463c4762a1bSJed Brown ierr = DMCreateGlobalVector(daf,&af);CHKERRQ(ierr); 464c4762a1bSJed Brown ierr = MatMult(INTERP,ac, af);CHKERRQ(ierr); 465c4762a1bSJed Brown 466c4762a1bSJed Brown { 467c4762a1bSJed Brown Vec afexact; 468c4762a1bSJed Brown PetscReal nrm; 469c4762a1bSJed Brown PetscInt N; 470c4762a1bSJed Brown 471c4762a1bSJed Brown ierr = DMCreateGlobalVector(daf,&afexact);CHKERRQ(ierr); 472c4762a1bSJed Brown ierr = VecZeroEntries(afexact);CHKERRQ(ierr); 473c4762a1bSJed Brown ierr = DADefineXLinearField2D(daf,afexact);CHKERRQ(ierr); 474c4762a1bSJed Brown ierr = VecAXPY(afexact,-1.0,af);CHKERRQ(ierr); /* af <= af - afinterp */ 475c4762a1bSJed Brown ierr = VecNorm(afexact,NORM_2,&nrm);CHKERRQ(ierr); 476c4762a1bSJed Brown ierr = VecGetSize(afexact,&N);CHKERRQ(ierr); 477*b5675b0fSBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"[%D x %D]=>[%D x %D], interp err = %1.4e\n",mx,my,Mx,My,(double)(nrm/PetscSqrtReal((PetscReal)N)));CHKERRQ(ierr); 478c4762a1bSJed Brown ierr = VecDestroy(&afexact);CHKERRQ(ierr); 479c4762a1bSJed Brown } 480c4762a1bSJed Brown 481*b5675b0fSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-output",&output,NULL);CHKERRQ(ierr); 482c4762a1bSJed Brown if (output) { 483c4762a1bSJed Brown ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD, "dac_2D.vtk", &vv);CHKERRQ(ierr); 484c4762a1bSJed Brown ierr = PetscViewerPushFormat(vv, PETSC_VIEWER_ASCII_VTK);CHKERRQ(ierr); 485c4762a1bSJed Brown ierr = DMView(dac, vv);CHKERRQ(ierr); 486c4762a1bSJed Brown ierr = VecView(ac, vv);CHKERRQ(ierr); 487c4762a1bSJed Brown ierr = PetscViewerPopFormat(vv);CHKERRQ(ierr); 488c4762a1bSJed Brown ierr = PetscViewerDestroy(&vv);CHKERRQ(ierr); 489c4762a1bSJed Brown 490c4762a1bSJed Brown ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD, "daf_2D.vtk", &vv);CHKERRQ(ierr); 491c4762a1bSJed Brown ierr = PetscViewerPushFormat(vv, PETSC_VIEWER_ASCII_VTK);CHKERRQ(ierr); 492c4762a1bSJed Brown ierr = DMView(daf, vv);CHKERRQ(ierr); 493c4762a1bSJed Brown ierr = VecView(af, vv);CHKERRQ(ierr); 494c4762a1bSJed Brown ierr = PetscViewerPopFormat(vv);CHKERRQ(ierr); 495c4762a1bSJed Brown ierr = PetscViewerDestroy(&vv);CHKERRQ(ierr); 496c4762a1bSJed Brown } 497c4762a1bSJed Brown 498c4762a1bSJed Brown ierr = MatDestroy(&INTERP);CHKERRQ(ierr); 499c4762a1bSJed Brown ierr = DMDestroy(&dac);CHKERRQ(ierr); 500c4762a1bSJed Brown ierr = DMDestroy(&daf);CHKERRQ(ierr); 501c4762a1bSJed Brown ierr = VecDestroy(&ac);CHKERRQ(ierr); 502c4762a1bSJed Brown ierr = VecDestroy(&af);CHKERRQ(ierr); 503c4762a1bSJed Brown PetscFunctionReturn(0); 504c4762a1bSJed Brown } 505c4762a1bSJed Brown 506c4762a1bSJed Brown PetscErrorCode da_test_RefineCoords3D(PetscInt mx,PetscInt my,PetscInt mz) 507c4762a1bSJed Brown { 508c4762a1bSJed Brown PetscErrorCode ierr; 509c4762a1bSJed Brown DM dac,daf; 510c4762a1bSJed Brown PetscViewer vv; 511c4762a1bSJed Brown Vec ac,af; 512c4762a1bSJed Brown PetscInt map_id,Mx,My,Mz; 513c4762a1bSJed Brown Mat II,INTERP; 514c4762a1bSJed Brown Vec scale; 515c4762a1bSJed Brown PetscBool output = PETSC_FALSE; 516c4762a1bSJed Brown 517c4762a1bSJed Brown PetscFunctionBeginUser; 518c4762a1bSJed 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 */ 519c4762a1bSJed Brown 1, /* stencil = 1 */NULL,NULL,NULL,&dac);CHKERRQ(ierr); 520c4762a1bSJed Brown ierr = DMSetFromOptions(dac);CHKERRQ(ierr); 521c4762a1bSJed Brown ierr = DMSetUp(dac);CHKERRQ(ierr); 522c4762a1bSJed Brown 523c4762a1bSJed Brown ierr = DMRefine(dac,MPI_COMM_NULL,&daf);CHKERRQ(ierr); 524c4762a1bSJed Brown ierr = DMDAGetInfo(daf,0,&Mx,&My,&Mz,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 525c4762a1bSJed Brown Mx--; My--; Mz--; 526c4762a1bSJed Brown 527c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(dac, -1.0,1.0, -1.0,1.0, -1.0,1.0);CHKERRQ(ierr); 528c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(daf, -1.0,1.0, -1.0,1.0, -1.0,1.0);CHKERRQ(ierr); 529c4762a1bSJed Brown 530c4762a1bSJed Brown /* apply trilinear mappings */ 531c4762a1bSJed Brown /*ierr = DAApplyTrilinearMapping(dac);CHKERRQ(ierr);*/ 532c4762a1bSJed Brown /* apply conformal mappings */ 533c4762a1bSJed Brown map_id = 0; 534c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-cmap", &map_id,NULL);CHKERRQ(ierr); 535c4762a1bSJed Brown if (map_id >= 1) { 536c4762a1bSJed Brown ierr = DAApplyConformalMapping(dac,map_id);CHKERRQ(ierr); 537c4762a1bSJed Brown } 538c4762a1bSJed Brown 539c4762a1bSJed Brown { 540c4762a1bSJed Brown DM cdaf,cdac; 541c4762a1bSJed Brown Vec coordsc,coordsf; 542c4762a1bSJed Brown 543c4762a1bSJed Brown ierr = DMGetCoordinateDM(dac,&cdac);CHKERRQ(ierr); 544c4762a1bSJed Brown ierr = DMGetCoordinateDM(daf,&cdaf);CHKERRQ(ierr); 545c4762a1bSJed Brown 546c4762a1bSJed Brown ierr = DMGetCoordinates(dac,&coordsc);CHKERRQ(ierr); 547c4762a1bSJed Brown ierr = DMGetCoordinates(daf,&coordsf);CHKERRQ(ierr); 548c4762a1bSJed Brown 549c4762a1bSJed Brown ierr = DMCreateInterpolation(cdac,cdaf,&II,&scale);CHKERRQ(ierr); 550c4762a1bSJed Brown ierr = MatInterpolate(II,coordsc,coordsf);CHKERRQ(ierr); 551c4762a1bSJed Brown ierr = MatDestroy(&II);CHKERRQ(ierr); 552c4762a1bSJed Brown ierr = VecDestroy(&scale);CHKERRQ(ierr); 553c4762a1bSJed Brown } 554c4762a1bSJed Brown 555c4762a1bSJed Brown ierr = DMCreateInterpolation(dac,daf,&INTERP,NULL);CHKERRQ(ierr); 556c4762a1bSJed Brown 557c4762a1bSJed Brown ierr = DMCreateGlobalVector(dac,&ac);CHKERRQ(ierr); 558c4762a1bSJed Brown ierr = VecZeroEntries(ac);CHKERRQ(ierr); 559c4762a1bSJed Brown ierr = DADefineXLinearField3D(dac,ac);CHKERRQ(ierr); 560c4762a1bSJed Brown 561c4762a1bSJed Brown ierr = DMCreateGlobalVector(daf,&af);CHKERRQ(ierr); 562c4762a1bSJed Brown ierr = VecZeroEntries(af);CHKERRQ(ierr); 563c4762a1bSJed Brown 564c4762a1bSJed Brown ierr = MatMult(INTERP,ac, af);CHKERRQ(ierr); 565c4762a1bSJed Brown 566c4762a1bSJed Brown { 567c4762a1bSJed Brown Vec afexact; 568c4762a1bSJed Brown PetscReal nrm; 569c4762a1bSJed Brown PetscInt N; 570c4762a1bSJed Brown 571c4762a1bSJed Brown ierr = DMCreateGlobalVector(daf,&afexact);CHKERRQ(ierr); 572c4762a1bSJed Brown ierr = VecZeroEntries(afexact);CHKERRQ(ierr); 573c4762a1bSJed Brown ierr = DADefineXLinearField3D(daf,afexact);CHKERRQ(ierr); 574c4762a1bSJed Brown ierr = VecAXPY(afexact,-1.0,af);CHKERRQ(ierr); /* af <= af - afinterp */ 575c4762a1bSJed Brown ierr = VecNorm(afexact,NORM_2,&nrm);CHKERRQ(ierr); 576c4762a1bSJed Brown ierr = VecGetSize(afexact,&N);CHKERRQ(ierr); 577*b5675b0fSBarry Smith ierr = 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)));CHKERRQ(ierr); 578c4762a1bSJed Brown ierr = VecDestroy(&afexact);CHKERRQ(ierr); 579c4762a1bSJed Brown } 580c4762a1bSJed Brown 581c4762a1bSJed Brown PetscOptionsGetBool(NULL,NULL,"-output",&output,NULL);CHKERRQ(ierr); 582c4762a1bSJed Brown if (output) { 583c4762a1bSJed Brown ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD, "dac_3D.vtk", &vv);CHKERRQ(ierr); 584c4762a1bSJed Brown ierr = PetscViewerPushFormat(vv, PETSC_VIEWER_ASCII_VTK);CHKERRQ(ierr); 585c4762a1bSJed Brown ierr = DMView(dac, vv);CHKERRQ(ierr); 586c4762a1bSJed Brown ierr = VecView(ac, vv);CHKERRQ(ierr); 587c4762a1bSJed Brown ierr = PetscViewerPopFormat(vv);CHKERRQ(ierr); 588c4762a1bSJed Brown ierr = PetscViewerDestroy(&vv);CHKERRQ(ierr); 589c4762a1bSJed Brown 590c4762a1bSJed Brown ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD, "daf_3D.vtk", &vv);CHKERRQ(ierr); 591c4762a1bSJed Brown ierr = PetscViewerPushFormat(vv, PETSC_VIEWER_ASCII_VTK);CHKERRQ(ierr); 592c4762a1bSJed Brown ierr = DMView(daf, vv);CHKERRQ(ierr); 593c4762a1bSJed Brown ierr = VecView(af, vv);CHKERRQ(ierr); 594c4762a1bSJed Brown ierr = PetscViewerPopFormat(vv);CHKERRQ(ierr); 595c4762a1bSJed Brown ierr = PetscViewerDestroy(&vv);CHKERRQ(ierr); 596c4762a1bSJed Brown } 597c4762a1bSJed Brown 598c4762a1bSJed Brown ierr = MatDestroy(&INTERP);CHKERRQ(ierr); 599c4762a1bSJed Brown ierr = DMDestroy(&dac);CHKERRQ(ierr); 600c4762a1bSJed Brown ierr = DMDestroy(&daf);CHKERRQ(ierr); 601c4762a1bSJed Brown ierr = VecDestroy(&ac);CHKERRQ(ierr); 602c4762a1bSJed Brown ierr = VecDestroy(&af);CHKERRQ(ierr); 603c4762a1bSJed Brown PetscFunctionReturn(0); 604c4762a1bSJed Brown } 605c4762a1bSJed Brown 606c4762a1bSJed Brown int main(int argc,char **argv) 607c4762a1bSJed Brown { 608c4762a1bSJed Brown PetscErrorCode ierr; 609c4762a1bSJed Brown PetscInt mx = 2,my = 2,mz = 2,l,nl,dim; 610c4762a1bSJed Brown 611c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr; 612*b5675b0fSBarry Smith ierr = PetscOptionsGetInt(NULL,NULL,"-mx", &mx, 0);CHKERRQ(ierr); 613*b5675b0fSBarry Smith ierr = PetscOptionsGetInt(NULL,NULL,"-my", &my, 0);CHKERRQ(ierr); 614*b5675b0fSBarry Smith ierr = PetscOptionsGetInt(NULL,NULL,"-mz", &mz, 0);CHKERRQ(ierr); 615c4762a1bSJed Brown nl = 1; 616*b5675b0fSBarry Smith ierr = PetscOptionsGetInt(NULL,NULL,"-nl", &nl, 0);CHKERRQ(ierr); 617c4762a1bSJed Brown dim = 2; 618*b5675b0fSBarry Smith ierr = PetscOptionsGetInt(NULL,NULL,"-dim", &dim, 0);CHKERRQ(ierr); 619c4762a1bSJed Brown 620c4762a1bSJed Brown for (l=0; l<nl; l++) { 621*b5675b0fSBarry Smith if (dim==1) { 622*b5675b0fSBarry Smith ierr = da_test_RefineCoords1D(mx);CHKERRQ(ierr); 623*b5675b0fSBarry Smith } else if (dim==2) { 624*b5675b0fSBarry Smith ierr = da_test_RefineCoords2D(mx,my);CHKERRQ(ierr); 625*b5675b0fSBarry Smith } else if (dim==3) { 626*b5675b0fSBarry Smith ierr = da_test_RefineCoords3D(mx,my,mz);CHKERRQ(ierr); 627*b5675b0fSBarry Smith } 628c4762a1bSJed Brown mx = mx * 2; 629c4762a1bSJed Brown my = my * 2; 630c4762a1bSJed Brown mz = mz * 2; 631c4762a1bSJed Brown } 632c4762a1bSJed Brown ierr = PetscFinalize(); 633c4762a1bSJed Brown return ierr; 634c4762a1bSJed Brown } 635c4762a1bSJed Brown 636c4762a1bSJed Brown 637c4762a1bSJed Brown /*TEST 638c4762a1bSJed Brown 639c4762a1bSJed Brown test: 640c4762a1bSJed Brown suffix: 1d 641c4762a1bSJed Brown args: -mx 10 -nl 6 -dim 1 642c4762a1bSJed Brown 643c4762a1bSJed Brown test: 644c4762a1bSJed Brown suffix: 2d 645c4762a1bSJed Brown 646c4762a1bSJed Brown test: 647c4762a1bSJed Brown args: -mx 10 -my 10 -nl 6 -dim 2 -cmap 0 648c4762a1bSJed Brown 649c4762a1bSJed Brown test: 650c4762a1bSJed Brown args: -mx 10 -my 10 -nl 6 -dim 2 -cmap 1 651c4762a1bSJed Brown 652c4762a1bSJed Brown test: 653c4762a1bSJed Brown args: -mx 10 -my 10 -nl 6 -dim 2 -cmap 2 654c4762a1bSJed Brown 655c4762a1bSJed Brown test: 656c4762a1bSJed Brown args: -mx 10 -my 10 -nl 6 -dim 2 -cmap 3 657c4762a1bSJed Brown 658c4762a1bSJed Brown test: 659c4762a1bSJed Brown suffix: 2dp1 660c4762a1bSJed Brown nsize: 8 661c4762a1bSJed Brown args: -mx 10 -my 10 -nl 4 -dim 2 -cmap 3 -da_refine_x 3 -da_refine_y 4 662c4762a1bSJed Brown timeoutfactor: 2 663c4762a1bSJed Brown 664c4762a1bSJed Brown test: 665c4762a1bSJed Brown suffix: 2dp2 666c4762a1bSJed Brown nsize: 8 667c4762a1bSJed Brown args: -mx 10 -my 10 -nl 4 -dim 2 -cmap 3 -da_refine_x 3 -da_refine_y 1 668c4762a1bSJed Brown timeoutfactor: 2 669c4762a1bSJed Brown 670c4762a1bSJed Brown test: 671c4762a1bSJed Brown suffix: 3d 672c4762a1bSJed Brown args: -mx 5 -my 5 -mz 5 -nl 4 -dim 3 -cmap 3 673c4762a1bSJed Brown 674c4762a1bSJed Brown test: 675c4762a1bSJed Brown suffix: 3dp1 676c4762a1bSJed Brown nsize: 32 677c4762a1bSJed 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 678c4762a1bSJed Brown 679c4762a1bSJed Brown TEST*/ 680