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