xref: /petsc/src/dm/tutorials/ex51.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
1 static char help[] = "Tests DMDAVecGetArrayDOF()\n";
2 #include <petscdm.h>
3 #include <petscdmda.h>
4 
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;
18   PetscErrorCode ierr;
19 
20   ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr;
21   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
22   /* Create 2D DMDA */
23   CHKERRQ(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   CHKERRQ(DMSetFromOptions(da));
25   CHKERRQ(DMSetUp(da));
26   /* Create 1D DMDAs along two directions. */
27   CHKERRQ(DMDAGetOwnershipRanges(da, &lx, &ly, NULL));
28   CHKERRQ(DMDAGetLocalInfo(da, &info));
29   /* Partitioning in the X direction makes a subcomm extending in the Y direction and vice-versa. */
30   CHKERRQ(DMDAGetProcessorSubsets(da, DM_X, &commY));
31   CHKERRQ(DMDAGetProcessorSubsets(da, DM_Y, &commX));
32   CHKERRMPI(MPI_Comm_size(commX, &subsize));
33   CHKERRMPI(MPI_Comm_rank(commX, &subrank));
34   CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]X subrank: %d subsize: %d\n", rank, subrank, subsize));
35   CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
36   CHKERRMPI(MPI_Comm_size(commY, &subsize));
37   CHKERRMPI(MPI_Comm_rank(commY, &subrank));
38   CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Y subrank: %d subsize: %d\n", rank, subrank, subsize));
39   CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
40   CHKERRQ(DMDACreate1d(commX, DM_BOUNDARY_NONE, info.mx, dof, 1, lx, &daX));
41   CHKERRQ(DMSetUp(daX));
42   CHKERRQ(DMDACreate1d(commY, DM_BOUNDARY_NONE, info.my, dof, 1, ly, &daY));
43   CHKERRQ(DMSetUp(daY));
44   /* Create 1D vectors for basis functions */
45   CHKERRQ(DMGetGlobalVector(daX, &basisX));
46   CHKERRQ(DMGetGlobalVector(daY, &basisY));
47   /* Extract basis functions */
48   CHKERRQ(DMDAVecGetArrayDOF(daX, basisX, &arrayX));
49   CHKERRQ(DMDAVecGetArrayDOF(daY, basisY, &arrayY));
50   /*arrayX[i][ndof]; */
51   /*arrayY[j][ndof]; */
52   CHKERRQ(DMDAVecRestoreArrayDOF(daX, basisX, &arrayX));
53   CHKERRQ(DMDAVecRestoreArrayDOF(daY, basisY, &arrayY));
54   /* Return basis vectors */
55   CHKERRQ(DMRestoreGlobalVector(daX, &basisX));
56   CHKERRQ(DMRestoreGlobalVector(daY, &basisY));
57   /* Cleanup */
58   CHKERRMPI(MPI_Comm_free(&commX));
59   CHKERRMPI(MPI_Comm_free(&commY));
60   CHKERRQ(DMDestroy(&daX));
61   CHKERRQ(DMDestroy(&daY));
62   CHKERRQ(DMDestroy(&da));
63   ierr = PetscFinalize();
64   return ierr;
65 }
66 
67 /*TEST
68 
69    test:
70       nsize: 2
71 
72 TEST*/
73