xref: /petsc/src/dm/tutorials/ex6.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*
5c4762a1bSJed Brown      Demonstrates using DM_BOUNDARY_GHOSTED how to handle a rotated boundary conditions where one edge
6c4762a1bSJed Brown     is connected to its immediate neighbor
7c4762a1bSJed Brown 
8c4762a1bSJed Brown     Consider the domain (with natural numbering)
9c4762a1bSJed Brown 
10c4762a1bSJed Brown      6   7   8
11c4762a1bSJed Brown      3   4   5
12c4762a1bSJed Brown      0   1   2
13c4762a1bSJed Brown 
14c4762a1bSJed Brown     The ghost points along the bottom (directly below the three columns above) should be 0 3 and 6
15c4762a1bSJed Brown     while the ghost points along the left side should be 0 1 2
16c4762a1bSJed Brown 
17c4762a1bSJed Brown     Note that the ghosted local vectors extend in both the x and y directions so, for example if we have a
18c4762a1bSJed Brown     single MPI process the ghosted vector has (in the original natural numbering)
19c4762a1bSJed Brown 
20c4762a1bSJed Brown      x  x  x  x  x
21c4762a1bSJed Brown      2  6  7  8  x
22c4762a1bSJed Brown      1  3  4  5  x
23c4762a1bSJed Brown      0  0  1  2  x
24c4762a1bSJed Brown      x  0  3  6  x
25c4762a1bSJed Brown 
26c4762a1bSJed Brown     where x indicates a location that is not updated by the communication and should be used.
27c4762a1bSJed Brown 
28c4762a1bSJed Brown     For this to make sense the number of grid points in the x and y directions must be the same
29c4762a1bSJed Brown 
30c4762a1bSJed Brown     This ghost point mapping was suggested by: Wenbo Zhao <zhaowenbo.npic@gmail.com>
31c4762a1bSJed Brown */
32c4762a1bSJed Brown 
33c4762a1bSJed Brown #include <petscdm.h>
34c4762a1bSJed Brown #include <petscdmda.h>
35c4762a1bSJed Brown 
36c4762a1bSJed Brown int main(int argc,char **argv)
37c4762a1bSJed Brown {
38c4762a1bSJed Brown   PetscInt         M = 6;
39c4762a1bSJed Brown   PetscErrorCode   ierr;
40c4762a1bSJed Brown   DM               da;
41c4762a1bSJed Brown   Vec              local,global,natural;
42c4762a1bSJed Brown   PetscInt         i,start,end,*ifrom,x,y,xm,ym;
43c4762a1bSJed Brown   PetscScalar      *xnatural;
44c4762a1bSJed Brown   IS               from,to;
45c4762a1bSJed Brown   AO               ao;
46c4762a1bSJed Brown   VecScatter       scatter1,scatter2;
47c4762a1bSJed Brown   PetscViewer      subviewer;
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   /* Create distributed array and get vectors */
52*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_GHOSTED,DM_BOUNDARY_GHOSTED,DMDA_STENCIL_STAR,M,M,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da));
53*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
54*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&global));
55*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateLocalVector(da,&local));
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   /* construct global to local scatter for the left side of the domain to the ghost on the bottom */
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&x,&y,NULL,&xm,&ym,NULL));
59c4762a1bSJed Brown   if (!y) { /* only processes on the bottom of the domain fill up the ghost locations */
60*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF,xm,1,1,&to));
61c4762a1bSJed Brown   } else {
62*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF,0,0,0,&to));
63c4762a1bSJed Brown   }
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(xm,&ifrom));
65c4762a1bSJed Brown   for (i=x;i<x+xm;i++) {
66c4762a1bSJed Brown     ifrom[i-x] = M*i;
67c4762a1bSJed Brown   }
68*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetAO(da,&ao));
69*5f80ce2aSJacob Faibussowitsch   CHKERRQ(AOApplicationToPetsc(ao,xm,ifrom));
70c4762a1bSJed Brown   if (!y) {
71*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,xm,ifrom,PETSC_OWN_POINTER,&from));
72c4762a1bSJed Brown   } else {
73*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(ifrom));
74*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,0,NULL,PETSC_COPY_VALUES,&from));
75c4762a1bSJed Brown   }
76*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(global,from,local,to,&scatter1));
77*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&to));
78*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&from));
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /* construct global to local scatter for the bottom side of the domain to the ghost on the right */
81c4762a1bSJed Brown   if (!x) { /* only processes on the left side of the domain fill up the ghost locations */
82*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF,ym,xm+2,xm+2,&to));
83c4762a1bSJed Brown   } else {
84*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF,0,0,0,&to));
85c4762a1bSJed Brown   }
86*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(ym,&ifrom));
87c4762a1bSJed Brown   for (i=y;i<y+ym;i++) {
88c4762a1bSJed Brown     ifrom[i-y] = i;
89c4762a1bSJed Brown   }
90*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetAO(da,&ao));
91*5f80ce2aSJacob Faibussowitsch   CHKERRQ(AOApplicationToPetsc(ao,ym,ifrom));
92c4762a1bSJed Brown   if (!x) {
93*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,ym,ifrom,PETSC_OWN_POINTER,&from));
94c4762a1bSJed Brown   } else {
95*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(ifrom));
96*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,0,NULL,PETSC_COPY_VALUES,&from));
97c4762a1bSJed Brown   }
98*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(global,from,local,to,&scatter2));
99*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&to));
100*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&from));
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   /*
103c4762a1bSJed Brown      fill the global vector with the natural global numbering for each local entry
104c4762a1bSJed Brown      this is only done for testing purposes since it is easy to see if the scatter worked correctly
105c4762a1bSJed Brown   */
106*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreateNaturalVector(da,&natural));
107*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetOwnershipRange(natural,&start,&end));
108*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(natural,&xnatural));
109c4762a1bSJed Brown   for (i=start; i<end; i++) {
110c4762a1bSJed Brown     xnatural[i-start] = i;
111c4762a1bSJed Brown   }
112*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(natural,&xnatural));
113*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDANaturalToGlobalBegin(da,natural,INSERT_VALUES,global));
114*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDANaturalToGlobalEnd(da,natural,INSERT_VALUES,global));
115*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&natural));
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   /* scatter from global to local */
118*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(scatter1,global,local,INSERT_VALUES,SCATTER_FORWARD));
119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(scatter1,global,local,INSERT_VALUES,SCATTER_FORWARD));
120*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(scatter2,global,local,INSERT_VALUES,SCATTER_FORWARD));
121*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(scatter2,global,local,INSERT_VALUES,SCATTER_FORWARD));
122c4762a1bSJed Brown   /*
123c4762a1bSJed Brown      normally here you would also call
124*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,global,INSERT_VALUES,local));
125*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,global,INSERT_VALUES,local));
126c4762a1bSJed Brown     to update all the interior ghost cells between neighboring processes.
127c4762a1bSJed Brown     We don't do it here since this is only a test of "special" ghost points.
128c4762a1bSJed Brown   */
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   /* view each local ghosted vector */
131*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&subviewer));
132*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(local,subviewer));
133*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&subviewer));
134c4762a1bSJed Brown 
135*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&scatter1));
136*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&scatter2));
137*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&local));
138*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&global));
139*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
140c4762a1bSJed Brown   ierr = PetscFinalize();
141c4762a1bSJed Brown   return ierr;
142c4762a1bSJed Brown }
143c4762a1bSJed Brown 
144c4762a1bSJed Brown /*TEST
145c4762a1bSJed Brown 
146c4762a1bSJed Brown    test:
147c4762a1bSJed Brown 
148c4762a1bSJed Brown    test:
149c4762a1bSJed Brown       suffix: 2
150c4762a1bSJed Brown       nsize: 2
151c4762a1bSJed Brown 
152c4762a1bSJed Brown    test:
153c4762a1bSJed Brown       suffix: 4
154c4762a1bSJed Brown       nsize: 4
155c4762a1bSJed Brown 
156c4762a1bSJed Brown    test:
157c4762a1bSJed Brown       suffix: 9
158c4762a1bSJed Brown       nsize: 9
159c4762a1bSJed Brown 
160c4762a1bSJed Brown TEST*/
161