Actual source code: ex51.c
 
   petsc-3.12.4 2020-02-04
   
  1: static char help[] = "Tests DMDAVecGetArrayDOF()\n";
  2:  #include <petscdm.h>
  3:  #include <petscdmda.h>
  5: int main(int argc, char *argv[])
  6: {
  7:   DM             da, daX, daY;
  8:   DMDALocalInfo  info;
  9:   MPI_Comm       commX, commY;
 10:   Vec            basisX, basisY;
 11:   PetscScalar    **arrayX, **arrayY;
 12:   const PetscInt *lx, *ly;
 13:   PetscInt       M     = 3, N = 3;
 14:   PetscInt       p     = 1;
 15:   PetscInt       numGP = 3;
 16:   PetscInt       dof   = 2*(p+1)*numGP;
 17:   PetscMPIInt    rank, subsize, subrank;
 20:   PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr;
 21:   MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
 22:   /* Create 2D DMDA */
 23:   DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,N,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);
 24:   DMSetFromOptions(da);
 25:   DMSetUp(da);
 26:   /* Create 1D DMDAs along two directions. */
 27:   DMDAGetOwnershipRanges(da, &lx, &ly, NULL);
 28:   DMDAGetLocalInfo(da, &info);
 29:   /* Partitioning in the X direction makes a subcomm extending in the Y direction and vice-versa. */
 30:   DMDAGetProcessorSubsets(da, DMDA_X, &commY);
 31:   DMDAGetProcessorSubsets(da, DMDA_Y, &commX);
 32:   MPI_Comm_size(commX, &subsize);
 33:   MPI_Comm_rank(commX, &subrank);
 34:   PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]X subrank: %d subsize: %d\n", rank, subrank, subsize);
 35:   PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);
 36:   MPI_Comm_size(commY, &subsize);
 37:   MPI_Comm_rank(commY, &subrank);
 38:   PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Y subrank: %d subsize: %d\n", rank, subrank, subsize);
 39:   PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);
 40:   DMDACreate1d(commX, DM_BOUNDARY_NONE, info.mx, dof, 1, lx, &daX);
 41:   DMSetUp(daX);
 42:   DMDACreate1d(commY, DM_BOUNDARY_NONE, info.my, dof, 1, ly, &daY);
 43:   DMSetUp(daY);
 44:   /* Create 1D vectors for basis functions */
 45:   DMGetGlobalVector(daX, &basisX);
 46:   DMGetGlobalVector(daY, &basisY);
 47:   /* Extract basis functions */
 48:   DMDAVecGetArrayDOF(daX, basisX, &arrayX);
 49:   DMDAVecGetArrayDOF(daY, basisY, &arrayY);
 50:   /*arrayX[i][ndof]; */
 51:   /*arrayY[j][ndof]; */
 52:   DMDAVecRestoreArrayDOF(daX, basisX, &arrayX);
 53:   DMDAVecRestoreArrayDOF(daY, basisY, &arrayY);
 54:   /* Return basis vectors */
 55:   DMRestoreGlobalVector(daX, &basisX);
 56:   DMRestoreGlobalVector(daY, &basisY);
 57:   /* Cleanup */
 58:   MPI_Comm_free(&commX);
 59:   MPI_Comm_free(&commY);
 60:   DMDestroy(&daX);
 61:   DMDestroy(&daY);
 62:   DMDestroy(&da);
 63:   PetscFinalize();
 64:   return ierr;
 65: }
 68: /*TEST
 70:    test:
 71:       nsize: 2
 73: TEST*/