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