xref: /petsc/src/vec/is/ao/impls/mapping/aomapping.c (revision 5b6bfdb9644f185dbf5e5a09b808ec241507e1e7)
1 
2 /*
3   These AO application ordering routines do not require that the input
4   be a permutation, but merely a 1-1 mapping. This implementation still
5   keeps the entire ordering on each processor.
6 */
7 
8 #include <../src/vec/is/ao/aoimpl.h>          /*I  "petscao.h" I*/
9 
10 typedef struct {
11   PetscInt N;
12   PetscInt *app;       /* app[i] is the partner for petsc[appPerm[i]] */
13   PetscInt *appPerm;
14   PetscInt *petsc;     /* petsc[j] is the partner for app[petscPerm[j]] */
15   PetscInt *petscPerm;
16 } AO_Mapping;
17 
18 PetscErrorCode AODestroy_Mapping(AO ao)
19 {
20   AO_Mapping     *aomap = (AO_Mapping*) ao->data;
21   PetscErrorCode ierr;
22 
23   PetscFunctionBegin;
24   ierr = PetscFree4(aomap->app,aomap->appPerm,aomap->petsc,aomap->petscPerm);CHKERRQ(ierr);
25   ierr = PetscFree(aomap);CHKERRQ(ierr);
26   PetscFunctionReturn(0);
27 }
28 
29 PetscErrorCode AOView_Mapping(AO ao, PetscViewer viewer)
30 {
31   AO_Mapping     *aomap = (AO_Mapping*) ao->data;
32   PetscMPIInt    rank;
33   PetscInt       i;
34   PetscBool      iascii;
35   PetscErrorCode ierr;
36 
37   PetscFunctionBegin;
38   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)ao), &rank);CHKERRQ(ierr);
39   if (rank) PetscFunctionReturn(0);
40   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
41   if (iascii) {
42     PetscViewerASCIIPrintf(viewer, "Number of elements in ordering %D\n", aomap->N);
43     PetscViewerASCIIPrintf(viewer, "   App.   PETSc\n");
44     for (i = 0; i < aomap->N; i++) {
45       PetscViewerASCIIPrintf(viewer, "%D   %D    %D\n", i, aomap->app[i], aomap->petsc[aomap->appPerm[i]]);
46     }
47   }
48   PetscFunctionReturn(0);
49 }
50 
51 PetscErrorCode AOPetscToApplication_Mapping(AO ao, PetscInt n, PetscInt *ia)
52 {
53   AO_Mapping *aomap = (AO_Mapping*) ao->data;
54   PetscInt   *app   = aomap->app;
55   PetscInt   *petsc = aomap->petsc;
56   PetscInt   *perm  = aomap->petscPerm;
57   PetscInt   N      = aomap->N;
58   PetscInt   low, high, mid=0;
59   PetscInt   idex;
60   PetscInt   i;
61 
62   /* It would be possible to use a single bisection search, which
63      recursively divided the indices to be converted, and searched
64      partitions which contained an index. This would result in
65      better running times if indices are clustered.
66   */
67   PetscFunctionBegin;
68   for (i = 0; i < n; i++) {
69     idex = ia[i];
70     if (idex < 0) continue;
71     /* Use bisection since the array is sorted */
72     low  = 0;
73     high = N - 1;
74     while (low <= high) {
75       mid = (low + high)/2;
76       if (idex == petsc[mid]) break;
77       else if (idex < petsc[mid]) high = mid - 1;
78       else low = mid + 1;
79     }
80     if (low > high) ia[i] = -1;
81     else            ia[i] = app[perm[mid]];
82   }
83   PetscFunctionReturn(0);
84 }
85 
86 PetscErrorCode AOApplicationToPetsc_Mapping(AO ao, PetscInt n, PetscInt *ia)
87 {
88   AO_Mapping *aomap = (AO_Mapping*) ao->data;
89   PetscInt   *app   = aomap->app;
90   PetscInt   *petsc = aomap->petsc;
91   PetscInt   *perm  = aomap->appPerm;
92   PetscInt   N      = aomap->N;
93   PetscInt   low, high, mid=0;
94   PetscInt   idex;
95   PetscInt   i;
96 
97   /* It would be possible to use a single bisection search, which
98      recursively divided the indices to be converted, and searched
99      partitions which contained an index. This would result in
100      better running times if indices are clustered.
101   */
102   PetscFunctionBegin;
103   for (i = 0; i < n; i++) {
104     idex = ia[i];
105     if (idex < 0) continue;
106     /* Use bisection since the array is sorted */
107     low  = 0;
108     high = N - 1;
109     while (low <= high) {
110       mid = (low + high)/2;
111       if (idex == app[mid]) break;
112       else if (idex < app[mid]) high = mid - 1;
113       else low = mid + 1;
114     }
115     if (low > high) ia[i] = -1;
116     else            ia[i] = petsc[perm[mid]];
117   }
118   PetscFunctionReturn(0);
119 }
120 
121 static struct _AOOps AOps = {AOView_Mapping,
122                              AODestroy_Mapping,
123                              AOPetscToApplication_Mapping,
124                              AOApplicationToPetsc_Mapping,
125                              NULL,
126                              NULL,
127                              NULL,
128                              NULL};
129 
130 /*@C
131   AOMappingHasApplicationIndex - Searches for the supplied application index.
132 
133   Input Parameters:
134 + ao       - The AOMapping
135 - index    - The application index
136 
137   Output Parameter:
138 . hasIndex - Flag is PETSC_TRUE if the index exists
139 
140   Level: intermediate
141 
142 .keywords: AO, index
143 .seealso: AOMappingHasPetscIndex(), AOCreateMapping()
144 @*/
145 PetscErrorCode  AOMappingHasApplicationIndex(AO ao, PetscInt idex, PetscBool  *hasIndex)
146 {
147   AO_Mapping *aomap;
148   PetscInt   *app;
149   PetscInt   low, high, mid;
150 
151   PetscFunctionBegin;
152   PetscValidHeaderSpecific(ao, AO_CLASSID,1);
153   PetscValidPointer(hasIndex,3);
154   aomap = (AO_Mapping*) ao->data;
155   app   = aomap->app;
156   /* Use bisection since the array is sorted */
157   low  = 0;
158   high = aomap->N - 1;
159   while (low <= high) {
160     mid = (low + high)/2;
161     if (idex == app[mid]) break;
162     else if (idex < app[mid]) high = mid - 1;
163     else low = mid + 1;
164   }
165   if (low > high) *hasIndex = PETSC_FALSE;
166   else *hasIndex = PETSC_TRUE;
167   PetscFunctionReturn(0);
168 }
169 
170 /*@
171   AOMappingHasPetscIndex - Searches for the supplied petsc index.
172 
173   Input Parameters:
174 + ao       - The AOMapping
175 - index    - The petsc index
176 
177   Output Parameter:
178 . hasIndex - Flag is PETSC_TRUE if the index exists
179 
180   Level: intermediate
181 
182 .keywords: AO, index
183 .seealso: AOMappingHasApplicationIndex(), AOCreateMapping()
184 @*/
185 PetscErrorCode  AOMappingHasPetscIndex(AO ao, PetscInt idex, PetscBool  *hasIndex)
186 {
187   AO_Mapping *aomap;
188   PetscInt   *petsc;
189   PetscInt   low, high, mid;
190 
191   PetscFunctionBegin;
192   PetscValidHeaderSpecific(ao, AO_CLASSID,1);
193   PetscValidPointer(hasIndex,3);
194   aomap = (AO_Mapping*) ao->data;
195   petsc = aomap->petsc;
196   /* Use bisection since the array is sorted */
197   low  = 0;
198   high = aomap->N - 1;
199   while (low <= high) {
200     mid = (low + high)/2;
201     if (idex == petsc[mid]) break;
202     else if (idex < petsc[mid]) high = mid - 1;
203     else low = mid + 1;
204   }
205   if (low > high) *hasIndex = PETSC_FALSE;
206   else *hasIndex = PETSC_TRUE;
207   PetscFunctionReturn(0);
208 }
209 
210 /*@C
211   AOCreateMapping - Creates a basic application mapping using two integer arrays.
212 
213   Input Parameters:
214 + comm    - MPI communicator that is to share AO
215 . napp    - size of integer arrays
216 . myapp   - integer array that defines an ordering
217 - mypetsc - integer array that defines another ordering (may be NULL to indicate the identity ordering)
218 
219   Output Parameter:
220 . aoout   - the new application mapping
221 
222   Options Database Key:
223 . -ao_view : call AOView() at the conclusion of AOCreateMapping()
224 
225   Level: beginner
226 
227     Notes: the arrays myapp and mypetsc need NOT contain the all the integers 0 to napp-1, that is there CAN be "holes"  in the indices.
228        Use AOCreateBasic() or AOCreateBasicIS() if they do not have holes for better performance.
229 
230 .keywords: AO, create
231 .seealso: AOCreateBasic(), AOCreateBasic(), AOCreateMappingIS(), AODestroy()
232 @*/
233 PetscErrorCode  AOCreateMapping(MPI_Comm comm,PetscInt napp,const PetscInt myapp[],const PetscInt mypetsc[],AO *aoout)
234 {
235   AO             ao;
236   AO_Mapping     *aomap;
237   PetscInt       *allpetsc,  *allapp;
238   PetscInt       *petscPerm, *appPerm;
239   PetscInt       *petsc;
240   PetscMPIInt    size, rank,*lens, *disp,nnapp;
241   PetscInt       N, start;
242   PetscInt       i;
243   PetscErrorCode ierr;
244 
245   PetscFunctionBegin;
246   PetscValidPointer(aoout,5);
247   *aoout = 0;
248   ierr = AOInitializePackage();CHKERRQ(ierr);
249 
250   ierr     = PetscHeaderCreate(ao, AO_CLASSID, "AO", "Application Ordering", "AO", comm, AODestroy, AOView);CHKERRQ(ierr);
251   ierr     = PetscNewLog(ao,&aomap);CHKERRQ(ierr);
252   ierr     = PetscMemcpy(ao->ops, &AOps, sizeof(AOps));CHKERRQ(ierr);
253   ao->data = (void*) aomap;
254 
255   /* transmit all lengths to all processors */
256   ierr  = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
257   ierr  = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
258   ierr  = PetscMalloc2(size, &lens,size,&disp);CHKERRQ(ierr);
259   nnapp = napp;
260   ierr  = MPI_Allgather(&nnapp, 1, MPI_INT, lens, 1, MPI_INT, comm);CHKERRQ(ierr);
261   N     = 0;
262   for (i = 0; i < size; i++) {
263     disp[i] = N;
264     N      += lens[i];
265   }
266   aomap->N = N;
267   ao->N    = N;
268   ao->n    = N;
269 
270   /* If mypetsc is 0 then use "natural" numbering */
271   if (!mypetsc) {
272     start = disp[rank];
273     ierr  = PetscMalloc1(napp+1, &petsc);CHKERRQ(ierr);
274     for (i = 0; i < napp; i++) petsc[i] = start + i;
275   } else {
276     petsc = (PetscInt*)mypetsc;
277   }
278 
279   /* get all indices on all processors */
280   ierr = PetscMalloc4(N, &allapp,N,&appPerm,N,&allpetsc,N,&petscPerm);CHKERRQ(ierr);
281   ierr = MPI_Allgatherv((void*)myapp, napp, MPIU_INT, allapp,   lens, disp, MPIU_INT, comm);CHKERRQ(ierr);
282   ierr = MPI_Allgatherv((void*)petsc, napp, MPIU_INT, allpetsc, lens, disp, MPIU_INT, comm);CHKERRQ(ierr);
283   ierr = PetscFree2(lens,disp);CHKERRQ(ierr);
284 
285   /* generate a list of application and PETSc node numbers */
286   ierr = PetscMalloc4(N, &aomap->app,N,&aomap->appPerm,N,&aomap->petsc,N,&aomap->petscPerm);CHKERRQ(ierr);
287   ierr = PetscLogObjectMemory((PetscObject)ao, 4*N * sizeof(PetscInt));CHKERRQ(ierr);
288   for (i = 0; i < N; i++) {
289     appPerm[i]   = i;
290     petscPerm[i] = i;
291   }
292   ierr = PetscSortIntWithPermutation(N, allpetsc, petscPerm);CHKERRQ(ierr);
293   ierr = PetscSortIntWithPermutation(N, allapp,   appPerm);CHKERRQ(ierr);
294   /* Form sorted arrays of indices */
295   for (i = 0; i < N; i++) {
296     aomap->app[i]   = allapp[appPerm[i]];
297     aomap->petsc[i] = allpetsc[petscPerm[i]];
298   }
299   /* Invert petscPerm[] into aomap->petscPerm[] */
300   for (i = 0; i < N; i++) aomap->petscPerm[petscPerm[i]] = i;
301 
302   /* Form map between aomap->app[] and aomap->petsc[] */
303   for (i = 0; i < N; i++) aomap->appPerm[i] = aomap->petscPerm[appPerm[i]];
304 
305   /* Invert appPerm[] into allapp[] */
306   for (i = 0; i < N; i++) allapp[appPerm[i]] = i;
307 
308   /* Form map between aomap->petsc[] and aomap->app[] */
309   for (i = 0; i < N; i++) aomap->petscPerm[i] = allapp[petscPerm[i]];
310 
311 #if defined(PETSC_USE_DEBUG)
312   /* Check that the permutations are complementary */
313   for (i = 0; i < N; i++) {
314     if (i != aomap->appPerm[aomap->petscPerm[i]]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Invalid ordering");
315   }
316 #endif
317   /* Cleanup */
318   if (!mypetsc) {
319     ierr = PetscFree(petsc);CHKERRQ(ierr);
320   }
321   ierr = PetscFree4(allapp,appPerm,allpetsc,petscPerm);CHKERRQ(ierr);
322 
323   ierr = AOViewFromOptions(ao,NULL,"-ao_view");CHKERRQ(ierr);
324 
325   *aoout = ao;
326   PetscFunctionReturn(0);
327 }
328 
329 /*@C
330   AOCreateMappingIS - Creates a basic application ordering using two index sets.
331 
332   Input Parameters:
333 + comm    - MPI communicator that is to share AO
334 . isapp   - index set that defines an ordering
335 - ispetsc - index set that defines another ordering, maybe NULL for identity IS
336 
337   Output Parameter:
338 . aoout   - the new application ordering
339 
340   Options Database Key:
341 . -ao_view : call AOView() at the conclusion of AOCreateMappingIS()
342 
343   Level: beginner
344 
345     Notes: the index sets isapp and ispetsc need NOT contain the all the integers 0 to N-1, that is there CAN be "holes"  in the indices.
346        Use AOCreateBasic() or AOCreateBasicIS() if they do not have holes for better performance.
347 
348 .keywords: AO, create
349 .seealso: AOCreateBasic(), AOCreateMapping(), AODestroy()
350 @*/
351 PetscErrorCode  AOCreateMappingIS(IS isapp, IS ispetsc, AO *aoout)
352 {
353   MPI_Comm       comm;
354   const PetscInt *mypetsc, *myapp;
355   PetscInt       napp, npetsc;
356   PetscErrorCode ierr;
357 
358   PetscFunctionBegin;
359   ierr = PetscObjectGetComm((PetscObject) isapp, &comm);CHKERRQ(ierr);
360   ierr = ISGetLocalSize(isapp, &napp);CHKERRQ(ierr);
361   if (ispetsc) {
362     ierr = ISGetLocalSize(ispetsc, &npetsc);CHKERRQ(ierr);
363     if (napp != npetsc) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "Local IS lengths must match");
364     ierr = ISGetIndices(ispetsc, &mypetsc);CHKERRQ(ierr);
365   } else {
366     mypetsc = NULL;
367   }
368   ierr = ISGetIndices(isapp, &myapp);CHKERRQ(ierr);
369 
370   ierr = AOCreateMapping(comm, napp, myapp, mypetsc, aoout);CHKERRQ(ierr);
371 
372   ierr = ISRestoreIndices(isapp, &myapp);CHKERRQ(ierr);
373   if (ispetsc) {
374     ierr = ISRestoreIndices(ispetsc, &mypetsc);CHKERRQ(ierr);
375   }
376   PetscFunctionReturn(0);
377 }
378