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