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