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