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