1 2 static char help[] = "Tests DMDAGlobalToNaturalAllCreate() using contour plotting for 2d DMDAs.\n\n"; 3 4 #include <petscdm.h> 5 #include <petscdmda.h> 6 #include <petscdraw.h> 7 8 int main(int argc,char **argv) 9 { 10 PetscInt i,j,M = 10,N = 8,m = PETSC_DECIDE,n = PETSC_DECIDE; 11 PetscMPIInt rank; 12 PetscErrorCode ierr; 13 PetscBool flg = PETSC_FALSE; 14 DM da; 15 PetscViewer viewer; 16 Vec localall,global; 17 PetscScalar value,*vlocal; 18 DMBoundaryType bx = DM_BOUNDARY_NONE,by = DM_BOUNDARY_NONE; 19 DMDAStencilType stype = DMDA_STENCIL_BOX; 20 VecScatter tolocalall,fromlocalall; 21 PetscInt start,end; 22 23 ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 24 CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",300,0,300,300,&viewer)); 25 26 /* Read options */ 27 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 28 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL)); 29 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 30 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 31 CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-star_stencil",&flg,NULL)); 32 if (flg) stype = DMDA_STENCIL_STAR; 33 34 /* Create distributed array and get vectors */ 35 CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,bx,by,stype,M,N,m,n,1,1,NULL,NULL,&da)); 36 CHKERRQ(DMSetFromOptions(da)); 37 CHKERRQ(DMSetUp(da)); 38 39 CHKERRQ(DMCreateGlobalVector(da,&global)); 40 CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,M*N,&localall)); 41 42 CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 43 CHKERRQ(VecGetOwnershipRange(global,&start,&end)); 44 for (i=start; i<end; i++) { 45 value = 5.0*rank; 46 CHKERRQ(VecSetValues(global,1,&i,&value,INSERT_VALUES)); 47 } 48 CHKERRQ(VecView(global,viewer)); 49 50 /* 51 Create Scatter from global DMDA parallel vector to local vector that 52 contains all entries 53 */ 54 CHKERRQ(DMDAGlobalToNaturalAllCreate(da,&tolocalall)); 55 CHKERRQ(DMDANaturalAllToGlobalCreate(da,&fromlocalall)); 56 57 CHKERRQ(VecScatterBegin(tolocalall,global,localall,INSERT_VALUES,SCATTER_FORWARD)); 58 CHKERRQ(VecScatterEnd(tolocalall,global,localall,INSERT_VALUES,SCATTER_FORWARD)); 59 60 CHKERRQ(VecGetArray(localall,&vlocal)); 61 for (j=0; j<N; j++) { 62 for (i=0; i<M; i++) { 63 *vlocal++ += i + j*M; 64 } 65 } 66 CHKERRQ(VecRestoreArray(localall,&vlocal)); 67 68 /* scatter back to global vector */ 69 CHKERRQ(VecScatterBegin(fromlocalall,localall,global,INSERT_VALUES,SCATTER_FORWARD)); 70 CHKERRQ(VecScatterEnd(fromlocalall,localall,global,INSERT_VALUES,SCATTER_FORWARD)); 71 72 CHKERRQ(VecView(global,viewer)); 73 74 /* Free memory */ 75 CHKERRQ(VecScatterDestroy(&tolocalall)); 76 CHKERRQ(VecScatterDestroy(&fromlocalall)); 77 CHKERRQ(PetscViewerDestroy(&viewer)); 78 CHKERRQ(VecDestroy(&localall)); 79 CHKERRQ(VecDestroy(&global)); 80 CHKERRQ(DMDestroy(&da)); 81 ierr = PetscFinalize(); 82 return ierr; 83 } 84 85 /*TEST 86 87 build: 88 requires: !complex 89 90 test: 91 nsize: 3 92 93 TEST*/ 94