xref: /petsc/src/vec/is/ao/impls/mapping/aomapping.c (revision 487a658c8b32ba712a1dc8280daad2fd70c1dcd9)
11447629fSBarry Smith 
21447629fSBarry Smith /*
31447629fSBarry Smith   These AO application ordering routines do not require that the input
41447629fSBarry Smith   be a permutation, but merely a 1-1 mapping. This implementation still
51447629fSBarry Smith   keeps the entire ordering on each processor.
61447629fSBarry Smith */
71447629fSBarry Smith 
81447629fSBarry Smith #include <../src/vec/is/ao/aoimpl.h>          /*I  "petscao.h" I*/
91447629fSBarry Smith 
101447629fSBarry Smith typedef struct {
111447629fSBarry Smith   PetscInt N;
121447629fSBarry Smith   PetscInt *app;       /* app[i] is the partner for petsc[appPerm[i]] */
131447629fSBarry Smith   PetscInt *appPerm;
141447629fSBarry Smith   PetscInt *petsc;     /* petsc[j] is the partner for app[petscPerm[j]] */
151447629fSBarry Smith   PetscInt *petscPerm;
161447629fSBarry Smith } AO_Mapping;
171447629fSBarry Smith 
181447629fSBarry Smith PetscErrorCode AODestroy_Mapping(AO ao)
191447629fSBarry Smith {
201447629fSBarry Smith   AO_Mapping     *aomap = (AO_Mapping*) ao->data;
211447629fSBarry Smith   PetscErrorCode ierr;
221447629fSBarry Smith 
231447629fSBarry Smith   PetscFunctionBegin;
241447629fSBarry Smith   ierr = PetscFree4(aomap->app,aomap->appPerm,aomap->petsc,aomap->petscPerm);CHKERRQ(ierr);
251447629fSBarry Smith   ierr = PetscFree(aomap);CHKERRQ(ierr);
261447629fSBarry Smith   PetscFunctionReturn(0);
271447629fSBarry Smith }
281447629fSBarry Smith 
291447629fSBarry Smith PetscErrorCode AOView_Mapping(AO ao, PetscViewer viewer)
301447629fSBarry Smith {
311447629fSBarry Smith   AO_Mapping     *aomap = (AO_Mapping*) ao->data;
321447629fSBarry Smith   PetscMPIInt    rank;
331447629fSBarry Smith   PetscInt       i;
341447629fSBarry Smith   PetscBool      iascii;
351447629fSBarry Smith   PetscErrorCode ierr;
361447629fSBarry Smith 
371447629fSBarry Smith   PetscFunctionBegin;
381447629fSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)ao), &rank);CHKERRQ(ierr);
391447629fSBarry Smith   if (rank) PetscFunctionReturn(0);
401447629fSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
411447629fSBarry Smith   if (iascii) {
421447629fSBarry Smith     PetscViewerASCIIPrintf(viewer, "Number of elements in ordering %D\n", aomap->N);
431447629fSBarry Smith     PetscViewerASCIIPrintf(viewer, "   App.   PETSc\n");
441447629fSBarry Smith     for (i = 0; i < aomap->N; i++) {
451447629fSBarry Smith       PetscViewerASCIIPrintf(viewer, "%D   %D    %D\n", i, aomap->app[i], aomap->petsc[aomap->appPerm[i]]);
461447629fSBarry Smith     }
471447629fSBarry Smith   }
481447629fSBarry Smith   PetscFunctionReturn(0);
491447629fSBarry Smith }
501447629fSBarry Smith 
511447629fSBarry Smith PetscErrorCode AOPetscToApplication_Mapping(AO ao, PetscInt n, PetscInt *ia)
521447629fSBarry Smith {
531447629fSBarry Smith   AO_Mapping *aomap = (AO_Mapping*) ao->data;
541447629fSBarry Smith   PetscInt   *app   = aomap->app;
551447629fSBarry Smith   PetscInt   *petsc = aomap->petsc;
561447629fSBarry Smith   PetscInt   *perm  = aomap->petscPerm;
571447629fSBarry Smith   PetscInt   N      = aomap->N;
581447629fSBarry Smith   PetscInt   low, high, mid=0;
591447629fSBarry Smith   PetscInt   idex;
601447629fSBarry Smith   PetscInt   i;
611447629fSBarry Smith 
621447629fSBarry Smith   /* It would be possible to use a single bisection search, which
631447629fSBarry Smith      recursively divided the indices to be converted, and searched
641447629fSBarry Smith      partitions which contained an index. This would result in
651447629fSBarry Smith      better running times if indices are clustered.
661447629fSBarry Smith   */
671447629fSBarry Smith   PetscFunctionBegin;
681447629fSBarry Smith   for (i = 0; i < n; i++) {
691447629fSBarry Smith     idex = ia[i];
701447629fSBarry Smith     if (idex < 0) continue;
711447629fSBarry Smith     /* Use bisection since the array is sorted */
721447629fSBarry Smith     low  = 0;
731447629fSBarry Smith     high = N - 1;
741447629fSBarry Smith     while (low <= high) {
751447629fSBarry Smith       mid = (low + high)/2;
761447629fSBarry Smith       if (idex == petsc[mid]) break;
771447629fSBarry Smith       else if (idex < petsc[mid]) high = mid - 1;
781447629fSBarry Smith       else low = mid + 1;
791447629fSBarry Smith     }
801447629fSBarry Smith     if (low > high) ia[i] = -1;
811447629fSBarry Smith     else            ia[i] = app[perm[mid]];
821447629fSBarry Smith   }
831447629fSBarry Smith   PetscFunctionReturn(0);
841447629fSBarry Smith }
851447629fSBarry Smith 
861447629fSBarry Smith PetscErrorCode AOApplicationToPetsc_Mapping(AO ao, PetscInt n, PetscInt *ia)
871447629fSBarry Smith {
881447629fSBarry Smith   AO_Mapping *aomap = (AO_Mapping*) ao->data;
891447629fSBarry Smith   PetscInt   *app   = aomap->app;
901447629fSBarry Smith   PetscInt   *petsc = aomap->petsc;
911447629fSBarry Smith   PetscInt   *perm  = aomap->appPerm;
921447629fSBarry Smith   PetscInt   N      = aomap->N;
931447629fSBarry Smith   PetscInt   low, high, mid=0;
941447629fSBarry Smith   PetscInt   idex;
951447629fSBarry Smith   PetscInt   i;
961447629fSBarry Smith 
971447629fSBarry Smith   /* It would be possible to use a single bisection search, which
981447629fSBarry Smith      recursively divided the indices to be converted, and searched
991447629fSBarry Smith      partitions which contained an index. This would result in
1001447629fSBarry Smith      better running times if indices are clustered.
1011447629fSBarry Smith   */
1021447629fSBarry Smith   PetscFunctionBegin;
1031447629fSBarry Smith   for (i = 0; i < n; i++) {
1041447629fSBarry Smith     idex = ia[i];
1051447629fSBarry Smith     if (idex < 0) continue;
1061447629fSBarry Smith     /* Use bisection since the array is sorted */
1071447629fSBarry Smith     low  = 0;
1081447629fSBarry Smith     high = N - 1;
1091447629fSBarry Smith     while (low <= high) {
1101447629fSBarry Smith       mid = (low + high)/2;
1111447629fSBarry Smith       if (idex == app[mid]) break;
1121447629fSBarry Smith       else if (idex < app[mid]) high = mid - 1;
1131447629fSBarry Smith       else low = mid + 1;
1141447629fSBarry Smith     }
1151447629fSBarry Smith     if (low > high) ia[i] = -1;
1161447629fSBarry Smith     else            ia[i] = petsc[perm[mid]];
1171447629fSBarry Smith   }
1181447629fSBarry Smith   PetscFunctionReturn(0);
1191447629fSBarry Smith }
1201447629fSBarry Smith 
1211447629fSBarry Smith static struct _AOOps AOps = {AOView_Mapping,
1221447629fSBarry Smith                              AODestroy_Mapping,
1231447629fSBarry Smith                              AOPetscToApplication_Mapping,
1241447629fSBarry Smith                              AOApplicationToPetsc_Mapping,
1251447629fSBarry Smith                              NULL,
1261447629fSBarry Smith                              NULL,
1271447629fSBarry Smith                              NULL,
1281447629fSBarry Smith                              NULL};
1291447629fSBarry Smith 
1301447629fSBarry Smith /*@C
1311447629fSBarry Smith   AOMappingHasApplicationIndex - Searches for the supplied application index.
1321447629fSBarry Smith 
1331447629fSBarry Smith   Input Parameters:
1341447629fSBarry Smith + ao       - The AOMapping
1351447629fSBarry Smith - index    - The application index
1361447629fSBarry Smith 
1371447629fSBarry Smith   Output Parameter:
1381447629fSBarry Smith . hasIndex - Flag is PETSC_TRUE if the index exists
1391447629fSBarry Smith 
1401447629fSBarry Smith   Level: intermediate
1411447629fSBarry Smith 
1421447629fSBarry Smith .keywords: AO, index
1431447629fSBarry Smith .seealso: AOMappingHasPetscIndex(), AOCreateMapping()
1441447629fSBarry Smith @*/
1451447629fSBarry Smith PetscErrorCode  AOMappingHasApplicationIndex(AO ao, PetscInt idex, PetscBool  *hasIndex)
1461447629fSBarry Smith {
1471447629fSBarry Smith   AO_Mapping *aomap;
1481447629fSBarry Smith   PetscInt   *app;
1491447629fSBarry Smith   PetscInt   low, high, mid;
1501447629fSBarry Smith 
1511447629fSBarry Smith   PetscFunctionBegin;
1521447629fSBarry Smith   PetscValidHeaderSpecific(ao, AO_CLASSID,1);
1531447629fSBarry Smith   PetscValidPointer(hasIndex,3);
1541447629fSBarry Smith   aomap = (AO_Mapping*) ao->data;
1551447629fSBarry Smith   app   = aomap->app;
1561447629fSBarry Smith   /* Use bisection since the array is sorted */
1571447629fSBarry Smith   low  = 0;
1581447629fSBarry Smith   high = aomap->N - 1;
1591447629fSBarry Smith   while (low <= high) {
1601447629fSBarry Smith     mid = (low + high)/2;
1611447629fSBarry Smith     if (idex == app[mid]) break;
1621447629fSBarry Smith     else if (idex < app[mid]) high = mid - 1;
1631447629fSBarry Smith     else low = mid + 1;
1641447629fSBarry Smith   }
1651447629fSBarry Smith   if (low > high) *hasIndex = PETSC_FALSE;
1661447629fSBarry Smith   else *hasIndex = PETSC_TRUE;
1671447629fSBarry Smith   PetscFunctionReturn(0);
1681447629fSBarry Smith }
1691447629fSBarry Smith 
1701447629fSBarry Smith /*@
1711447629fSBarry Smith   AOMappingHasPetscIndex - Searches for the supplied petsc index.
1721447629fSBarry Smith 
1731447629fSBarry Smith   Input Parameters:
1741447629fSBarry Smith + ao       - The AOMapping
1751447629fSBarry Smith - index    - The petsc index
1761447629fSBarry Smith 
1771447629fSBarry Smith   Output Parameter:
1781447629fSBarry Smith . hasIndex - Flag is PETSC_TRUE if the index exists
1791447629fSBarry Smith 
1801447629fSBarry Smith   Level: intermediate
1811447629fSBarry Smith 
1821447629fSBarry Smith .keywords: AO, index
1831447629fSBarry Smith .seealso: AOMappingHasApplicationIndex(), AOCreateMapping()
1841447629fSBarry Smith @*/
1851447629fSBarry Smith PetscErrorCode  AOMappingHasPetscIndex(AO ao, PetscInt idex, PetscBool  *hasIndex)
1861447629fSBarry Smith {
1871447629fSBarry Smith   AO_Mapping *aomap;
1881447629fSBarry Smith   PetscInt   *petsc;
1891447629fSBarry Smith   PetscInt   low, high, mid;
1901447629fSBarry Smith 
1911447629fSBarry Smith   PetscFunctionBegin;
1921447629fSBarry Smith   PetscValidHeaderSpecific(ao, AO_CLASSID,1);
1931447629fSBarry Smith   PetscValidPointer(hasIndex,3);
1941447629fSBarry Smith   aomap = (AO_Mapping*) ao->data;
1951447629fSBarry Smith   petsc = aomap->petsc;
1961447629fSBarry Smith   /* Use bisection since the array is sorted */
1971447629fSBarry Smith   low  = 0;
1981447629fSBarry Smith   high = aomap->N - 1;
1991447629fSBarry Smith   while (low <= high) {
2001447629fSBarry Smith     mid = (low + high)/2;
2011447629fSBarry Smith     if (idex == petsc[mid]) break;
2021447629fSBarry Smith     else if (idex < petsc[mid]) high = mid - 1;
2031447629fSBarry Smith     else low = mid + 1;
2041447629fSBarry Smith   }
2051447629fSBarry Smith   if (low > high) *hasIndex = PETSC_FALSE;
2061447629fSBarry Smith   else *hasIndex = PETSC_TRUE;
2071447629fSBarry Smith   PetscFunctionReturn(0);
2081447629fSBarry Smith }
2091447629fSBarry Smith 
2101447629fSBarry Smith /*@C
2111447629fSBarry Smith   AOCreateMapping - Creates a basic application mapping using two integer arrays.
2121447629fSBarry Smith 
2131447629fSBarry Smith   Input Parameters:
2141447629fSBarry Smith + comm    - MPI communicator that is to share AO
2151447629fSBarry Smith . napp    - size of integer arrays
2161447629fSBarry Smith . myapp   - integer array that defines an ordering
2171447629fSBarry Smith - mypetsc - integer array that defines another ordering (may be NULL to indicate the identity ordering)
2181447629fSBarry Smith 
2191447629fSBarry Smith   Output Parameter:
2201447629fSBarry Smith . aoout   - the new application mapping
2211447629fSBarry Smith 
2221447629fSBarry Smith   Options Database Key:
223e1bc860dSBarry Smith . -ao_view : call AOView() at the conclusion of AOCreateMapping()
2241447629fSBarry Smith 
2251447629fSBarry Smith   Level: beginner
2261447629fSBarry Smith 
22795452b02SPatrick Sanan     Notes:
22895452b02SPatrick Sanan     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.
2291447629fSBarry Smith        Use AOCreateBasic() or AOCreateBasicIS() if they do not have holes for better performance.
2301447629fSBarry Smith 
2311447629fSBarry Smith .keywords: AO, create
2321447629fSBarry Smith .seealso: AOCreateBasic(), AOCreateBasic(), AOCreateMappingIS(), AODestroy()
2331447629fSBarry Smith @*/
2341447629fSBarry Smith PetscErrorCode  AOCreateMapping(MPI_Comm comm,PetscInt napp,const PetscInt myapp[],const PetscInt mypetsc[],AO *aoout)
2351447629fSBarry Smith {
2361447629fSBarry Smith   AO             ao;
2371447629fSBarry Smith   AO_Mapping     *aomap;
2381447629fSBarry Smith   PetscInt       *allpetsc,  *allapp;
2391447629fSBarry Smith   PetscInt       *petscPerm, *appPerm;
2401447629fSBarry Smith   PetscInt       *petsc;
2411447629fSBarry Smith   PetscMPIInt    size, rank,*lens, *disp,nnapp;
2421447629fSBarry Smith   PetscInt       N, start;
2431447629fSBarry Smith   PetscInt       i;
2441447629fSBarry Smith   PetscErrorCode ierr;
2451447629fSBarry Smith 
2461447629fSBarry Smith   PetscFunctionBegin;
2471447629fSBarry Smith   PetscValidPointer(aoout,5);
2481447629fSBarry Smith   *aoout = 0;
249607a6623SBarry Smith   ierr = AOInitializePackage();CHKERRQ(ierr);
2501447629fSBarry Smith 
25173107ff1SLisandro Dalcin   ierr     = PetscHeaderCreate(ao, AO_CLASSID, "AO", "Application Ordering", "AO", comm, AODestroy, AOView);CHKERRQ(ierr);
252b00a9115SJed Brown   ierr     = PetscNewLog(ao,&aomap);CHKERRQ(ierr);
2531447629fSBarry Smith   ierr     = PetscMemcpy(ao->ops, &AOps, sizeof(AOps));CHKERRQ(ierr);
2541447629fSBarry Smith   ao->data = (void*) aomap;
2551447629fSBarry Smith 
2561447629fSBarry Smith   /* transmit all lengths to all processors */
2571447629fSBarry Smith   ierr  = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
2581447629fSBarry Smith   ierr  = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
259dcca6d9dSJed Brown   ierr  = PetscMalloc2(size, &lens,size,&disp);CHKERRQ(ierr);
2601447629fSBarry Smith   nnapp = napp;
2611447629fSBarry Smith   ierr  = MPI_Allgather(&nnapp, 1, MPI_INT, lens, 1, MPI_INT, comm);CHKERRQ(ierr);
2621447629fSBarry Smith   N     = 0;
2631447629fSBarry Smith   for (i = 0; i < size; i++) {
2641447629fSBarry Smith     disp[i] = N;
2651447629fSBarry Smith     N      += lens[i];
2661447629fSBarry Smith   }
2671447629fSBarry Smith   aomap->N = N;
2681447629fSBarry Smith   ao->N    = N;
2691447629fSBarry Smith   ao->n    = N;
2701447629fSBarry Smith 
2711447629fSBarry Smith   /* If mypetsc is 0 then use "natural" numbering */
2721447629fSBarry Smith   if (!mypetsc) {
2731447629fSBarry Smith     start = disp[rank];
274854ce69bSBarry Smith     ierr  = PetscMalloc1(napp+1, &petsc);CHKERRQ(ierr);
2751447629fSBarry Smith     for (i = 0; i < napp; i++) petsc[i] = start + i;
2761447629fSBarry Smith   } else {
2771447629fSBarry Smith     petsc = (PetscInt*)mypetsc;
2781447629fSBarry Smith   }
2791447629fSBarry Smith 
2801447629fSBarry Smith   /* get all indices on all processors */
281dcca6d9dSJed Brown   ierr = PetscMalloc4(N, &allapp,N,&appPerm,N,&allpetsc,N,&petscPerm);CHKERRQ(ierr);
2821447629fSBarry Smith   ierr = MPI_Allgatherv((void*)myapp, napp, MPIU_INT, allapp,   lens, disp, MPIU_INT, comm);CHKERRQ(ierr);
2831447629fSBarry Smith   ierr = MPI_Allgatherv((void*)petsc, napp, MPIU_INT, allpetsc, lens, disp, MPIU_INT, comm);CHKERRQ(ierr);
2841447629fSBarry Smith   ierr = PetscFree2(lens,disp);CHKERRQ(ierr);
2851447629fSBarry Smith 
2861447629fSBarry Smith   /* generate a list of application and PETSc node numbers */
287dcca6d9dSJed Brown   ierr = PetscMalloc4(N, &aomap->app,N,&aomap->appPerm,N,&aomap->petsc,N,&aomap->petscPerm);CHKERRQ(ierr);
2883bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)ao, 4*N * sizeof(PetscInt));CHKERRQ(ierr);
2891447629fSBarry Smith   for (i = 0; i < N; i++) {
2901447629fSBarry Smith     appPerm[i]   = i;
2911447629fSBarry Smith     petscPerm[i] = i;
2921447629fSBarry Smith   }
2931447629fSBarry Smith   ierr = PetscSortIntWithPermutation(N, allpetsc, petscPerm);CHKERRQ(ierr);
2941447629fSBarry Smith   ierr = PetscSortIntWithPermutation(N, allapp,   appPerm);CHKERRQ(ierr);
2951447629fSBarry Smith   /* Form sorted arrays of indices */
2961447629fSBarry Smith   for (i = 0; i < N; i++) {
2971447629fSBarry Smith     aomap->app[i]   = allapp[appPerm[i]];
2981447629fSBarry Smith     aomap->petsc[i] = allpetsc[petscPerm[i]];
2991447629fSBarry Smith   }
3001447629fSBarry Smith   /* Invert petscPerm[] into aomap->petscPerm[] */
3011447629fSBarry Smith   for (i = 0; i < N; i++) aomap->petscPerm[petscPerm[i]] = i;
3021447629fSBarry Smith 
3031447629fSBarry Smith   /* Form map between aomap->app[] and aomap->petsc[] */
3041447629fSBarry Smith   for (i = 0; i < N; i++) aomap->appPerm[i] = aomap->petscPerm[appPerm[i]];
3051447629fSBarry Smith 
3061447629fSBarry Smith   /* Invert appPerm[] into allapp[] */
3071447629fSBarry Smith   for (i = 0; i < N; i++) allapp[appPerm[i]] = i;
3081447629fSBarry Smith 
3091447629fSBarry Smith   /* Form map between aomap->petsc[] and aomap->app[] */
3101447629fSBarry Smith   for (i = 0; i < N; i++) aomap->petscPerm[i] = allapp[petscPerm[i]];
3111447629fSBarry Smith 
3121447629fSBarry Smith #if defined(PETSC_USE_DEBUG)
3131447629fSBarry Smith   /* Check that the permutations are complementary */
3141447629fSBarry Smith   for (i = 0; i < N; i++) {
3151447629fSBarry Smith     if (i != aomap->appPerm[aomap->petscPerm[i]]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Invalid ordering");
3161447629fSBarry Smith   }
3171447629fSBarry Smith #endif
3181447629fSBarry Smith   /* Cleanup */
3191447629fSBarry Smith   if (!mypetsc) {
3201447629fSBarry Smith     ierr = PetscFree(petsc);CHKERRQ(ierr);
3211447629fSBarry Smith   }
3221447629fSBarry Smith   ierr = PetscFree4(allapp,appPerm,allpetsc,petscPerm);CHKERRQ(ierr);
3231447629fSBarry Smith 
324817ea411SJed Brown   ierr = AOViewFromOptions(ao,NULL,"-ao_view");CHKERRQ(ierr);
3251447629fSBarry Smith 
3261447629fSBarry Smith   *aoout = ao;
3271447629fSBarry Smith   PetscFunctionReturn(0);
3281447629fSBarry Smith }
3291447629fSBarry Smith 
330*487a658cSBarry Smith /*@
3311447629fSBarry Smith   AOCreateMappingIS - Creates a basic application ordering using two index sets.
3321447629fSBarry Smith 
3331447629fSBarry Smith   Input Parameters:
3341447629fSBarry Smith + comm    - MPI communicator that is to share AO
3351447629fSBarry Smith . isapp   - index set that defines an ordering
3361447629fSBarry Smith - ispetsc - index set that defines another ordering, maybe NULL for identity IS
3371447629fSBarry Smith 
3381447629fSBarry Smith   Output Parameter:
3391447629fSBarry Smith . aoout   - the new application ordering
3401447629fSBarry Smith 
3411447629fSBarry Smith   Options Database Key:
342e1bc860dSBarry Smith . -ao_view : call AOView() at the conclusion of AOCreateMappingIS()
3431447629fSBarry Smith 
3441447629fSBarry Smith   Level: beginner
3451447629fSBarry Smith 
34695452b02SPatrick Sanan     Notes:
34795452b02SPatrick Sanan     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.
3481447629fSBarry Smith        Use AOCreateBasic() or AOCreateBasicIS() if they do not have holes for better performance.
3491447629fSBarry Smith 
3501447629fSBarry Smith .keywords: AO, create
3511447629fSBarry Smith .seealso: AOCreateBasic(), AOCreateMapping(), AODestroy()
3521447629fSBarry Smith @*/
3531447629fSBarry Smith PetscErrorCode  AOCreateMappingIS(IS isapp, IS ispetsc, AO *aoout)
3541447629fSBarry Smith {
3551447629fSBarry Smith   MPI_Comm       comm;
3561447629fSBarry Smith   const PetscInt *mypetsc, *myapp;
3571447629fSBarry Smith   PetscInt       napp, npetsc;
3581447629fSBarry Smith   PetscErrorCode ierr;
3591447629fSBarry Smith 
3601447629fSBarry Smith   PetscFunctionBegin;
3611447629fSBarry Smith   ierr = PetscObjectGetComm((PetscObject) isapp, &comm);CHKERRQ(ierr);
3621447629fSBarry Smith   ierr = ISGetLocalSize(isapp, &napp);CHKERRQ(ierr);
3631447629fSBarry Smith   if (ispetsc) {
3641447629fSBarry Smith     ierr = ISGetLocalSize(ispetsc, &npetsc);CHKERRQ(ierr);
3651447629fSBarry Smith     if (napp != npetsc) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "Local IS lengths must match");
3661447629fSBarry Smith     ierr = ISGetIndices(ispetsc, &mypetsc);CHKERRQ(ierr);
3671447629fSBarry Smith   } else {
3681447629fSBarry Smith     mypetsc = NULL;
3691447629fSBarry Smith   }
3701447629fSBarry Smith   ierr = ISGetIndices(isapp, &myapp);CHKERRQ(ierr);
3711447629fSBarry Smith 
3721447629fSBarry Smith   ierr = AOCreateMapping(comm, napp, myapp, mypetsc, aoout);CHKERRQ(ierr);
3731447629fSBarry Smith 
3741447629fSBarry Smith   ierr = ISRestoreIndices(isapp, &myapp);CHKERRQ(ierr);
3751447629fSBarry Smith   if (ispetsc) {
3761447629fSBarry Smith     ierr = ISRestoreIndices(ispetsc, &mypetsc);CHKERRQ(ierr);
3771447629fSBarry Smith   }
3781447629fSBarry Smith   PetscFunctionReturn(0);
3791447629fSBarry Smith }
380