1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests various DMComposite routines.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscdm.h> 5c4762a1bSJed Brown #include <petscdmda.h> 6c4762a1bSJed Brown #include <petscdmcomposite.h> 7c4762a1bSJed Brown 8c4762a1bSJed Brown int main(int argc,char **argv) 9c4762a1bSJed Brown { 10c4762a1bSJed Brown PetscMPIInt rank; 11c4762a1bSJed Brown DM da1,da2,packer; 12c4762a1bSJed Brown Vec local,global,globals[2],buffer; 13c4762a1bSJed Brown PetscScalar value; 14c4762a1bSJed Brown PetscViewer viewer; 15c4762a1bSJed Brown 16*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 17c4762a1bSJed Brown 185f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeCreate(PETSC_COMM_WORLD,&packer)); 195f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,8,1,1,NULL,&da1)); 205f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da1)); 215f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da1)); 225f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeAddDM(packer,da1)); 235f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,6,1,1,NULL,&da2)); 245f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da2)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da2)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeAddDM(packer,da2)); 27c4762a1bSJed Brown 285f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(packer,&global)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLocalVector(packer,&local)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLocalVector(packer,&buffer)); 31c4762a1bSJed Brown 325f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeGetAccessArray(packer,global,2,NULL,globals)); 33c4762a1bSJed Brown value = 1; 345f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(globals[0], value)); 35c4762a1bSJed Brown value = -1; 365f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(globals[1], value)); 375f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 38c4762a1bSJed Brown value = rank + 1; 395f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(globals[0], value)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(globals[1], value)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(DMCompositeRestoreAccessArray(packer,global,2,NULL,globals)); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* Test GlobalToLocal in insert mode */ 445f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(packer,global,INSERT_VALUES,local)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(packer,global,INSERT_VALUES,local)); 46c4762a1bSJed Brown 475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"\nLocal Vector: processor %d\n",rank)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&viewer)); 505f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(local,viewer)); 515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&viewer)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* Test LocalToGlobal in insert mode */ 565f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalBegin(packer,local,INSERT_VALUES,global)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalEnd(packer,local,INSERT_VALUES,global)); 58c4762a1bSJed Brown 595f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(global,PETSC_VIEWER_STDOUT_WORLD)); 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* Test LocalToLocal in insert mode */ 625f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToLocalBegin(packer,local,INSERT_VALUES,buffer)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToLocalEnd(packer,local,INSERT_VALUES,buffer)); 64c4762a1bSJed Brown 655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"\nLocal Vector: processor %d\n",rank)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&viewer)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(buffer,viewer)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&viewer)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 72c4762a1bSJed Brown 735f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&buffer)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&local)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&global)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&packer)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da2)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da1)); 79c4762a1bSJed Brown 80*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 81*b122ec5aSJacob Faibussowitsch return 0; 82c4762a1bSJed Brown } 83c4762a1bSJed Brown 84c4762a1bSJed Brown /*TEST 85c4762a1bSJed Brown 86c4762a1bSJed Brown test: 87c4762a1bSJed Brown nsize: 3 88c4762a1bSJed Brown 89c4762a1bSJed Brown TEST*/ 90