1 2 static char help[] = "Tests various DM routines.\n\n"; 3 4 #include <petscdm.h> 5 #include <petscdmda.h> 6 7 int main(int argc,char **argv) 8 { 9 PetscMPIInt rank; 10 PetscInt M = 10,N = 8,m = PETSC_DECIDE,n = PETSC_DECIDE; 11 DM da; 12 PetscViewer viewer; 13 Vec local,global; 14 PetscScalar value; 15 16 PetscFunctionBeginUser; 17 PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 18 PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",300,0,300,300,&viewer)); 19 20 /* Read options */ 21 PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 22 PetscCall(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL)); 23 PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 24 PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 25 26 /* Create distributed array and get vectors */ 27 PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,m,n,1,1,NULL,NULL,&da)); 28 PetscCall(DMSetFromOptions(da)); 29 PetscCall(DMSetUp(da)); 30 PetscCall(DMCreateGlobalVector(da,&global)); 31 PetscCall(DMCreateLocalVector(da,&local)); 32 33 value = -3.0; 34 PetscCall(VecSet(global,value)); 35 PetscCall(DMGlobalToLocalBegin(da,global,INSERT_VALUES,local)); 36 PetscCall(DMGlobalToLocalEnd(da,global,INSERT_VALUES,local)); 37 38 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 39 value = rank+1; 40 PetscCall(VecScale(local,value)); 41 PetscCall(DMLocalToGlobalBegin(da,local,ADD_VALUES,global)); 42 PetscCall(DMLocalToGlobalEnd(da,local,ADD_VALUES,global)); 43 44 PetscCall(VecView(global,PETSC_VIEWER_STDOUT_WORLD)); 45 PetscCall(DMView(da,viewer)); 46 47 /* Free memory */ 48 PetscCall(PetscViewerDestroy(&viewer)); 49 PetscCall(VecDestroy(&local)); 50 PetscCall(VecDestroy(&global)); 51 PetscCall(DMDestroy(&da)); 52 PetscCall(PetscFinalize()); 53 return 0; 54 } 55 56 /*TEST 57 58 test: 59 nsize: 2 60 args: -nox 61 filter: grep -v -i Object 62 63 test: 64 suffix: cuda1 65 requires: cuda 66 args: -dm_vec_type cuda -nox 67 filter: grep -v -i Object 68 69 test: 70 suffix: cuda2 71 nsize: 2 72 requires: cuda 73 args: -dm_vec_type cuda -nox 74 filter: grep -v -i Object 75 76 TEST*/ 77