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