xref: /petsc/src/dm/impls/da/daindex.c (revision 021e0df48735593a4b10a8d433fd5ac359e13c2e)
1 
2 /*
3   Code for manipulating distributed regular arrays in parallel.
4 */
5 
6 #include <petsc-private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
7 
8 #undef __FUNCT__
9 #define __FUNCT__ "DMDAGetGlobalIndices"
10 /*@C
11    DMDAGetGlobalIndices - Returns the global node number of all local nodes,
12    including ghost nodes.
13 
14    Not Collective
15 
16    Input Parameter:
17 .  da - the distributed array
18 
19    Output Parameters:
20 +  n - the number of local elements, including ghost nodes (or NULL)
21 -  idx - the global indices
22 
23    Level: intermediate
24 
25    Note:
26    For DMDA_STENCIL_STAR stencils the inactive corner ghost nodes are also included
27    in the list of local indices (even though those nodes are not updated
28    during calls to DMDAXXXToXXX().
29 
30    Essentially the same data is returned in the form of a local-to-global mapping
31    with the routine DMDAGetISLocalToGlobalMapping(), that is the recommended interface.
32 
33    You must call DMDARestoreGlobalIndices() after you are finished using the indices
34 
35    Fortran Note:
36    This routine is used differently from Fortran
37 .vb
38         DM          da
39         integer     n,da_array(1)
40         PetscOffset i_da
41         integer     ierr
42         call DMDAGetGlobalIndices(da,n,da_array,i_da,ierr)
43 
44    C Access first local entry in list
45         value = da_array(i_da + 1)
46 .ve
47 
48    See Users-Manual: ch_fortran for details.
49 
50 .keywords: distributed array, get, global, indices, local-to-global
51 
52 .seealso: DMDACreate2d(), DMDAGetGhostCorners(), DMDAGetCorners(), DMLocalToGlobalBegin(), DMDARestoreGlobalIndices()
53           DMGlobalToLocalBegin(), DMGlobalToLocalEnd(), DMLocalToLocalBegin(), DMDAGetAO(), DMDAGetGlobalIndicesF90()
54           DMDAGetISLocalToGlobalMapping(), DMDACreate3d(), DMDACreate1d(), DMLocalToLocalEnd(), DMDAGetOwnershipRanges()
55 @*/
56 PetscErrorCode  DMDAGetGlobalIndices(DM da,PetscInt *n,const PetscInt *idx[])
57 {
58   PetscErrorCode ierr;
59 
60   PetscFunctionBegin;
61   PetscValidHeaderSpecific(da,DM_CLASSID,1);
62   if (n) {
63     PetscInt bs;
64     ierr = ISLocalToGlobalMappingGetSize(da->ltogmap,n);CHKERRQ(ierr);
65     ierr = ISLocalToGlobalMappingGetBlockSize(da->ltogmap,&bs);CHKERRQ(ierr);
66     *n = (*n)*bs;
67   }
68   if (idx) {
69     ierr = ISLocalToGlobalMappingGetIndices(da->ltogmap,idx);CHKERRQ(ierr);
70   }
71   PetscFunctionReturn(0);
72 }
73 
74 #undef __FUNCT__
75 #define __FUNCT__ "DMDAGetNatural_Private"
76 /*
77    Gets the natural number for each global number on the process.
78 
79    Used by DMDAGetAO() and DMDAGlobalToNatural_Create()
80 */
81 PetscErrorCode DMDAGetNatural_Private(DM da,PetscInt *outNlocal,IS *isnatural)
82 {
83   PetscErrorCode ierr;
84   PetscInt       Nlocal,i,j,k,*lidx,lict = 0;
85   DM_DA          *dd = (DM_DA*)da->data;
86 
87   PetscFunctionBegin;
88   Nlocal = (dd->xe-dd->xs);
89   if (dd->dim > 1) Nlocal *= (dd->ye-dd->ys);
90   if (dd->dim > 2) Nlocal *= (dd->ze-dd->zs);
91 
92   ierr = PetscMalloc1(Nlocal,&lidx);CHKERRQ(ierr);
93 
94   if (dd->dim == 1) {
95     for (i=dd->xs; i<dd->xe; i++) {
96       /*  global number in natural ordering */
97       lidx[lict++] = i;
98     }
99   } else if (dd->dim == 2) {
100     for (j=dd->ys; j<dd->ye; j++) {
101       for (i=dd->xs; i<dd->xe; i++) {
102         /*  global number in natural ordering */
103         lidx[lict++] = i + j*dd->M*dd->w;
104       }
105     }
106   } else if (dd->dim == 3) {
107     for (k=dd->zs; k<dd->ze; k++) {
108       for (j=dd->ys; j<dd->ye; j++) {
109         for (i=dd->xs; i<dd->xe; i++) {
110           lidx[lict++] = i + j*dd->M*dd->w + k*dd->M*dd->N*dd->w;
111         }
112       }
113     }
114   }
115   *outNlocal = Nlocal;
116   ierr       = ISCreateGeneral(PetscObjectComm((PetscObject)da),Nlocal,lidx,PETSC_OWN_POINTER,isnatural);CHKERRQ(ierr);
117   PetscFunctionReturn(0);
118 }
119 
120 #undef __FUNCT__
121 #define __FUNCT__ "DMDARestoreGlobalIndices"
122 /*@C
123    DMDARestoreGlobalIndices - Restores the global node number of all local nodes,   including ghost nodes.
124 
125    Not Collective
126 
127    Input Parameter:
128 .  da - the distributed array
129 
130    Output Parameters:
131 +  n - the number of local elements, including ghost nodes (or NULL)
132 -  idx - the global indices
133 
134    Level: intermediate
135 
136    Note:
137    For DMDA_STENCIL_STAR stencils the inactive corner ghost nodes are also included
138    in the list of local indices (even though those nodes are not updated
139    during calls to DMDAXXXToXXX().
140 
141    Essentially the same data is returned in the form of a local-to-global mapping
142    with the routine DMDAGetISLocalToGlobalMapping();
143 
144    Fortran Note:
145    This routine is used differently from Fortran
146 .vb
147         DM          da
148         integer     n,da_array(1)
149         PetscOffset i_da
150         integer     ierr
151         call DMDAGetGlobalIndices(da,n,da_array,i_da,ierr)
152 
153    C Access first local entry in list
154         value = da_array(i_da + 1)
155 .ve
156 
157    See Users-Manual: ch_fortran for details.
158 
159 .keywords: distributed array, get, global, indices, local-to-global
160 
161 .seealso: DMDACreate2d(), DMDAGetGhostCorners(), DMDAGetCorners(), DMLocalToGlobalBegin()
162           DMGlobalToLocalBegin(), DMGlobalToLocalEnd(), DMLocalToLocalBegin(), DMDAGetAO(), DMDAGetGlobalIndicesF90()
163           DMDAGetISLocalToGlobalMapping(), DMDACreate3d(), DMDACreate1d(), DMLocalToLocalEnd(), DMDAGetOwnershipRanges()
164 @*/
165 PetscErrorCode  DMDARestoreGlobalIndices(DM da,PetscInt *n,const PetscInt *idx[])
166 {
167   PetscErrorCode ierr;
168 
169   PetscFunctionBegin;
170   PetscValidHeaderSpecific(da,DM_CLASSID,1);
171   if (idx) {
172     ierr = ISLocalToGlobalMappingRestoreIndices(da->ltogmap,idx);CHKERRQ(ierr);
173   }
174   PetscFunctionReturn(0);
175 }
176 
177 #undef __FUNCT__
178 #define __FUNCT__ "DMDAGetAO"
179 /*@
180    DMDAGetAO - Gets the application ordering context for a distributed array.
181 
182    Collective on DMDA
183 
184    Input Parameter:
185 .  da - the distributed array
186 
187    Output Parameters:
188 .  ao - the application ordering context for DMDAs
189 
190    Level: intermediate
191 
192    Notes:
193    In this case, the AO maps to the natural grid ordering that would be used
194    for the DMDA if only 1 processor were employed (ordering most rapidly in the
195    x-direction, then y, then z).  Multiple degrees of freedom are numbered
196    for each node (rather than 1 component for the whole grid, then the next
197    component, etc.)
198 
199 .keywords: distributed array, get, global, indices, local-to-global
200 
201 .seealso: DMDACreate2d(), DMDAGetGhostCorners(), DMDAGetCorners(), DMDALocalToGlocal()
202           DMGlobalToLocalBegin(), DMGlobalToLocalEnd(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDAGetGlobalIndices(), DMDAGetOwnershipRanges(),
203           AO, AOPetscToApplication(), AOApplicationToPetsc()
204 @*/
205 PetscErrorCode  DMDAGetAO(DM da,AO *ao)
206 {
207   DM_DA *dd = (DM_DA*)da->data;
208 
209   PetscFunctionBegin;
210   PetscValidHeaderSpecific(da,DM_CLASSID,1);
211   PetscValidPointer(ao,2);
212 
213   /*
214      Build the natural ordering to PETSc ordering mappings.
215   */
216   if (!dd->ao) {
217     IS             ispetsc,isnatural;
218     PetscErrorCode ierr;
219     PetscInt       Nlocal;
220 
221     ierr = DMDAGetNatural_Private(da,&Nlocal,&isnatural);CHKERRQ(ierr);
222     ierr = ISCreateStride(PetscObjectComm((PetscObject)da),Nlocal,dd->base,1,&ispetsc);CHKERRQ(ierr);
223     ierr = AOCreateBasicIS(isnatural,ispetsc,&dd->ao);CHKERRQ(ierr);
224     ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)dd->ao);CHKERRQ(ierr);
225     ierr = ISDestroy(&ispetsc);CHKERRQ(ierr);
226     ierr = ISDestroy(&isnatural);CHKERRQ(ierr);
227   }
228   *ao = dd->ao;
229   PetscFunctionReturn(0);
230 }
231 
232 /*MC
233     DMDAGetGlobalIndicesF90 - Returns a Fortran90 pointer to the list of
234     global indices (global node number of all local nodes, including
235     ghost nodes).
236 
237     Synopsis:
238     DMDAGetGlobalIndicesF90(DM da,integer n,{integer, pointer :: idx(:)},integer ierr)
239 
240     Not Collective
241 
242     Input Parameter:
243 .   da - the distributed array
244 
245     Output Parameters:
246 +   n - the number of local elements, including ghost nodes (or NULL)
247 .   idx - the Fortran90 pointer to the global indices
248 -   ierr - error code
249 
250     Level: intermediate
251 
252 .keywords: distributed array, get, global, indices, local-to-global, f90
253 
254 .seealso: DMDAGetGlobalIndices(), DMDARestoreGlobalIndicesF90(), DMDARestoreGlobalIndices()
255 M*/
256 
257 /*MC
258     DMDARestoreGlobalIndicesF90 - Returns a Fortran90 pointer to the list of
259     global indices (global node number of all local nodes, including
260     ghost nodes).
261 
262     Synopsis:
263     DMDARestoreGlobalIndicesF90(DM da,integer n,{integer, pointer :: idx(:)},integer ierr)
264 
265     Not Collective
266 
267     Input Parameter:
268 .   da - the distributed array
269 
270     Output Parameters:
271 +   n - the number of local elements, including ghost nodes (or NULL)
272 .   idx - the Fortran90 pointer to the global indices
273 -   ierr - error code
274 
275     Level: intermediate
276 
277 .keywords: distributed array, get, global, indices, local-to-global, f90
278 
279 .seealso: DMDARestoreGlobalIndices(), DMDAGetGlobalIndicesF90(), DMDAGetGlobalIndices()
280 M*/
281 
282