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 */ 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)); 193c4762a1bSJed Brown PetscFunctionReturn(0); 194c4762a1bSJed Brown } 195c4762a1bSJed Brown 196c4762a1bSJed Brown PetscErrorCode DAApplyTrilinearMapping(DM da) 197c4762a1bSJed Brown { 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)); 248c4762a1bSJed Brown PetscFunctionReturn(0); 249c4762a1bSJed Brown } 250c4762a1bSJed Brown 251c4762a1bSJed Brown PetscErrorCode DADefineXLinearField2D(DM da,Vec field) 252c4762a1bSJed Brown { 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++) { 270c4762a1bSJed Brown for (j=sy; j<sy+ny; j++) { 271c4762a1bSJed 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; 272c4762a1bSJed Brown } 273c4762a1bSJed Brown } 274c4762a1bSJed Brown 2759566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,field,&FF)); 2769566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda,Gcoords,&XX)); 277c4762a1bSJed Brown PetscFunctionReturn(0); 278c4762a1bSJed Brown } 279c4762a1bSJed Brown 280c4762a1bSJed Brown PetscErrorCode DADefineXLinearField3D(DM da,Vec field) 281c4762a1bSJed Brown { 282c4762a1bSJed Brown PetscInt i,j,k; 283c4762a1bSJed Brown PetscInt sx,nx,sy,ny,sz,nz; 284c4762a1bSJed Brown Vec Gcoords; 285c4762a1bSJed Brown DMDACoor3d ***XX; 286c4762a1bSJed Brown PetscScalar ***FF; 287c4762a1bSJed Brown DM cda; 288c4762a1bSJed Brown 289c4762a1bSJed Brown PetscFunctionBeginUser; 2909566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da,&cda)); 2919566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da,&Gcoords)); 292c4762a1bSJed Brown 2939566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda,Gcoords,&XX)); 2949566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,field,&FF)); 295c4762a1bSJed Brown 2969566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&sx,&sy,&sz,&nx,&ny,&nz)); 297c4762a1bSJed Brown 298c4762a1bSJed Brown for (k=sz; k<sz+nz; k++) { 299c4762a1bSJed Brown for (j=sy; j<sy+ny; j++) { 300c4762a1bSJed Brown for (i=sx; i<sx+nx; i++) { 301c4762a1bSJed Brown FF[k][j][i] = 10.0 302c4762a1bSJed Brown + 4.05 * XX[k][j][i].x 303c4762a1bSJed Brown + 5.50 * XX[k][j][i].y 304c4762a1bSJed Brown + 1.33 * XX[k][j][i].z 305c4762a1bSJed Brown + 2.03 * XX[k][j][i].x * XX[k][j][i].y 306c4762a1bSJed Brown + 0.03 * XX[k][j][i].x * XX[k][j][i].z 307c4762a1bSJed Brown + 0.83 * XX[k][j][i].y * XX[k][j][i].z 308c4762a1bSJed Brown + 3.79 * XX[k][j][i].x * XX[k][j][i].y * XX[k][j][i].z; 309c4762a1bSJed Brown } 310c4762a1bSJed Brown } 311c4762a1bSJed Brown } 312c4762a1bSJed Brown 3139566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,field,&FF)); 3149566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda,Gcoords,&XX)); 315c4762a1bSJed Brown PetscFunctionReturn(0); 316c4762a1bSJed Brown } 317c4762a1bSJed Brown 318c4762a1bSJed Brown PetscErrorCode da_test_RefineCoords1D(PetscInt mx) 319c4762a1bSJed Brown { 320c4762a1bSJed Brown DM dac,daf; 321c4762a1bSJed Brown PetscViewer vv; 322c4762a1bSJed Brown Vec ac,af; 323c4762a1bSJed Brown PetscInt Mx; 324c4762a1bSJed Brown Mat II,INTERP; 325c4762a1bSJed Brown Vec scale; 326c4762a1bSJed Brown PetscBool output = PETSC_FALSE; 327c4762a1bSJed Brown 328c4762a1bSJed Brown PetscFunctionBeginUser; 3299566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE,mx+1,1, /* 1 dof */1, /* stencil = 1 */NULL,&dac)); 3309566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dac)); 3319566063dSJacob Faibussowitsch PetscCall(DMSetUp(dac)); 332c4762a1bSJed Brown 3339566063dSJacob Faibussowitsch PetscCall(DMRefine(dac,MPI_COMM_NULL,&daf)); 3349566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf,0,&Mx,0,0,0,0,0,0,0,0,0,0,0)); 335c4762a1bSJed Brown Mx--; 336c4762a1bSJed Brown 3379566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(dac, -1.0,1.0, PETSC_DECIDE,PETSC_DECIDE, PETSC_DECIDE,PETSC_DECIDE)); 3389566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(daf, -1.0,1.0, PETSC_DECIDE,PETSC_DECIDE, PETSC_DECIDE,PETSC_DECIDE)); 339c4762a1bSJed Brown 340c4762a1bSJed Brown { 341c4762a1bSJed Brown DM cdaf,cdac; 342c4762a1bSJed Brown Vec coordsc,coordsf; 343c4762a1bSJed Brown 3449566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dac,&cdac)); 3459566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(daf,&cdaf)); 346c4762a1bSJed Brown 3479566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dac,&coordsc)); 3489566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(daf,&coordsf)); 349c4762a1bSJed Brown 3509566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(cdac,cdaf,&II,&scale)); 3519566063dSJacob Faibussowitsch PetscCall(MatInterpolate(II,coordsc,coordsf)); 3529566063dSJacob Faibussowitsch PetscCall(MatDestroy(&II)); 3539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&scale)); 354c4762a1bSJed Brown } 355c4762a1bSJed Brown 3569566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dac,daf,&INTERP,NULL)); 357c4762a1bSJed Brown 3589566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dac,&ac)); 3599566063dSJacob Faibussowitsch PetscCall(VecSet(ac,66.99)); 360c4762a1bSJed Brown 3619566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(daf,&af)); 362c4762a1bSJed Brown 3639566063dSJacob Faibussowitsch PetscCall(MatMult(INTERP,ac, af)); 364c4762a1bSJed Brown 365c4762a1bSJed Brown { 366c4762a1bSJed Brown Vec afexact; 367c4762a1bSJed Brown PetscReal nrm; 368c4762a1bSJed Brown PetscInt N; 369c4762a1bSJed Brown 3709566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(daf,&afexact)); 3719566063dSJacob Faibussowitsch PetscCall(VecSet(afexact,66.99)); 3729566063dSJacob Faibussowitsch PetscCall(VecAXPY(afexact,-1.0,af)); /* af <= af - afinterp */ 3739566063dSJacob Faibussowitsch PetscCall(VecNorm(afexact,NORM_2,&nrm)); 3749566063dSJacob Faibussowitsch PetscCall(VecGetSize(afexact,&N)); 37563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%" PetscInt_FMT "=>%" PetscInt_FMT ", interp err = %1.4e\n",mx,Mx,(double)(nrm/PetscSqrtReal((PetscReal)N)))); 3769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&afexact)); 377c4762a1bSJed Brown } 378c4762a1bSJed Brown 3799566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-output",&output,NULL)); 380c4762a1bSJed Brown if (output) { 3819566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "dac_1D.vtr", &vv)); 3829566063dSJacob Faibussowitsch PetscCall(VecView(ac, vv)); 3839566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vv)); 384c4762a1bSJed Brown 3859566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "daf_1D.vtr", &vv)); 3869566063dSJacob Faibussowitsch PetscCall(VecView(af, vv)); 3879566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vv)); 388c4762a1bSJed Brown } 389c4762a1bSJed Brown 3909566063dSJacob Faibussowitsch PetscCall(MatDestroy(&INTERP)); 3919566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dac)); 3929566063dSJacob Faibussowitsch PetscCall(DMDestroy(&daf)); 3939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ac)); 3949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&af)); 395c4762a1bSJed Brown PetscFunctionReturn(0); 396c4762a1bSJed Brown } 397c4762a1bSJed Brown 398c4762a1bSJed Brown PetscErrorCode da_test_RefineCoords2D(PetscInt mx,PetscInt my) 399c4762a1bSJed Brown { 400c4762a1bSJed Brown DM dac,daf; 401c4762a1bSJed Brown PetscViewer vv; 402c4762a1bSJed Brown Vec ac,af; 403c4762a1bSJed Brown PetscInt map_id,Mx,My; 404c4762a1bSJed Brown Mat II,INTERP; 405c4762a1bSJed Brown Vec scale; 406c4762a1bSJed Brown PetscBool output = PETSC_FALSE; 407c4762a1bSJed Brown 408c4762a1bSJed Brown PetscFunctionBeginUser; 4099566063dSJacob 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)); 4109566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dac)); 4119566063dSJacob Faibussowitsch PetscCall(DMSetUp(dac)); 412c4762a1bSJed Brown 4139566063dSJacob Faibussowitsch PetscCall(DMRefine(dac,MPI_COMM_NULL,&daf)); 4149566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf,0,&Mx,&My,0,0,0,0,0,0,0,0,0,0)); 415c4762a1bSJed Brown Mx--; My--; 416c4762a1bSJed Brown 4179566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(dac, -1.0,1.0, -1.0,1.0, PETSC_DECIDE,PETSC_DECIDE)); 4189566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(daf, -1.0,1.0, -1.0,1.0, PETSC_DECIDE,PETSC_DECIDE)); 419c4762a1bSJed Brown 420c4762a1bSJed Brown /* apply conformal mappings */ 421c4762a1bSJed Brown map_id = 0; 4229566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-cmap", &map_id,NULL)); 423c4762a1bSJed Brown if (map_id >= 1) { 4249566063dSJacob Faibussowitsch PetscCall(DAApplyConformalMapping(dac,map_id)); 425c4762a1bSJed Brown } 426c4762a1bSJed Brown 427c4762a1bSJed Brown { 428c4762a1bSJed Brown DM cdaf,cdac; 429c4762a1bSJed Brown Vec coordsc,coordsf; 430c4762a1bSJed Brown 4319566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dac,&cdac)); 4329566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(daf,&cdaf)); 433c4762a1bSJed Brown 4349566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dac,&coordsc)); 4359566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(daf,&coordsf)); 436c4762a1bSJed Brown 4379566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(cdac,cdaf,&II,&scale)); 4389566063dSJacob Faibussowitsch PetscCall(MatInterpolate(II,coordsc,coordsf)); 4399566063dSJacob Faibussowitsch PetscCall(MatDestroy(&II)); 4409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&scale)); 441c4762a1bSJed Brown } 442c4762a1bSJed Brown 4439566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dac,daf,&INTERP,NULL)); 444c4762a1bSJed Brown 4459566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dac,&ac)); 4469566063dSJacob Faibussowitsch PetscCall(DADefineXLinearField2D(dac,ac)); 447c4762a1bSJed Brown 4489566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(daf,&af)); 4499566063dSJacob Faibussowitsch PetscCall(MatMult(INTERP,ac, af)); 450c4762a1bSJed Brown 451c4762a1bSJed Brown { 452c4762a1bSJed Brown Vec afexact; 453c4762a1bSJed Brown PetscReal nrm; 454c4762a1bSJed Brown PetscInt N; 455c4762a1bSJed Brown 4569566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(daf,&afexact)); 4579566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(afexact)); 4589566063dSJacob Faibussowitsch PetscCall(DADefineXLinearField2D(daf,afexact)); 4599566063dSJacob Faibussowitsch PetscCall(VecAXPY(afexact,-1.0,af)); /* af <= af - afinterp */ 4609566063dSJacob Faibussowitsch PetscCall(VecNorm(afexact,NORM_2,&nrm)); 4619566063dSJacob Faibussowitsch PetscCall(VecGetSize(afexact,&N)); 46263a3b9bcSJacob 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)))); 4639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&afexact)); 464c4762a1bSJed Brown } 465c4762a1bSJed Brown 4669566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-output",&output,NULL)); 467c4762a1bSJed Brown if (output) { 4689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "dac_2D.vtr", &vv)); 4699566063dSJacob Faibussowitsch PetscCall(VecView(ac, vv)); 4709566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vv)); 471c4762a1bSJed Brown 4729566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "daf_2D.vtr", &vv)); 4739566063dSJacob Faibussowitsch PetscCall(VecView(af, vv)); 4749566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vv)); 475c4762a1bSJed Brown } 476c4762a1bSJed Brown 4779566063dSJacob Faibussowitsch PetscCall(MatDestroy(&INTERP)); 4789566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dac)); 4799566063dSJacob Faibussowitsch PetscCall(DMDestroy(&daf)); 4809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ac)); 4819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&af)); 482c4762a1bSJed Brown PetscFunctionReturn(0); 483c4762a1bSJed Brown } 484c4762a1bSJed Brown 485c4762a1bSJed Brown PetscErrorCode da_test_RefineCoords3D(PetscInt mx,PetscInt my,PetscInt mz) 486c4762a1bSJed Brown { 487c4762a1bSJed Brown DM dac,daf; 488c4762a1bSJed Brown PetscViewer vv; 489c4762a1bSJed Brown Vec ac,af; 490c4762a1bSJed Brown PetscInt map_id,Mx,My,Mz; 491c4762a1bSJed Brown Mat II,INTERP; 492c4762a1bSJed Brown Vec scale; 493c4762a1bSJed Brown PetscBool output = PETSC_FALSE; 494c4762a1bSJed Brown 495c4762a1bSJed Brown PetscFunctionBeginUser; 496d0609cedSBarry 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 */ 497d0609cedSBarry Smith 1, /* stencil = 1 */NULL,NULL,NULL,&dac)); 4989566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dac)); 4999566063dSJacob Faibussowitsch PetscCall(DMSetUp(dac)); 500c4762a1bSJed Brown 5019566063dSJacob Faibussowitsch PetscCall(DMRefine(dac,MPI_COMM_NULL,&daf)); 5029566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daf,0,&Mx,&My,&Mz,0,0,0,0,0,0,0,0,0)); 503c4762a1bSJed Brown Mx--; My--; Mz--; 504c4762a1bSJed Brown 5059566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(dac, -1.0,1.0, -1.0,1.0, -1.0,1.0)); 5069566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(daf, -1.0,1.0, -1.0,1.0, -1.0,1.0)); 507c4762a1bSJed Brown 508c4762a1bSJed Brown /* apply trilinear mappings */ 5099566063dSJacob Faibussowitsch /*PetscCall(DAApplyTrilinearMapping(dac));*/ 510c4762a1bSJed Brown /* apply conformal mappings */ 511c4762a1bSJed Brown map_id = 0; 5129566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-cmap", &map_id,NULL)); 513c4762a1bSJed Brown if (map_id >= 1) { 5149566063dSJacob Faibussowitsch PetscCall(DAApplyConformalMapping(dac,map_id)); 515c4762a1bSJed Brown } 516c4762a1bSJed Brown 517c4762a1bSJed Brown { 518c4762a1bSJed Brown DM cdaf,cdac; 519c4762a1bSJed Brown Vec coordsc,coordsf; 520c4762a1bSJed Brown 5219566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dac,&cdac)); 5229566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(daf,&cdaf)); 523c4762a1bSJed Brown 5249566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dac,&coordsc)); 5259566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(daf,&coordsf)); 526c4762a1bSJed Brown 5279566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(cdac,cdaf,&II,&scale)); 5289566063dSJacob Faibussowitsch PetscCall(MatInterpolate(II,coordsc,coordsf)); 5299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&II)); 5309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&scale)); 531c4762a1bSJed Brown } 532c4762a1bSJed Brown 5339566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dac,daf,&INTERP,NULL)); 534c4762a1bSJed Brown 5359566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dac,&ac)); 5369566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(ac)); 5379566063dSJacob Faibussowitsch PetscCall(DADefineXLinearField3D(dac,ac)); 538c4762a1bSJed Brown 5399566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(daf,&af)); 5409566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(af)); 541c4762a1bSJed Brown 5429566063dSJacob Faibussowitsch PetscCall(MatMult(INTERP,ac, af)); 543c4762a1bSJed Brown 544c4762a1bSJed Brown { 545c4762a1bSJed Brown Vec afexact; 546c4762a1bSJed Brown PetscReal nrm; 547c4762a1bSJed Brown PetscInt N; 548c4762a1bSJed Brown 5499566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(daf,&afexact)); 5509566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(afexact)); 5519566063dSJacob Faibussowitsch PetscCall(DADefineXLinearField3D(daf,afexact)); 5529566063dSJacob Faibussowitsch PetscCall(VecAXPY(afexact,-1.0,af)); /* af <= af - afinterp */ 5539566063dSJacob Faibussowitsch PetscCall(VecNorm(afexact,NORM_2,&nrm)); 5549566063dSJacob Faibussowitsch PetscCall(VecGetSize(afexact,&N)); 55563a3b9bcSJacob 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)))); 5569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&afexact)); 557c4762a1bSJed Brown } 558c4762a1bSJed Brown 5599566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-output",&output,NULL)); 560c4762a1bSJed Brown if (output) { 5619566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "dac_3D.vtr", &vv)); 5629566063dSJacob Faibussowitsch PetscCall(VecView(ac, vv)); 5639566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vv)); 564c4762a1bSJed Brown 5659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "daf_3D.vtr", &vv)); 5669566063dSJacob Faibussowitsch PetscCall(VecView(af, vv)); 5679566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vv)); 568c4762a1bSJed Brown } 569c4762a1bSJed Brown 5709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&INTERP)); 5719566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dac)); 5729566063dSJacob Faibussowitsch PetscCall(DMDestroy(&daf)); 5739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ac)); 5749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&af)); 575c4762a1bSJed Brown PetscFunctionReturn(0); 576c4762a1bSJed Brown } 577c4762a1bSJed Brown 578c4762a1bSJed Brown int main(int argc,char **argv) 579c4762a1bSJed Brown { 580c4762a1bSJed Brown PetscInt mx = 2,my = 2,mz = 2,l,nl,dim; 581c4762a1bSJed Brown 582*327415f7SBarry Smith PetscFunctionBeginUser; 5839566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,0,help)); 5849566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mx", &mx, 0)); 5859566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-my", &my, 0)); 5869566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mz", &mz, 0)); 587c4762a1bSJed Brown nl = 1; 5889566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-nl", &nl, 0)); 589c4762a1bSJed Brown dim = 2; 5909566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-dim", &dim, 0)); 591c4762a1bSJed Brown 592c4762a1bSJed Brown for (l=0; l<nl; l++) { 5931baa6e33SBarry Smith if (dim==1) PetscCall(da_test_RefineCoords1D(mx)); 5941baa6e33SBarry Smith else if (dim==2) PetscCall(da_test_RefineCoords2D(mx,my)); 5951baa6e33SBarry Smith else if (dim==3) PetscCall(da_test_RefineCoords3D(mx,my,mz)); 596c4762a1bSJed Brown mx = mx * 2; 597c4762a1bSJed Brown my = my * 2; 598c4762a1bSJed Brown mz = mz * 2; 599c4762a1bSJed Brown } 6009566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 601b122ec5aSJacob Faibussowitsch return 0; 602c4762a1bSJed Brown } 603c4762a1bSJed Brown 604c4762a1bSJed Brown /*TEST 605c4762a1bSJed Brown 606c4762a1bSJed Brown test: 607c4762a1bSJed Brown suffix: 1d 608c4762a1bSJed Brown args: -mx 10 -nl 6 -dim 1 609c4762a1bSJed Brown 610c4762a1bSJed Brown test: 611c4762a1bSJed Brown suffix: 2d 61263037964SStefano Zampini output_file: output/ex36_2d.out 61363037964SStefano Zampini args: -mx 10 -my 10 -nl 6 -dim 2 -cmap {{0 1 2 3}} 614c4762a1bSJed Brown 615c4762a1bSJed Brown test: 616c4762a1bSJed Brown suffix: 2dp1 617c4762a1bSJed Brown nsize: 8 618c4762a1bSJed Brown args: -mx 10 -my 10 -nl 4 -dim 2 -cmap 3 -da_refine_x 3 -da_refine_y 4 619c4762a1bSJed Brown timeoutfactor: 2 620c4762a1bSJed Brown 621c4762a1bSJed Brown test: 622c4762a1bSJed Brown suffix: 2dp2 623c4762a1bSJed Brown nsize: 8 624c4762a1bSJed Brown args: -mx 10 -my 10 -nl 4 -dim 2 -cmap 3 -da_refine_x 3 -da_refine_y 1 625c4762a1bSJed Brown timeoutfactor: 2 626c4762a1bSJed Brown 627c4762a1bSJed Brown test: 628c4762a1bSJed Brown suffix: 3d 629c4762a1bSJed Brown args: -mx 5 -my 5 -mz 5 -nl 4 -dim 3 -cmap 3 630c4762a1bSJed Brown 631c4762a1bSJed Brown test: 632c4762a1bSJed Brown suffix: 3dp1 633c4762a1bSJed Brown nsize: 32 634c4762a1bSJed 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 635c4762a1bSJed Brown 636c4762a1bSJed Brown TEST*/ 637