Actual source code: ex44.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

  2: static char help[] = "Tests various DMComposite routines.\n\n";

  4:  #include <petscdm.h>
  5:  #include <petscdmda.h>
  6:  #include <petscdmcomposite.h>

  8: int main(int argc,char **argv)
  9: {
 10:   PetscMPIInt            rank;
 11:   PetscErrorCode         ierr;
 12:   DM                     da1,da2,packer;
 13:   Vec                    local,global,globals[2],buffer;
 14:   PetscScalar            value;
 15:   PetscViewer            viewer;

 17:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;

 19:   DMCompositeCreate(PETSC_COMM_WORLD,&packer);
 20:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,8,1,1,NULL,&da1);
 21:   DMSetFromOptions(da1);
 22:   DMSetUp(da1);
 23:   DMCompositeAddDM(packer,da1);
 24:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,6,1,1,NULL,&da2);
 25:   DMSetFromOptions(da2);
 26:   DMSetUp(da2);
 27:   DMCompositeAddDM(packer,da2);

 29:   DMCreateGlobalVector(packer,&global);
 30:   DMCreateLocalVector(packer,&local);
 31:   DMCreateLocalVector(packer,&buffer);

 33:   DMCompositeGetAccessArray(packer,global,2,NULL,globals);
 34:   value = 1;
 35:   VecSet(globals[0], value);
 36:   value = -1;
 37:   VecSet(globals[1], value);
 38:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 39:   value = rank + 1;
 40:   VecScale(globals[0], value);
 41:   VecScale(globals[1], value);
 42:   DMCompositeRestoreAccessArray(packer,global,2,NULL,globals);

 44:   /* Test GlobalToLocal in insert mode */
 45:   DMGlobalToLocalBegin(packer,global,INSERT_VALUES,local);
 46:   DMGlobalToLocalEnd(packer,global,INSERT_VALUES,local);

 48:   PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);
 49:   PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"\nLocal Vector: processor %d\n",rank);
 50:   PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&viewer);
 51:   VecView(local,viewer);
 52:   PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&viewer);
 53:   PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);
 54:   PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);

 56:   /* Test LocalToGlobal in insert mode */
 57:   DMLocalToGlobalBegin(packer,local,INSERT_VALUES,global);
 58:   DMLocalToGlobalEnd(packer,local,INSERT_VALUES,global);

 60:   VecView(global,PETSC_VIEWER_STDOUT_WORLD);

 62:   /* Test LocalToLocal in insert mode */
 63:   DMLocalToLocalBegin(packer,local,INSERT_VALUES,buffer);
 64:   DMLocalToLocalEnd(packer,local,INSERT_VALUES,buffer);

 66:   PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);
 67:   PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"\nLocal Vector: processor %d\n",rank);
 68:   PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&viewer);
 69:   VecView(buffer,viewer);
 70:   PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&viewer);
 71:   PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);
 72:   PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);

 74:   VecDestroy(&buffer);
 75:   VecDestroy(&local);
 76:   VecDestroy(&global);
 77:   DMDestroy(&packer);
 78:   DMDestroy(&da2);
 79:   DMDestroy(&da1);

 81:   PetscFinalize();

 83:   return ierr;
 84: }


 87: /*TEST

 89:    test:
 90:       nsize: 3

 92: TEST*/