xref: /petsc/src/dm/impls/da/daindex.c (revision 2205254efee3a00a594e5e2a3a70f74dcb40bc03)
1 
2 /*
3   Code for manipulating distributed regular arrays in parallel.
4 */
5 
6 #include <petsc-private/daimpl.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 PETSC_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();
32 
33    Fortran Note:
34    This routine is used differently from Fortran
35 .vb
36         DM          da
37         integer     n,da_array(1)
38         PetscOffset i_da
39         integer     ierr
40         call DMDAGetGlobalIndices(da,n,da_array,i_da,ierr)
41 
42    C Access first local entry in list
43         value = da_array(i_da + 1)
44 .ve
45 
46    See the <A href="../../docs/manual.pdf#nameddest=ch_fortran">Fortran chapter</A> of the users manual for details.
47 
48 .keywords: distributed array, get, global, indices, local-to-global
49 
50 .seealso: DMDACreate2d(), DMDAGetGhostCorners(), DMDAGetCorners(), DMLocalToGlobalBegin()
51           DMGlobalToLocalBegin(), DMGlobalToLocalEnd(), DMDALocalToLocalBegin(), DMDAGetAO(), DMDAGetGlobalIndicesF90()
52           DMDAGetISLocalToGlobalMapping(), DMDACreate3d(), DMDACreate1d(), DMDALocalToLocalEnd(), DMDAGetOwnershipRanges()
53 @*/
54 PetscErrorCode  DMDAGetGlobalIndices(DM da,PetscInt *n,PetscInt **idx)
55 {
56   DM_DA *dd = (DM_DA*)da->data;
57 
58   PetscFunctionBegin;
59   PetscValidHeaderSpecific(da,DM_CLASSID,1);
60   if (n) *n = dd->Nl;
61   if (idx) *idx = dd->idx;
62   PetscFunctionReturn(0);
63 }
64 
65 #undef __FUNCT__
66 #define __FUNCT__ "DMDAGetNatural_Private"
67 /*
68    Gets the natural number for each global number on the process.
69 
70    Used by DMDAGetAO() and DMDAGlobalToNatural_Create()
71 */
72 PetscErrorCode DMDAGetNatural_Private(DM da,PetscInt *outNlocal,IS *isnatural)
73 {
74   PetscErrorCode ierr;
75   PetscInt       Nlocal,i,j,k,*lidx,lict = 0;
76   DM_DA          *dd = (DM_DA*)da->data;
77 
78   PetscFunctionBegin;
79   Nlocal = (dd->xe-dd->xs);
80   if (dd->dim > 1) Nlocal *= (dd->ye-dd->ys);
81   if (dd->dim > 2) Nlocal *= (dd->ze-dd->zs);
82 
83   ierr = PetscMalloc(Nlocal*sizeof(PetscInt),&lidx);CHKERRQ(ierr);
84 
85   if (dd->dim == 1) {
86     for (i=dd->xs; i<dd->xe; i++) {
87       /*  global number in natural ordering */
88       lidx[lict++] = i;
89     }
90   } else if (dd->dim == 2) {
91     for (j=dd->ys; j<dd->ye; j++) {
92       for (i=dd->xs; i<dd->xe; i++) {
93         /*  global number in natural ordering */
94         lidx[lict++] = i + j*dd->M*dd->w;
95       }
96     }
97   } else if (dd->dim == 3) {
98     for (k=dd->zs; k<dd->ze; k++) {
99       for (j=dd->ys; j<dd->ye; j++) {
100         for (i=dd->xs; i<dd->xe; i++) {
101           lidx[lict++] = i + j*dd->M*dd->w + k*dd->M*dd->N*dd->w;
102         }
103       }
104     }
105   }
106   *outNlocal = Nlocal;
107   ierr       = ISCreateGeneral(((PetscObject)da)->comm,Nlocal,lidx,PETSC_OWN_POINTER,isnatural);CHKERRQ(ierr);
108   PetscFunctionReturn(0);
109 }
110 
111 #undef __FUNCT__
112 #define __FUNCT__ "DMDAGetAO"
113 /*@
114    DMDAGetAO - Gets the application ordering context for a distributed array.
115 
116    Collective on DMDA
117 
118    Input Parameter:
119 .  da - the distributed array
120 
121    Output Parameters:
122 .  ao - the application ordering context for DMDAs
123 
124    Level: intermediate
125 
126    Notes:
127    In this case, the AO maps to the natural grid ordering that would be used
128    for the DMDA if only 1 processor were employed (ordering most rapidly in the
129    x-direction, then y, then z).  Multiple degrees of freedom are numbered
130    for each node (rather than 1 component for the whole grid, then the next
131    component, etc.)
132 
133 .keywords: distributed array, get, global, indices, local-to-global
134 
135 .seealso: DMDACreate2d(), DMDAGetGhostCorners(), DMDAGetCorners(), DMDALocalToGlocal()
136           DMGlobalToLocalBegin(), DMGlobalToLocalEnd(), DMDALocalToLocalBegin(), DMDALocalToLocalEnd(), DMDAGetGlobalIndices(), DMDAGetOwnershipRanges(),
137           AO, AOPetscToApplication(), AOApplicationToPetsc()
138 @*/
139 PetscErrorCode  DMDAGetAO(DM da,AO *ao)
140 {
141   DM_DA *dd = (DM_DA*)da->data;
142 
143   PetscFunctionBegin;
144   PetscValidHeaderSpecific(da,DM_CLASSID,1);
145   PetscValidPointer(ao,2);
146 
147   /*
148      Build the natural ordering to PETSc ordering mappings.
149   */
150   if (!dd->ao) {
151     IS             ispetsc,isnatural;
152     PetscErrorCode ierr;
153     PetscInt       Nlocal;
154 
155     ierr = DMDAGetNatural_Private(da,&Nlocal,&isnatural);CHKERRQ(ierr);
156     ierr = ISCreateStride(((PetscObject)da)->comm,Nlocal,dd->base,1,&ispetsc);CHKERRQ(ierr);
157     ierr = AOCreateBasicIS(isnatural,ispetsc,&dd->ao);CHKERRQ(ierr);
158     ierr = PetscLogObjectParent(da,dd->ao);CHKERRQ(ierr);
159     ierr = ISDestroy(&ispetsc);CHKERRQ(ierr);
160     ierr = ISDestroy(&isnatural);CHKERRQ(ierr);
161   }
162   *ao = dd->ao;
163   PetscFunctionReturn(0);
164 }
165 
166 /*MC
167     DMDAGetGlobalIndicesF90 - Returns a Fortran90 pointer to the list of
168     global indices (global node number of all local nodes, including
169     ghost nodes).
170 
171     Synopsis:
172     DMDAGetGlobalIndicesF90(DM da,integer n,{integer, pointer :: idx(:)},integer ierr)
173 
174     Not Collective
175 
176     Input Parameter:
177 .   da - the distributed array
178 
179     Output Parameters:
180 +   n - the number of local elements, including ghost nodes (or PETSC_NULL)
181 .   idx - the Fortran90 pointer to the global indices
182 -   ierr - error code
183 
184     Level: intermediate
185 
186     Notes:
187      Not yet supported for all F90 compilers
188 
189 .keywords: distributed array, get, global, indices, local-to-global, f90
190 
191 .seealso: DMDAGetGlobalIndices()
192 M*/
193