xref: /petsc/src/vec/is/ao/impls/mapping/aomapping.c (revision fa2bb9fe3cba50bac3b093e83d6f95e65b517df9)
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:
228     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.
229        Use AOCreateBasic() or AOCreateBasicIS() if they do not have holes for better performance.
230 
231 .keywords: AO, create
232 .seealso: AOCreateBasic(), AOCreateBasic(), AOCreateMappingIS(), AODestroy()
233 @*/
234 PetscErrorCode  AOCreateMapping(MPI_Comm comm,PetscInt napp,const PetscInt myapp[],const PetscInt mypetsc[],AO *aoout)
235 {
236   AO             ao;
237   AO_Mapping     *aomap;
238   PetscInt       *allpetsc,  *allapp;
239   PetscInt       *petscPerm, *appPerm;
240   PetscInt       *petsc;
241   PetscMPIInt    size, rank,*lens, *disp,nnapp;
242   PetscInt       N, start;
243   PetscInt       i;
244   PetscErrorCode ierr;
245 
246   PetscFunctionBegin;
247   PetscValidPointer(aoout,5);
248   *aoout = 0;
249   ierr = AOInitializePackage();CHKERRQ(ierr);
250 
251   ierr     = PetscHeaderCreate(ao, AO_CLASSID, "AO", "Application Ordering", "AO", comm, AODestroy, AOView);CHKERRQ(ierr);
252   ierr     = PetscNewLog(ao,&aomap);CHKERRQ(ierr);
253   ierr     = PetscMemcpy(ao->ops, &AOps, sizeof(AOps));CHKERRQ(ierr);
254   ao->data = (void*) aomap;
255 
256   /* transmit all lengths to all processors */
257   ierr  = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
258   ierr  = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
259   ierr  = PetscMalloc2(size, &lens,size,&disp);CHKERRQ(ierr);
260   nnapp = napp;
261   ierr  = MPI_Allgather(&nnapp, 1, MPI_INT, lens, 1, MPI_INT, comm);CHKERRQ(ierr);
262   N     = 0;
263   for (i = 0; i < size; i++) {
264     disp[i] = N;
265     N      += lens[i];
266   }
267   aomap->N = N;
268   ao->N    = N;
269   ao->n    = N;
270 
271   /* If mypetsc is 0 then use "natural" numbering */
272   if (!mypetsc) {
273     start = disp[rank];
274     ierr  = PetscMalloc1(napp+1, &petsc);CHKERRQ(ierr);
275     for (i = 0; i < napp; i++) petsc[i] = start + i;
276   } else {
277     petsc = (PetscInt*)mypetsc;
278   }
279 
280   /* get all indices on all processors */
281   ierr = PetscMalloc4(N, &allapp,N,&appPerm,N,&allpetsc,N,&petscPerm);CHKERRQ(ierr);
282   ierr = MPI_Allgatherv((void*)myapp, napp, MPIU_INT, allapp,   lens, disp, MPIU_INT, comm);CHKERRQ(ierr);
283   ierr = MPI_Allgatherv((void*)petsc, napp, MPIU_INT, allpetsc, lens, disp, MPIU_INT, comm);CHKERRQ(ierr);
284   ierr = PetscFree2(lens,disp);CHKERRQ(ierr);
285 
286   /* generate a list of application and PETSc node numbers */
287   ierr = PetscMalloc4(N, &aomap->app,N,&aomap->appPerm,N,&aomap->petsc,N,&aomap->petscPerm);CHKERRQ(ierr);
288   ierr = PetscLogObjectMemory((PetscObject)ao, 4*N * sizeof(PetscInt));CHKERRQ(ierr);
289   for (i = 0; i < N; i++) {
290     appPerm[i]   = i;
291     petscPerm[i] = i;
292   }
293   ierr = PetscSortIntWithPermutation(N, allpetsc, petscPerm);CHKERRQ(ierr);
294   ierr = PetscSortIntWithPermutation(N, allapp,   appPerm);CHKERRQ(ierr);
295   /* Form sorted arrays of indices */
296   for (i = 0; i < N; i++) {
297     aomap->app[i]   = allapp[appPerm[i]];
298     aomap->petsc[i] = allpetsc[petscPerm[i]];
299   }
300   /* Invert petscPerm[] into aomap->petscPerm[] */
301   for (i = 0; i < N; i++) aomap->petscPerm[petscPerm[i]] = i;
302 
303   /* Form map between aomap->app[] and aomap->petsc[] */
304   for (i = 0; i < N; i++) aomap->appPerm[i] = aomap->petscPerm[appPerm[i]];
305 
306   /* Invert appPerm[] into allapp[] */
307   for (i = 0; i < N; i++) allapp[appPerm[i]] = i;
308 
309   /* Form map between aomap->petsc[] and aomap->app[] */
310   for (i = 0; i < N; i++) aomap->petscPerm[i] = allapp[petscPerm[i]];
311 
312 #if defined(PETSC_USE_DEBUG)
313   /* Check that the permutations are complementary */
314   for (i = 0; i < N; i++) {
315     if (i != aomap->appPerm[aomap->petscPerm[i]]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Invalid ordering");
316   }
317 #endif
318   /* Cleanup */
319   if (!mypetsc) {
320     ierr = PetscFree(petsc);CHKERRQ(ierr);
321   }
322   ierr = PetscFree4(allapp,appPerm,allpetsc,petscPerm);CHKERRQ(ierr);
323 
324   ierr = AOViewFromOptions(ao,NULL,"-ao_view");CHKERRQ(ierr);
325 
326   *aoout = ao;
327   PetscFunctionReturn(0);
328 }
329 
330 /*@C
331   AOCreateMappingIS - Creates a basic application ordering using two index sets.
332 
333   Input Parameters:
334 + comm    - MPI communicator that is to share AO
335 . isapp   - index set that defines an ordering
336 - ispetsc - index set that defines another ordering, maybe NULL for identity IS
337 
338   Output Parameter:
339 . aoout   - the new application ordering
340 
341   Options Database Key:
342 . -ao_view : call AOView() at the conclusion of AOCreateMappingIS()
343 
344   Level: beginner
345 
346     Notes:
347     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.
348        Use AOCreateBasic() or AOCreateBasicIS() if they do not have holes for better performance.
349 
350 .keywords: AO, create
351 .seealso: AOCreateBasic(), AOCreateMapping(), AODestroy()
352 @*/
353 PetscErrorCode  AOCreateMappingIS(IS isapp, IS ispetsc, AO *aoout)
354 {
355   MPI_Comm       comm;
356   const PetscInt *mypetsc, *myapp;
357   PetscInt       napp, npetsc;
358   PetscErrorCode ierr;
359 
360   PetscFunctionBegin;
361   ierr = PetscObjectGetComm((PetscObject) isapp, &comm);CHKERRQ(ierr);
362   ierr = ISGetLocalSize(isapp, &napp);CHKERRQ(ierr);
363   if (ispetsc) {
364     ierr = ISGetLocalSize(ispetsc, &npetsc);CHKERRQ(ierr);
365     if (napp != npetsc) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "Local IS lengths must match");
366     ierr = ISGetIndices(ispetsc, &mypetsc);CHKERRQ(ierr);
367   } else {
368     mypetsc = NULL;
369   }
370   ierr = ISGetIndices(isapp, &myapp);CHKERRQ(ierr);
371 
372   ierr = AOCreateMapping(comm, napp, myapp, mypetsc, aoout);CHKERRQ(ierr);
373 
374   ierr = ISRestoreIndices(isapp, &myapp);CHKERRQ(ierr);
375   if (ispetsc) {
376     ierr = ISRestoreIndices(ispetsc, &mypetsc);CHKERRQ(ierr);
377   }
378   PetscFunctionReturn(0);
379 }
380