Actual source code: ex51.c
petsc-3.4.5 2014-06-29
1: static char help[] = "Tests DMDAVecGetArrayDOF()\n";
2: #include <petscdmda.h>
4: int main(int argc, char *argv[])
5: {
6: DM da, daX, daY;
7: DMDALocalInfo info;
8: MPI_Comm commX, commY;
9: Vec basisX, basisY;
10: PetscScalar **arrayX, **arrayY;
11: const PetscInt *lx, *ly;
12: PetscInt M = 3, N = 3;
13: PetscInt p = 1;
14: PetscInt numGP = 3;
15: PetscInt dof = 2*(p+1)*numGP;
16: PetscMPIInt rank, subsize, subrank;
19: PetscInitialize(&argc,&argv,0,help);
20: MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
21: /* Create 2D DMDA */
22: DMDACreate2d(PETSC_COMM_WORLD, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE, DMDA_STENCIL_STAR, M, N, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da);
23: /* Create 1D DMDAs along two directions */
24: DMDAGetOwnershipRanges(da, &lx, &ly, NULL);
25: DMDAGetLocalInfo(da, &info);
26: DMDAGetProcessorSubsets(da, DMDA_X, &commX);
27: DMDAGetProcessorSubsets(da, DMDA_Y, &commY);
28: MPI_Comm_size(commX, &subsize);
29: MPI_Comm_rank(commX, &subrank);
30: PetscPrintf(PETSC_COMM_SELF, "[%d]X subrank: %d subsize: %d\n", rank, subrank, subsize);
31: MPI_Comm_size(commY, &subsize);
32: MPI_Comm_rank(commY, &subrank);
33: PetscPrintf(PETSC_COMM_SELF, "[%d]Y subrank: %d subsize: %d\n", rank, subrank, subsize);
34: DMDACreate1d(commX, DMDA_BOUNDARY_NONE, M, dof, 1, lx, &daX);
35: DMDACreate1d(commY, DMDA_BOUNDARY_NONE, N, dof, 1, ly, &daY);
36: /* Create 1D vectors for basis functions */
37: DMGetGlobalVector(daX, &basisX);
38: DMGetGlobalVector(daY, &basisY);
39: /* Extract basis functions */
40: DMDAVecGetArrayDOF(daX, basisX, &arrayX);
41: DMDAVecGetArrayDOF(daY, basisY, &arrayY);
42: /*arrayX[i][ndof]; */
43: /*arrayY[j][ndof]; */
44: DMDAVecRestoreArrayDOF(daX, basisX, &arrayX);
45: DMDAVecRestoreArrayDOF(daY, basisY, &arrayY);
46: /* Return basis vectors */
47: DMRestoreGlobalVector(daX, &basisX);
48: DMRestoreGlobalVector(daY, &basisY);
49: /* Cleanup */
50: DMDestroy(&daX);
51: DMDestroy(&daY);
52: DMDestroy(&da);
53: PetscFinalize();
54: return 0;
55: }