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