1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Tests various DMComposite routines.\n\n"; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown #include <petscdm.h> 5*c4762a1bSJed Brown #include <petscdmda.h> 6*c4762a1bSJed Brown #include <petscdmcomposite.h> 7*c4762a1bSJed Brown 8*c4762a1bSJed Brown int main(int argc,char **argv) 9*c4762a1bSJed Brown { 10*c4762a1bSJed Brown PetscMPIInt rank; 11*c4762a1bSJed Brown PetscErrorCode ierr; 12*c4762a1bSJed Brown DM da1,da2,packer; 13*c4762a1bSJed Brown Vec local,global,globals[2],buffer; 14*c4762a1bSJed Brown PetscScalar value; 15*c4762a1bSJed Brown PetscViewer viewer; 16*c4762a1bSJed Brown 17*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 18*c4762a1bSJed Brown 19*c4762a1bSJed Brown ierr = DMCompositeCreate(PETSC_COMM_WORLD,&packer);CHKERRQ(ierr); 20*c4762a1bSJed Brown ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,8,1,1,NULL,&da1);CHKERRQ(ierr); 21*c4762a1bSJed Brown ierr = DMSetFromOptions(da1);CHKERRQ(ierr); 22*c4762a1bSJed Brown ierr = DMSetUp(da1);CHKERRQ(ierr); 23*c4762a1bSJed Brown ierr = DMCompositeAddDM(packer,da1);CHKERRQ(ierr); 24*c4762a1bSJed Brown ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,6,1,1,NULL,&da2);CHKERRQ(ierr); 25*c4762a1bSJed Brown ierr = DMSetFromOptions(da2);CHKERRQ(ierr); 26*c4762a1bSJed Brown ierr = DMSetUp(da2);CHKERRQ(ierr); 27*c4762a1bSJed Brown ierr = DMCompositeAddDM(packer,da2);CHKERRQ(ierr); 28*c4762a1bSJed Brown 29*c4762a1bSJed Brown ierr = DMCreateGlobalVector(packer,&global);CHKERRQ(ierr); 30*c4762a1bSJed Brown ierr = DMCreateLocalVector(packer,&local);CHKERRQ(ierr); 31*c4762a1bSJed Brown ierr = DMCreateLocalVector(packer,&buffer);CHKERRQ(ierr); 32*c4762a1bSJed Brown 33*c4762a1bSJed Brown ierr = DMCompositeGetAccessArray(packer,global,2,NULL,globals);CHKERRQ(ierr); 34*c4762a1bSJed Brown value = 1; 35*c4762a1bSJed Brown ierr = VecSet(globals[0], value);CHKERRQ(ierr); 36*c4762a1bSJed Brown value = -1; 37*c4762a1bSJed Brown ierr = VecSet(globals[1], value);CHKERRQ(ierr); 38*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 39*c4762a1bSJed Brown value = rank + 1; 40*c4762a1bSJed Brown ierr = VecScale(globals[0], value);CHKERRQ(ierr); 41*c4762a1bSJed Brown ierr = VecScale(globals[1], value);CHKERRQ(ierr); 42*c4762a1bSJed Brown ierr = DMCompositeRestoreAccessArray(packer,global,2,NULL,globals);CHKERRQ(ierr); 43*c4762a1bSJed Brown 44*c4762a1bSJed Brown /* Test GlobalToLocal in insert mode */ 45*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(packer,global,INSERT_VALUES,local);CHKERRQ(ierr); 46*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(packer,global,INSERT_VALUES,local);CHKERRQ(ierr); 47*c4762a1bSJed Brown 48*c4762a1bSJed Brown ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 49*c4762a1bSJed Brown ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"\nLocal Vector: processor %d\n",rank);CHKERRQ(ierr); 50*c4762a1bSJed Brown ierr = PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&viewer);CHKERRQ(ierr); 51*c4762a1bSJed Brown ierr = VecView(local,viewer);CHKERRQ(ierr); 52*c4762a1bSJed Brown ierr = PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&viewer);CHKERRQ(ierr); 53*c4762a1bSJed Brown ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 54*c4762a1bSJed Brown ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 55*c4762a1bSJed Brown 56*c4762a1bSJed Brown /* Test LocalToGlobal in insert mode */ 57*c4762a1bSJed Brown ierr = DMLocalToGlobalBegin(packer,local,INSERT_VALUES,global);CHKERRQ(ierr); 58*c4762a1bSJed Brown ierr = DMLocalToGlobalEnd(packer,local,INSERT_VALUES,global);CHKERRQ(ierr); 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown ierr = VecView(global,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 61*c4762a1bSJed Brown 62*c4762a1bSJed Brown /* Test LocalToLocal in insert mode */ 63*c4762a1bSJed Brown ierr = DMLocalToLocalBegin(packer,local,INSERT_VALUES,buffer);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = DMLocalToLocalEnd(packer,local,INSERT_VALUES,buffer);CHKERRQ(ierr); 65*c4762a1bSJed Brown 66*c4762a1bSJed Brown ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 67*c4762a1bSJed Brown ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"\nLocal Vector: processor %d\n",rank);CHKERRQ(ierr); 68*c4762a1bSJed Brown ierr = PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&viewer);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = VecView(buffer,viewer);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&viewer);CHKERRQ(ierr); 71*c4762a1bSJed Brown ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 72*c4762a1bSJed Brown ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 73*c4762a1bSJed Brown 74*c4762a1bSJed Brown ierr = VecDestroy(&buffer);CHKERRQ(ierr); 75*c4762a1bSJed Brown ierr = VecDestroy(&local);CHKERRQ(ierr); 76*c4762a1bSJed Brown ierr = VecDestroy(&global);CHKERRQ(ierr); 77*c4762a1bSJed Brown ierr = DMDestroy(&packer);CHKERRQ(ierr); 78*c4762a1bSJed Brown ierr = DMDestroy(&da2);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = DMDestroy(&da1);CHKERRQ(ierr); 80*c4762a1bSJed Brown 81*c4762a1bSJed Brown ierr = PetscFinalize(); 82*c4762a1bSJed Brown 83*c4762a1bSJed Brown return ierr; 84*c4762a1bSJed Brown } 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown 87*c4762a1bSJed Brown /*TEST 88*c4762a1bSJed Brown 89*c4762a1bSJed Brown test: 90*c4762a1bSJed Brown nsize: 3 91*c4762a1bSJed Brown 92*c4762a1bSJed Brown TEST*/ 93