1*c4762a1bSJed Brown static char help[] = "Tests DMDA ghost coordinates\n\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscdm.h> 4*c4762a1bSJed Brown #include <petscdmda.h> 5*c4762a1bSJed Brown 6*c4762a1bSJed Brown static PetscErrorCode CompareGhostedCoords(Vec gc1,Vec gc2) 7*c4762a1bSJed Brown { 8*c4762a1bSJed Brown PetscErrorCode ierr; 9*c4762a1bSJed Brown PetscReal nrm,gnrm; 10*c4762a1bSJed Brown Vec tmp; 11*c4762a1bSJed Brown 12*c4762a1bSJed Brown PetscFunctionBeginUser; 13*c4762a1bSJed Brown ierr = VecDuplicate(gc1,&tmp);CHKERRQ(ierr); 14*c4762a1bSJed Brown ierr = VecWAXPY(tmp,-1.0,gc1,gc2);CHKERRQ(ierr); 15*c4762a1bSJed Brown ierr = VecNorm(tmp,NORM_INFINITY,&nrm);CHKERRQ(ierr); 16*c4762a1bSJed Brown ierr = MPI_Allreduce(&nrm,&gnrm,1,MPIU_REAL,MPIU_MAX,PETSC_COMM_WORLD);CHKERRQ(ierr); 17*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"norm of difference of ghosted coordinates %8.2e\n",gnrm);CHKERRQ(ierr); 18*c4762a1bSJed Brown ierr = VecDestroy(&tmp);CHKERRQ(ierr); 19*c4762a1bSJed Brown PetscFunctionReturn(0); 20*c4762a1bSJed Brown } 21*c4762a1bSJed Brown 22*c4762a1bSJed Brown static PetscErrorCode TestQ2Q1DA(void) 23*c4762a1bSJed Brown { 24*c4762a1bSJed Brown DM Q2_da,Q1_da,cda; 25*c4762a1bSJed Brown PetscInt mx,my,mz; 26*c4762a1bSJed Brown Vec coords,gcoords,gcoords2; 27*c4762a1bSJed Brown PetscErrorCode ierr; 28*c4762a1bSJed Brown 29*c4762a1bSJed Brown mx = 7; 30*c4762a1bSJed Brown my = 11; 31*c4762a1bSJed Brown mz = 13; 32*c4762a1bSJed Brown ierr = DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,mx,my,mz,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,3,2,0,0,0,&Q2_da);CHKERRQ(ierr); 33*c4762a1bSJed Brown ierr = DMSetFromOptions(Q2_da);CHKERRQ(ierr); 34*c4762a1bSJed Brown ierr = DMSetUp(Q2_da);CHKERRQ(ierr); 35*c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(Q2_da,-1.0,1.0,-2.0,2.0,-3.0,3.0);CHKERRQ(ierr); 36*c4762a1bSJed Brown ierr = DMGetCoordinates(Q2_da,&coords);CHKERRQ(ierr); 37*c4762a1bSJed Brown ierr = DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,mx,my,mz,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,3,1,0,0,0,&Q1_da);CHKERRQ(ierr); 38*c4762a1bSJed Brown ierr = DMSetFromOptions(Q1_da);CHKERRQ(ierr); 39*c4762a1bSJed Brown ierr = DMSetUp(Q1_da);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = DMSetCoordinates(Q1_da,coords);CHKERRQ(ierr); 41*c4762a1bSJed Brown 42*c4762a1bSJed Brown /* Get ghost coordinates one way */ 43*c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(Q1_da,&gcoords);CHKERRQ(ierr); 44*c4762a1bSJed Brown 45*c4762a1bSJed Brown /* And another */ 46*c4762a1bSJed Brown ierr = DMGetCoordinates(Q1_da,&coords);CHKERRQ(ierr); 47*c4762a1bSJed Brown ierr = DMGetCoordinateDM(Q1_da,&cda);CHKERRQ(ierr); 48*c4762a1bSJed Brown ierr = DMGetLocalVector(cda,&gcoords2);CHKERRQ(ierr); 49*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(cda,coords,INSERT_VALUES,gcoords2);CHKERRQ(ierr); 50*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(cda,coords,INSERT_VALUES,gcoords2);CHKERRQ(ierr); 51*c4762a1bSJed Brown 52*c4762a1bSJed Brown ierr = CompareGhostedCoords(gcoords,gcoords2);CHKERRQ(ierr); 53*c4762a1bSJed Brown ierr = DMRestoreLocalVector(cda,&gcoords2);CHKERRQ(ierr); 54*c4762a1bSJed Brown 55*c4762a1bSJed Brown ierr = VecScale(coords,10.0);CHKERRQ(ierr); 56*c4762a1bSJed Brown ierr = VecScale(gcoords,10.0);CHKERRQ(ierr); 57*c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(Q1_da,&gcoords2);CHKERRQ(ierr); 58*c4762a1bSJed Brown ierr = CompareGhostedCoords(gcoords,gcoords2);CHKERRQ(ierr); 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown ierr = DMDestroy(&Q2_da);CHKERRQ(ierr); 61*c4762a1bSJed Brown ierr = DMDestroy(&Q1_da);CHKERRQ(ierr); 62*c4762a1bSJed Brown PetscFunctionReturn(0); 63*c4762a1bSJed Brown } 64*c4762a1bSJed Brown 65*c4762a1bSJed Brown int main(int argc,char **argv) 66*c4762a1bSJed Brown { 67*c4762a1bSJed Brown PetscErrorCode ierr; 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr; 70*c4762a1bSJed Brown ierr = TestQ2Q1DA();CHKERRQ(ierr); 71*c4762a1bSJed Brown ierr = PetscFinalize(); 72*c4762a1bSJed Brown return ierr; 73*c4762a1bSJed Brown } 74*c4762a1bSJed Brown 75*c4762a1bSJed Brown 76*c4762a1bSJed Brown /*TEST 77*c4762a1bSJed Brown 78*c4762a1bSJed Brown test: 79*c4762a1bSJed Brown nsize: 2 80*c4762a1bSJed Brown 81*c4762a1bSJed Brown TEST*/ 82