xref: /petsc/src/dm/tutorials/ex25.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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 
1360c22052SBarry Smith int main(int argc,char **argv)
1460c22052SBarry Smith {
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 
23*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
2460c22052SBarry Smith 
25a5b23f4aSJose E. Roman   /* create the large DMDA and set coordinates (which we will copy down to the small DA). */
265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,m,n,PETSC_DECIDE,PETSC_DECIDE,dof,1,0,0,&da));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,1.0));
3060c22052SBarry Smith   /* Just as a simple example we use the coordinates as the variables in the vectors we wish to examine. */
315f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinates(da,&xy));
3260c22052SBarry 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 */
335f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(xy,0));
3460c22052SBarry Smith 
355f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
36dd400576SPatrick Sanan   if (rank == 0) comm = MPI_COMM_SELF;
3760c22052SBarry Smith   else comm = MPI_COMM_NULL;
3860c22052SBarry Smith 
395f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPatchZoom(da,lower,upper,comm,&sda, NULL,&sf));
40dd400576SPatrick Sanan   if (rank == 0) {
415f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCreateGlobalVector(sda,&sxy));
4260c22052SBarry Smith   } else {
435f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,0,&sxy));
4460c22052SBarry Smith   }
4560c22052SBarry Smith   /*  A PetscSF can also be used as a VecScatter context */
465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(sf,xy,sxy,INSERT_VALUES,SCATTER_FORWARD));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(sf,xy,sxy,INSERT_VALUES,SCATTER_FORWARD));
4860c22052SBarry Smith   /* Only rank == 0 has the entries of the patch, so run code only at that rank */
49dd400576SPatrick Sanan   if (rank == 0) {
5060c22052SBarry Smith     Field         **vars;
5160c22052SBarry Smith     DMDALocalInfo info;
5260c22052SBarry Smith     PetscInt      i,j;
5360c22052SBarry Smith     PetscScalar   sum = 0;
5460c22052SBarry Smith 
5560c22052SBarry 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 */
565f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(sxy,PETSC_VIEWER_STDOUT_SELF));
5760c22052SBarry Smith     /* Compute some trivial statistic of the coordinates */
585f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAGetLocalInfo(sda,&info));
595f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArray(sda,sxy,&vars));
6060c22052SBarry Smith     /* Loop over the patch of the entire domain */
6160c22052SBarry Smith     for (j=info.ys; j<info.ys+info.ym; j++) {
6260c22052SBarry Smith       for (i=info.xs; i<info.xs+info.xm; i++) {
6360c22052SBarry Smith         sum += vars[j][i].x;
6460c22052SBarry Smith       }
6560c22052SBarry Smith     }
665f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"The sum of the x coordinates is %g\n",(double)PetscRealPart(sum)));
675f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArray(sda,sxy,&vars));
6860c22052SBarry Smith   }
6960c22052SBarry Smith 
705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&sxy));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFDestroy(&sf));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&sda));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
74*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
75*b122ec5aSJacob Faibussowitsch   return 0;
7660c22052SBarry Smith }
7760c22052SBarry Smith 
7860c22052SBarry Smith /*TEST
7960c22052SBarry Smith 
8060c22052SBarry Smith    test:
8160c22052SBarry Smith 
8260c22052SBarry Smith    test:
8360c22052SBarry Smith      nsize: 4
8460c22052SBarry Smith      suffix: 2
8560c22052SBarry Smith 
8660c22052SBarry Smith TEST*/
87