Actual source code: aomapping.c
1: #ifdef PETSC_RCS_HEADER
2: static char vcid[] = "$Id: aomapping.c,v 1.3 2000/01/10 03:15:40 knepley Exp $";
3: #endif
5: /*
6: These AO application ordering routines do not require that the input
7: be a permutation, but merely a 1-1 mapping. This implementation still
8: keeps the entire ordering on each processor.
9: */
11: #include src/dm/ao/aoimpl.h
12: #include petscsys.h
14: typedef struct {
15: int N;
16: int *app; /* app[i] is the partner for petsc[appPerm[i]] */
17: int *appPerm;
18: int *petsc; /* petsc[j] is the partner for app[petscPerm[j]] */
19: int *petscPerm;
20: } AO_Mapping;
24: int AODestroy_Mapping(AO ao)
25: {
26: AO_Mapping *aomap = (AO_Mapping *) ao->data;
27: int ierr;
30: PetscFree(aomap->app);
31: PetscFree(ao->data);
32: return(0);
33: }
37: int AOView_Mapping(AO ao, PetscViewer viewer)
38: {
39: AO_Mapping *aomap = (AO_Mapping *) ao->data;
40: int rank, i;
41: PetscTruth isascii;
42: int ierr;
45: MPI_Comm_rank(ao->comm, &rank);
46: if (rank) return(0);
48: if (!viewer) {
49: viewer = PETSC_VIEWER_STDOUT_SELF;
50: }
52: PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_ASCII, &isascii);
53: if (isascii == PETSC_TRUE) {
54: PetscViewerASCIIPrintf(viewer, "Number of elements in ordering %d\n", aomap->N);
55: PetscViewerASCIIPrintf(viewer, " App. PETSc\n");
56: for(i = 0; i < aomap->N; i++) {
57: PetscViewerASCIIPrintf(viewer, "%d %d %d\n", i, aomap->app[i], aomap->petsc[aomap->appPerm[i]]);
58: }
59: }
60: return(0);
61: }
65: int AOPetscToApplication_Mapping(AO ao, int n, int *ia)
66: {
67: AO_Mapping *aomap = (AO_Mapping *) ao->data;
68: int *app = aomap->app;
69: int *petsc = aomap->petsc;
70: int *perm = aomap->petscPerm;
71: int N = aomap->N;
72: int low, high, mid=0;
73: int idex;
74: int i;
76: /* It would be possible to use a single bisection search, which
77: recursively divided the indices to be converted, and searched
78: partitions which contained an index. This would result in
79: better running times if indices are clustered.
80: */
82: for(i = 0; i < n; i++) {
83: idex = ia[i];
84: if (idex < 0) continue;
85: /* Use bisection since the array is sorted */
86: low = 0;
87: high = N - 1;
88: while (low <= high) {
89: mid = (low + high)/2;
90: if (idex == petsc[mid]) {
91: break;
92: } else if (idex < petsc[mid]) {
93: high = mid - 1;
94: } else {
95: low = mid + 1;
96: }
97: }
98: if (low > high) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Invalid input index %d", idex);
99: ia[i] = app[perm[mid]];
100: }
101: return(0);
102: }
106: int AOApplicationToPetsc_Mapping(AO ao, int n, int *ia)
107: {
108: AO_Mapping *aomap = (AO_Mapping *) ao->data;
109: int *app = aomap->app;
110: int *petsc = aomap->petsc;
111: int *perm = aomap->appPerm;
112: int N = aomap->N;
113: int low, high, mid=0;
114: int idex;
115: int i;
117: /* It would be possible to use a single bisection search, which
118: recursively divided the indices to be converted, and searched
119: partitions which contained an index. This would result in
120: better running times if indices are clustered.
121: */
123: for(i = 0; i < n; i++) {
124: idex = ia[i];
125: if (idex < 0) continue;
126: /* Use bisection since the array is sorted */
127: low = 0;
128: high = N - 1;
129: while (low <= high) {
130: mid = (low + high)/2;
131: if (idex == app[mid]) {
132: break;
133: } else if (idex < app[mid]) {
134: high = mid - 1;
135: } else {
136: low = mid + 1;
137: }
138: }
139: if (low > high) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Invalid input index %d", idex);
140: ia[i] = petsc[perm[mid]];
141: }
142: return(0);
143: }
145: static struct _AOOps AOps = {AOView_Mapping,
146: AODestroy_Mapping,
147: AOPetscToApplication_Mapping,
148: AOApplicationToPetsc_Mapping,
149: PETSC_NULL,
150: PETSC_NULL,
151: PETSC_NULL,
152: PETSC_NULL};
156: /*@C
157: AOMappingHasApplicationIndex - Searches for the supplied application index.
159: Input Parameters:
160: + ao - The AOMapping
161: - index - The application index
163: Output Parameter:
164: . hasIndex - Flag is PETSC_TRUE if the index exists
166: Level: intermediate
168: .keywords: AO, index
169: .seealso: AOMappingHasPetscIndex(), AOCreateMapping()
170: @*/
171: int AOMappingHasApplicationIndex(AO ao, int idex, PetscTruth *hasIndex)
172: {
173: AO_Mapping *aomap;
174: int *app;
175: int low, high, mid;
180: aomap = (AO_Mapping *) ao->data;
181: app = aomap->app;
182: /* Use bisection since the array is sorted */
183: low = 0;
184: high = aomap->N - 1;
185: while (low <= high) {
186: mid = (low + high)/2;
187: if (idex == app[mid]) {
188: break;
189: } else if (idex < app[mid]) {
190: high = mid - 1;
191: } else {
192: low = mid + 1;
193: }
194: }
195: if (low > high) {
196: *hasIndex = PETSC_FALSE;
197: } else {
198: *hasIndex = PETSC_TRUE;
199: }
200: return(0);
201: }
205: /*@C
206: AOMappingHasPetscIndex - Searches for the supplied petsc index.
208: Input Parameters:
209: + ao - The AOMapping
210: - index - The petsc index
212: Output Parameter:
213: . hasIndex - Flag is PETSC_TRUE if the index exists
215: Level: intermediate
217: .keywords: AO, index
218: .seealso: AOMappingHasApplicationIndex(), AOCreateMapping()
219: @*/
220: int AOMappingHasPetscIndex(AO ao, int idex, PetscTruth *hasIndex)
221: {
222: AO_Mapping *aomap;
223: int *petsc;
224: int low, high, mid;
229: aomap = (AO_Mapping *) ao->data;
230: petsc = aomap->petsc;
231: /* Use bisection since the array is sorted */
232: low = 0;
233: high = aomap->N - 1;
234: while (low <= high) {
235: mid = (low + high)/2;
236: if (idex == petsc[mid]) {
237: break;
238: } else if (idex < petsc[mid]) {
239: high = mid - 1;
240: } else {
241: low = mid + 1;
242: }
243: }
244: if (low > high) {
245: *hasIndex = PETSC_FALSE;
246: } else {
247: *hasIndex = PETSC_TRUE;
248: }
249: return(0);
250: }
254: /*@C
255: AOCreateMapping - Creates a basic application mapping using two integer arrays.
257: Input Parameters:
258: + comm - MPI communicator that is to share AO
259: . napp - size of integer arrays
260: . myapp - integer array that defines an ordering
261: - mypetsc - integer array that defines another ordering
263: Output Parameter:
264: . aoout - the new application mapping
266: Options Database Key:
267: $ -ao_view : call AOView() at the conclusion of AOCreateMapping()
269: Level: beginner
271: .keywords: AO, create
272: .seealso: AOCreateDebug(), AOCreateBasic(), AOCreateMappingIS(), AODestroy()
273: @*/
274: int AOCreateMapping(MPI_Comm comm,int napp,const int myapp[],const int mypetsc[],AO *aoout)
275: {
276: AO ao;
277: AO_Mapping *aomap;
278: int *lens, *disp;
279: int *allpetsc, *allapp;
280: int *petscPerm, *appPerm;
281: int *petsc;
282: int size, rank;
283: int N, start;
284: int i;
285: PetscTruth opt;
286: int ierr;
290: *aoout = 0;
291: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
292: DMInitializePackage(PETSC_NULL);
293: #endif
295: PetscHeaderCreate(ao, _p_AO, struct _AOOps, AO_COOKIE, AO_MAPPING, "AO", comm, AODestroy, AOView);
296: PetscLogObjectCreate(ao);
297: PetscNew(AO_Mapping, &aomap);
298: PetscLogObjectMemory(ao, sizeof(struct _p_AO) + sizeof(AO_Mapping));
299: PetscMemcpy(ao->ops, &AOps, sizeof(AOps));
300: ao->data = (void *) aomap;
302: /* transmit all lengths to all processors */
303: MPI_Comm_size(comm, &size);
304: MPI_Comm_rank(comm, &rank);
305: PetscMalloc(2*size * sizeof(int), &lens);
306: disp = lens + size;
307: MPI_Allgather(&napp, 1, MPI_INT, lens, 1, MPI_INT, comm);
308: N = 0;
309: for(i = 0; i < size; i++) {
310: disp[i] = N;
311: N += lens[i];
312: }
313: aomap->N = N;
314: ao->N = N;
315: ao->n = N;
317: /* If mypetsc is 0 then use "natural" numbering */
318: if (!mypetsc) {
319: start = disp[rank];
320: PetscMalloc((napp+1) * sizeof(int), &petsc);
321: for(i = 0; i < napp; i++) {
322: petsc[i] = start + i;
323: }
324: } else {
325: petsc = (int *)mypetsc;
326: }
328: /* get all indices on all processors */
329: PetscMalloc(N*4 * sizeof(int), &allapp);
330: appPerm = allapp + N;
331: allpetsc = appPerm + N;
332: petscPerm = allpetsc + N;
333: MPI_Allgatherv((void*)myapp, napp, MPI_INT, allapp, lens, disp, MPI_INT, comm);
334: MPI_Allgatherv((void*)mypetsc, napp, MPI_INT, allpetsc, lens, disp, MPI_INT, comm);
335: PetscFree(lens);
337: /* generate a list of application and PETSc node numbers */
338: PetscMalloc(N*4 * sizeof(int), &aomap->app);
339: PetscLogObjectMemory(ao, 4*N * sizeof(int));
340: aomap->appPerm = aomap->app + N;
341: aomap->petsc = aomap->appPerm + N;
342: aomap->petscPerm = aomap->petsc + N;
343: for(i = 0; i < N; i++) {
344: appPerm[i] = i;
345: petscPerm[i] = i;
346: }
347: PetscSortIntWithPermutation(N, allpetsc, petscPerm);
348: PetscSortIntWithPermutation(N, allapp, appPerm);
349: /* Form sorted arrays of indices */
350: for(i = 0; i < N; i++) {
351: aomap->app[i] = allapp[appPerm[i]];
352: aomap->petsc[i] = allpetsc[petscPerm[i]];
353: }
354: /* Invert petscPerm[] into aomap->petscPerm[] */
355: for(i = 0; i < N; i++) {
356: aomap->petscPerm[petscPerm[i]] = i;
357: }
358: /* Form map between aomap->app[] and aomap->petsc[] */
359: for(i = 0; i < N; i++) {
360: aomap->appPerm[i] = aomap->petscPerm[appPerm[i]];
361: }
362: /* Invert appPerm[] into allapp[] */
363: for(i = 0; i < N; i++) {
364: allapp[appPerm[i]] = i;
365: }
366: /* Form map between aomap->petsc[] and aomap->app[] */
367: for(i = 0; i < N; i++) {
368: aomap->petscPerm[i] = allapp[petscPerm[i]];
369: }
370: #ifdef PETSC_USE_BOPT_g
371: /* Check that the permutations are complementary */
372: for(i = 0; i < N; i++) {
373: if (i != aomap->appPerm[aomap->petscPerm[i]])
374: SETERRQ(PETSC_ERR_PLIB, "Invalid ordering");
375: }
376: #endif
377: /* Cleanup */
378: if (!mypetsc) {
379: PetscFree(petsc);
380: }
381: PetscFree(allapp);
383: PetscOptionsHasName(PETSC_NULL, "-ao_view", &opt);
384: if (opt == PETSC_TRUE) {
385: AOView(ao, PETSC_VIEWER_STDOUT_SELF);
386: }
388: *aoout = ao;
389: return(0);
390: }
394: /*@C
395: AOCreateMappingIS - Creates a basic application ordering using two index sets.
397: Input Parameters:
398: + comm - MPI communicator that is to share AO
399: . isapp - index set that defines an ordering
400: - ispetsc - index set that defines another ordering
402: Output Parameter:
403: . aoout - the new application ordering
405: Options Database Key:
406: $ -ao_view : call AOView() at the conclusion of AOCreateMappingIS()
408: Level: beginner
410: .keywords: AO, create
411: .seealso: AOCreateBasic(), AOCreateMapping(), AODestroy()
412: @*/
413: int AOCreateMappingIS(IS isapp, IS ispetsc, AO *aoout)
414: {
415: MPI_Comm comm;
416: int *mypetsc, *myapp;
417: int napp, npetsc;
418: int ierr;
421: PetscObjectGetComm((PetscObject) isapp, &comm);
422: ISGetSize(isapp, &napp);
423: if (ispetsc) {
424: ISGetSize(ispetsc, &npetsc);
425: if (napp != npetsc) SETERRQ(PETSC_ERR_ARG_SIZ, "Local IS lengths must match");
426: ISGetIndices(ispetsc, &mypetsc);
427: }
428: ISGetIndices(isapp, &myapp);
430: AOCreateMapping(comm, napp, myapp, mypetsc, aoout);
432: ISRestoreIndices(isapp, &myapp);
433: if (ispetsc) {
434: ISRestoreIndices(ispetsc, &mypetsc);
435: }
436: return(0);
437: }