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