xref: /petsc/src/dm/tutorials/ex6.c (revision 9566063d113dddea24716c546802770db7481bc0)
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   DM               da;
40c4762a1bSJed Brown   Vec              local,global,natural;
41c4762a1bSJed Brown   PetscInt         i,start,end,*ifrom,x,y,xm,ym;
42c4762a1bSJed Brown   PetscScalar      *xnatural;
43c4762a1bSJed Brown   IS               from,to;
44c4762a1bSJed Brown   AO               ao;
45c4762a1bSJed Brown   VecScatter       scatter1,scatter2;
46c4762a1bSJed Brown   PetscViewer      subviewer;
47c4762a1bSJed Brown 
48*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   /* Create distributed array and get vectors */
51*9566063dSJacob Faibussowitsch   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));
52*9566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
53*9566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da,&global));
54*9566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(da,&local));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   /* construct global to local scatter for the left side of the domain to the ghost on the bottom */
57*9566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&x,&y,NULL,&xm,&ym,NULL));
58c4762a1bSJed Brown   if (!y) { /* only processes on the bottom of the domain fill up the ghost locations */
59*9566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,xm,1,1,&to));
60c4762a1bSJed Brown   } else {
61*9566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,0,0,0,&to));
62c4762a1bSJed Brown   }
63*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(xm,&ifrom));
64c4762a1bSJed Brown   for (i=x;i<x+xm;i++) {
65c4762a1bSJed Brown     ifrom[i-x] = M*i;
66c4762a1bSJed Brown   }
67*9566063dSJacob Faibussowitsch   PetscCall(DMDAGetAO(da,&ao));
68*9566063dSJacob Faibussowitsch   PetscCall(AOApplicationToPetsc(ao,xm,ifrom));
69c4762a1bSJed Brown   if (!y) {
70*9566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,xm,ifrom,PETSC_OWN_POINTER,&from));
71c4762a1bSJed Brown   } else {
72*9566063dSJacob Faibussowitsch     PetscCall(PetscFree(ifrom));
73*9566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,0,NULL,PETSC_COPY_VALUES,&from));
74c4762a1bSJed Brown   }
75*9566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(global,from,local,to,&scatter1));
76*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&to));
77*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&from));
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   /* construct global to local scatter for the bottom side of the domain to the ghost on the right */
80c4762a1bSJed Brown   if (!x) { /* only processes on the left side of the domain fill up the ghost locations */
81*9566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,ym,xm+2,xm+2,&to));
82c4762a1bSJed Brown   } else {
83*9566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,0,0,0,&to));
84c4762a1bSJed Brown   }
85*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ym,&ifrom));
86c4762a1bSJed Brown   for (i=y;i<y+ym;i++) {
87c4762a1bSJed Brown     ifrom[i-y] = i;
88c4762a1bSJed Brown   }
89*9566063dSJacob Faibussowitsch   PetscCall(DMDAGetAO(da,&ao));
90*9566063dSJacob Faibussowitsch   PetscCall(AOApplicationToPetsc(ao,ym,ifrom));
91c4762a1bSJed Brown   if (!x) {
92*9566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,ym,ifrom,PETSC_OWN_POINTER,&from));
93c4762a1bSJed Brown   } else {
94*9566063dSJacob Faibussowitsch     PetscCall(PetscFree(ifrom));
95*9566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,0,NULL,PETSC_COPY_VALUES,&from));
96c4762a1bSJed Brown   }
97*9566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(global,from,local,to,&scatter2));
98*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&to));
99*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&from));
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   /*
102c4762a1bSJed Brown      fill the global vector with the natural global numbering for each local entry
103c4762a1bSJed Brown      this is only done for testing purposes since it is easy to see if the scatter worked correctly
104c4762a1bSJed Brown   */
105*9566063dSJacob Faibussowitsch   PetscCall(DMDACreateNaturalVector(da,&natural));
106*9566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(natural,&start,&end));
107*9566063dSJacob Faibussowitsch   PetscCall(VecGetArray(natural,&xnatural));
108c4762a1bSJed Brown   for (i=start; i<end; i++) {
109c4762a1bSJed Brown     xnatural[i-start] = i;
110c4762a1bSJed Brown   }
111*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(natural,&xnatural));
112*9566063dSJacob Faibussowitsch   PetscCall(DMDANaturalToGlobalBegin(da,natural,INSERT_VALUES,global));
113*9566063dSJacob Faibussowitsch   PetscCall(DMDANaturalToGlobalEnd(da,natural,INSERT_VALUES,global));
114*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&natural));
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   /* scatter from global to local */
117*9566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter1,global,local,INSERT_VALUES,SCATTER_FORWARD));
118*9566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter1,global,local,INSERT_VALUES,SCATTER_FORWARD));
119*9566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter2,global,local,INSERT_VALUES,SCATTER_FORWARD));
120*9566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter2,global,local,INSERT_VALUES,SCATTER_FORWARD));
121c4762a1bSJed Brown   /*
122c4762a1bSJed Brown      normally here you would also call
123*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,global,INSERT_VALUES,local));
124*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,global,INSERT_VALUES,local));
125c4762a1bSJed Brown     to update all the interior ghost cells between neighboring processes.
126c4762a1bSJed Brown     We don't do it here since this is only a test of "special" ghost points.
127c4762a1bSJed Brown   */
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /* view each local ghosted vector */
130*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&subviewer));
131*9566063dSJacob Faibussowitsch   PetscCall(VecView(local,subviewer));
132*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&subviewer));
133c4762a1bSJed Brown 
134*9566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scatter1));
135*9566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scatter2));
136*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&local));
137*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&global));
138*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
139*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
140b122ec5aSJacob Faibussowitsch   return 0;
141c4762a1bSJed Brown }
142c4762a1bSJed Brown 
143c4762a1bSJed Brown /*TEST
144c4762a1bSJed Brown 
145c4762a1bSJed Brown    test:
146c4762a1bSJed Brown 
147c4762a1bSJed Brown    test:
148c4762a1bSJed Brown       suffix: 2
149c4762a1bSJed Brown       nsize: 2
150c4762a1bSJed Brown 
151c4762a1bSJed Brown    test:
152c4762a1bSJed Brown       suffix: 4
153c4762a1bSJed Brown       nsize: 4
154c4762a1bSJed Brown 
155c4762a1bSJed Brown    test:
156c4762a1bSJed Brown       suffix: 9
157c4762a1bSJed Brown       nsize: 9
158c4762a1bSJed Brown 
159c4762a1bSJed Brown TEST*/
160