xref: /petsc/src/dm/tutorials/ex25.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
160c22052SBarry Smith 
260c22052SBarry Smith static char help[] = "Takes a patch of a large DMDA vector to one process.\n\n";
360c22052SBarry Smith 
460c22052SBarry Smith #include <petscdm.h>
560c22052SBarry Smith #include <petscdmda.h>
660c22052SBarry Smith #include <petscdmpatch.h>
760c22052SBarry Smith #include <petscsf.h>
860c22052SBarry Smith 
960c22052SBarry Smith typedef struct {
1060c22052SBarry Smith   PetscScalar x, y;
1160c22052SBarry Smith } Field;
1260c22052SBarry Smith 
13*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
14*d71ae5a4SJacob Faibussowitsch {
1560c22052SBarry Smith   Vec         xy, sxy;
1660c22052SBarry Smith   DM          da, sda = NULL;
1760c22052SBarry Smith   PetscSF     sf;
1860c22052SBarry Smith   PetscInt    m = 10, n = 10, dof = 2;
1960c22052SBarry Smith   MatStencil  lower = {0, 3, 2, 0}, upper = {0, 7, 8, 0}; /* These are in the order of the z, y, x, logical coordinates, the fourth entry is ignored */
2060c22052SBarry Smith   MPI_Comm    comm;
2160c22052SBarry Smith   PetscMPIInt rank;
2260c22052SBarry Smith 
23327415f7SBarry Smith   PetscFunctionBeginUser;
249566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
2560c22052SBarry Smith 
26a5b23f4aSJose E. Roman   /* create the large DMDA and set coordinates (which we will copy down to the small DA). */
279566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, m, n, PETSC_DECIDE, PETSC_DECIDE, dof, 1, 0, 0, &da));
289566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
299566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
309566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
3160c22052SBarry Smith   /* Just as a simple example we use the coordinates as the variables in the vectors we wish to examine. */
329566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(da, &xy));
3360c22052SBarry Smith   /* The vector entries are displayed in the "natural" ordering on the two dimensional grid; interlaced x and y with with the x variable increasing more rapidly than the y */
349566063dSJacob Faibussowitsch   PetscCall(VecView(xy, 0));
3560c22052SBarry Smith 
369566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
37dd400576SPatrick Sanan   if (rank == 0) comm = MPI_COMM_SELF;
3860c22052SBarry Smith   else comm = MPI_COMM_NULL;
3960c22052SBarry Smith 
409566063dSJacob Faibussowitsch   PetscCall(DMPatchZoom(da, lower, upper, comm, &sda, NULL, &sf));
41dd400576SPatrick Sanan   if (rank == 0) {
429566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(sda, &sxy));
4360c22052SBarry Smith   } else {
449566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF, 0, &sxy));
4560c22052SBarry Smith   }
4660c22052SBarry Smith   /*  A PetscSF can also be used as a VecScatter context */
479566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(sf, xy, sxy, INSERT_VALUES, SCATTER_FORWARD));
489566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(sf, xy, sxy, INSERT_VALUES, SCATTER_FORWARD));
4960c22052SBarry Smith   /* Only rank == 0 has the entries of the patch, so run code only at that rank */
50dd400576SPatrick Sanan   if (rank == 0) {
5160c22052SBarry Smith     Field       **vars;
5260c22052SBarry Smith     DMDALocalInfo info;
5360c22052SBarry Smith     PetscInt      i, j;
5460c22052SBarry Smith     PetscScalar   sum = 0;
5560c22052SBarry Smith 
5660c22052SBarry Smith     /* The vector entries of the patch are displayed in the "natural" ordering on the two grid; interlaced x and y with with the x variable increasing more rapidly */
579566063dSJacob Faibussowitsch     PetscCall(VecView(sxy, PETSC_VIEWER_STDOUT_SELF));
5860c22052SBarry Smith     /* Compute some trivial statistic of the coordinates */
599566063dSJacob Faibussowitsch     PetscCall(DMDAGetLocalInfo(sda, &info));
609566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(sda, sxy, &vars));
6160c22052SBarry Smith     /* Loop over the patch of the entire domain */
6260c22052SBarry Smith     for (j = info.ys; j < info.ys + info.ym; j++) {
63ad540459SPierre Jolivet       for (i = info.xs; i < info.xs + info.xm; i++) sum += vars[j][i].x;
6460c22052SBarry Smith     }
659566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "The sum of the x coordinates is %g\n", (double)PetscRealPart(sum)));
669566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(sda, sxy, &vars));
6760c22052SBarry Smith   }
6860c22052SBarry Smith 
699566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sxy));
709566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sf));
719566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&sda));
729566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
739566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
74b122ec5aSJacob Faibussowitsch   return 0;
7560c22052SBarry Smith }
7660c22052SBarry Smith 
7760c22052SBarry Smith /*TEST
7860c22052SBarry Smith 
7960c22052SBarry Smith    test:
8060c22052SBarry Smith 
8160c22052SBarry Smith    test:
8260c22052SBarry Smith      nsize: 4
8360c22052SBarry Smith      suffix: 2
8460c22052SBarry Smith 
8560c22052SBarry Smith TEST*/
86