Actual source code: ex16.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  2: static char help[] = "Tests DMComposite routines.\n\n";

  4: #include <petscdmredundant.h>
  5: #include <petscdm.h>
  6: #include <petscdmda.h>
  7: #include <petscdmcomposite.h>
  8: #include <petscpf.h>

 12: int main(int argc,char **argv)
 13: {
 14:   PetscErrorCode         ierr;
 15:   PetscInt               nredundant1 = 5,nredundant2 = 2,i;
 16:   ISLocalToGlobalMapping *ltog;
 17:   PetscMPIInt            rank,size;
 18:   DM                     packer;
 19:   Vec                    global,local1,local2,redundant1,redundant2;
 20:   PF                     pf;
 21:   DM                     da1,da2,dmred1,dmred2;
 22:   PetscScalar            *redundant1a,*redundant2a;
 23:   PetscViewer            sviewer;
 24:   PetscBool              gather_add = PETSC_FALSE;

 26:   PetscInitialize(&argc,&argv,(char*)0,help);
 27:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 28:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 30:   PetscOptionsGetBool(NULL,NULL,"-gather_add",&gather_add,NULL);

 32:   DMCompositeCreate(PETSC_COMM_WORLD,&packer);

 34:   DMRedundantCreate(PETSC_COMM_WORLD,0,nredundant1,&dmred1);
 35:   DMCreateLocalVector(dmred1,&redundant1);
 36:   DMCompositeAddDM(packer,dmred1);

 38:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,8,1,1,NULL,&da1);
 39:   DMCreateLocalVector(da1,&local1);
 40:   DMCompositeAddDM(packer,da1);

 42:   DMRedundantCreate(PETSC_COMM_WORLD,1%size,nredundant2,&dmred2);
 43:   DMCreateLocalVector(dmred2,&redundant2);
 44:   DMCompositeAddDM(packer,dmred2);

 46:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,6,1,1,NULL,&da2);
 47:   DMCreateLocalVector(da2,&local2);
 48:   DMCompositeAddDM(packer,da2);

 50:   DMCreateGlobalVector(packer,&global);
 51:   PFCreate(PETSC_COMM_WORLD,1,1,&pf);
 52:   PFSetType(pf,PFIDENTITY,NULL);
 53:   PFApplyVec(pf,NULL,global);
 54:   PFDestroy(&pf);
 55:   VecView(global,PETSC_VIEWER_STDOUT_WORLD);

 57:   DMCompositeScatter(packer,global,redundant1,local1,redundant2,local2);
 58:   PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);
 59:   PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"[%d] My part of redundant1 vector\n",rank);
 60:   PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);
 61:   VecView(redundant1,sviewer);
 62:   PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);
 63:   PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);
 64:   PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"[%d] My part of da1 vector\n",rank);
 65:   PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);
 66:   VecView(local1,sviewer);
 67:   PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);
 68:   PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);
 69:   PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"[%d] My part of redundant2 vector\n",rank);
 70:   PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);
 71:   VecView(redundant2,sviewer);
 72:   PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);
 73:   PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);
 74:   PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"[%d] My part of da2 vector\n",rank);
 75:   PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);
 76:   VecView(local2,sviewer);
 77:   PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);
 78:   PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);
 79:   PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);

 81:   VecGetArray(redundant1,&redundant1a);
 82:   VecGetArray(redundant2,&redundant2a);
 83:   for (i=0; i<nredundant1; i++) redundant1a[i] = (rank+2)*i;
 84:   for (i=0; i<nredundant2; i++) redundant2a[i] = (rank+10)*i;
 85:   VecRestoreArray(redundant1,&redundant1a);
 86:   VecRestoreArray(redundant2,&redundant2a);

 88:   DMCompositeGather(packer,global,gather_add ? ADD_VALUES : INSERT_VALUES,redundant1,local1,redundant2,local2);
 89:   VecView(global,PETSC_VIEWER_STDOUT_WORLD);

 91:   /* get the global numbering for each subvector element */
 92:   DMCompositeGetISLocalToGlobalMappings(packer,&ltog);

 94:   PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"Local to global mapping of redundant1 vector\n");
 95:   ISLocalToGlobalMappingView(ltog[0],PETSC_VIEWER_STDOUT_WORLD);
 96:   PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"Local to global mapping of local1 vector\n");
 97:   ISLocalToGlobalMappingView(ltog[1],PETSC_VIEWER_STDOUT_WORLD);
 98:   PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"Local to global mapping of redundant2 vector\n");
 99:   ISLocalToGlobalMappingView(ltog[2],PETSC_VIEWER_STDOUT_WORLD);
100:   PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"Local to global mapping of local2 vector\n");
101:   ISLocalToGlobalMappingView(ltog[3],PETSC_VIEWER_STDOUT_WORLD);

103:   for (i=0; i<4; i++) {ISLocalToGlobalMappingDestroy(&ltog[i]);}
104:   PetscFree(ltog);

106:   DMDestroy(&da1);
107:   DMDestroy(&dmred1);
108:   DMDestroy(&dmred2);
109:   DMDestroy(&da2);
110:   VecDestroy(&redundant1);
111:   VecDestroy(&redundant2);
112:   VecDestroy(&local1);
113:   VecDestroy(&local2);
114:   VecDestroy(&global);
115:   DMDestroy(&packer);
116:   PetscFinalize();
117:   return 0;
118: }