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