xref: /petsc/src/dm/tutorials/ex51.c (revision 061e922f3926be00487707c73b78dd3d40309129)
1c4762a1bSJed Brown static char help[] = "Tests DMDAVecGetArrayDOF()\n";
2c4762a1bSJed Brown #include <petscdm.h>
3c4762a1bSJed Brown #include <petscdmda.h>
4c4762a1bSJed Brown 
main(int argc,char * argv[])5*d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[])
6*d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown   DM              da, daX, daY;
8c4762a1bSJed Brown   DMDALocalInfo   info;
9c4762a1bSJed Brown   MPI_Comm        commX, commY;
10c4762a1bSJed Brown   Vec             basisX, basisY;
11c4762a1bSJed Brown   PetscScalar   **arrayX, **arrayY;
12c4762a1bSJed Brown   const PetscInt *lx, *ly;
13c4762a1bSJed Brown   PetscInt        M = 3, N = 3;
14c4762a1bSJed Brown   PetscInt        p     = 1;
15c4762a1bSJed Brown   PetscInt        numGP = 3;
16c4762a1bSJed Brown   PetscInt        dof   = 2 * (p + 1) * numGP;
17c4762a1bSJed Brown   PetscMPIInt     rank, subsize, subrank;
18c4762a1bSJed Brown 
19327415f7SBarry Smith   PetscFunctionBeginUser;
209566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, 0, help));
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
22c4762a1bSJed Brown   /* Create 2D DMDA */
239566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, M, N, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da));
249566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
259566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
26c4762a1bSJed Brown   /* Create 1D DMDAs along two directions. */
279566063dSJacob Faibussowitsch   PetscCall(DMDAGetOwnershipRanges(da, &lx, &ly, NULL));
289566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
29c4762a1bSJed Brown   /* Partitioning in the X direction makes a subcomm extending in the Y direction and vice-versa. */
309566063dSJacob Faibussowitsch   PetscCall(DMDAGetProcessorSubsets(da, DM_X, &commY));
319566063dSJacob Faibussowitsch   PetscCall(DMDAGetProcessorSubsets(da, DM_Y, &commX));
329566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(commX, &subsize));
339566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(commX, &subrank));
349566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]X subrank: %d subsize: %d\n", rank, subrank, subsize));
359566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
369566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(commY, &subsize));
379566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(commY, &subrank));
389566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Y subrank: %d subsize: %d\n", rank, subrank, subsize));
399566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
409566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(commX, DM_BOUNDARY_NONE, info.mx, dof, 1, lx, &daX));
419566063dSJacob Faibussowitsch   PetscCall(DMSetUp(daX));
429566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(commY, DM_BOUNDARY_NONE, info.my, dof, 1, ly, &daY));
439566063dSJacob Faibussowitsch   PetscCall(DMSetUp(daY));
44c4762a1bSJed Brown   /* Create 1D vectors for basis functions */
459566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(daX, &basisX));
469566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(daY, &basisY));
47c4762a1bSJed Brown   /* Extract basis functions */
489566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayDOF(daX, basisX, &arrayX));
499566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayDOF(daY, basisY, &arrayY));
50c4762a1bSJed Brown   /*arrayX[i][ndof]; */
51c4762a1bSJed Brown   /*arrayY[j][ndof]; */
529566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayDOF(daX, basisX, &arrayX));
539566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayDOF(daY, basisY, &arrayY));
54c4762a1bSJed Brown   /* Return basis vectors */
559566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(daX, &basisX));
569566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(daY, &basisY));
57c4762a1bSJed Brown   /* Cleanup */
589566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free(&commX));
599566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free(&commY));
609566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&daX));
619566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&daY));
629566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
639566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
64b122ec5aSJacob Faibussowitsch   return 0;
65c4762a1bSJed Brown }
66c4762a1bSJed Brown 
67c4762a1bSJed Brown /*TEST
68c4762a1bSJed Brown 
69c4762a1bSJed Brown    test:
70c4762a1bSJed Brown       nsize: 2
71c4762a1bSJed Brown 
72c4762a1bSJed Brown TEST*/
73