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