xref: /petsc/src/vec/is/ao/impls/mapping/aomapping.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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 
18*d71ae5a4SJacob Faibussowitsch PetscErrorCode AODestroy_Mapping(AO ao)
19*d71ae5a4SJacob Faibussowitsch {
201447629fSBarry Smith   AO_Mapping *aomap = (AO_Mapping *)ao->data;
211447629fSBarry Smith 
221447629fSBarry Smith   PetscFunctionBegin;
239566063dSJacob Faibussowitsch   PetscCall(PetscFree4(aomap->app, aomap->appPerm, aomap->petsc, aomap->petscPerm));
249566063dSJacob Faibussowitsch   PetscCall(PetscFree(aomap));
251447629fSBarry Smith   PetscFunctionReturn(0);
261447629fSBarry Smith }
271447629fSBarry Smith 
28*d71ae5a4SJacob Faibussowitsch PetscErrorCode AOView_Mapping(AO ao, PetscViewer viewer)
29*d71ae5a4SJacob Faibussowitsch {
301447629fSBarry Smith   AO_Mapping *aomap = (AO_Mapping *)ao->data;
311447629fSBarry Smith   PetscMPIInt rank;
321447629fSBarry Smith   PetscInt    i;
331447629fSBarry Smith   PetscBool   iascii;
341447629fSBarry Smith 
351447629fSBarry Smith   PetscFunctionBegin;
369566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)ao), &rank));
371447629fSBarry Smith   if (rank) PetscFunctionReturn(0);
389566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
391447629fSBarry Smith   if (iascii) {
402abc8c78SJacob Faibussowitsch     PetscViewerASCIIPrintf(viewer, "Number of elements in ordering %" PetscInt_FMT "\n", aomap->N);
411447629fSBarry Smith     PetscViewerASCIIPrintf(viewer, "   App.   PETSc\n");
42ad540459SPierre Jolivet     for (i = 0; i < aomap->N; i++) PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "   %" PetscInt_FMT "    %" PetscInt_FMT "\n", i, aomap->app[i], aomap->petsc[aomap->appPerm[i]]);
431447629fSBarry Smith   }
441447629fSBarry Smith   PetscFunctionReturn(0);
451447629fSBarry Smith }
461447629fSBarry Smith 
47*d71ae5a4SJacob Faibussowitsch PetscErrorCode AOPetscToApplication_Mapping(AO ao, PetscInt n, PetscInt *ia)
48*d71ae5a4SJacob Faibussowitsch {
491447629fSBarry Smith   AO_Mapping *aomap = (AO_Mapping *)ao->data;
501447629fSBarry Smith   PetscInt   *app   = aomap->app;
511447629fSBarry Smith   PetscInt   *petsc = aomap->petsc;
521447629fSBarry Smith   PetscInt   *perm  = aomap->petscPerm;
531447629fSBarry Smith   PetscInt    N     = aomap->N;
541447629fSBarry Smith   PetscInt    low, high, mid = 0;
551447629fSBarry Smith   PetscInt    idex;
561447629fSBarry Smith   PetscInt    i;
571447629fSBarry Smith 
581447629fSBarry Smith   /* It would be possible to use a single bisection search, which
591447629fSBarry Smith      recursively divided the indices to be converted, and searched
601447629fSBarry Smith      partitions which contained an index. This would result in
611447629fSBarry Smith      better running times if indices are clustered.
621447629fSBarry Smith   */
631447629fSBarry Smith   PetscFunctionBegin;
641447629fSBarry Smith   for (i = 0; i < n; i++) {
651447629fSBarry Smith     idex = ia[i];
661447629fSBarry Smith     if (idex < 0) continue;
671447629fSBarry Smith     /* Use bisection since the array is sorted */
681447629fSBarry Smith     low  = 0;
691447629fSBarry Smith     high = N - 1;
701447629fSBarry Smith     while (low <= high) {
711447629fSBarry Smith       mid = (low + high) / 2;
721447629fSBarry Smith       if (idex == petsc[mid]) break;
731447629fSBarry Smith       else if (idex < petsc[mid]) high = mid - 1;
741447629fSBarry Smith       else low = mid + 1;
751447629fSBarry Smith     }
761447629fSBarry Smith     if (low > high) ia[i] = -1;
771447629fSBarry Smith     else ia[i] = app[perm[mid]];
781447629fSBarry Smith   }
791447629fSBarry Smith   PetscFunctionReturn(0);
801447629fSBarry Smith }
811447629fSBarry Smith 
82*d71ae5a4SJacob Faibussowitsch PetscErrorCode AOApplicationToPetsc_Mapping(AO ao, PetscInt n, PetscInt *ia)
83*d71ae5a4SJacob Faibussowitsch {
841447629fSBarry Smith   AO_Mapping *aomap = (AO_Mapping *)ao->data;
851447629fSBarry Smith   PetscInt   *app   = aomap->app;
861447629fSBarry Smith   PetscInt   *petsc = aomap->petsc;
871447629fSBarry Smith   PetscInt   *perm  = aomap->appPerm;
881447629fSBarry Smith   PetscInt    N     = aomap->N;
891447629fSBarry Smith   PetscInt    low, high, mid = 0;
901447629fSBarry Smith   PetscInt    idex;
911447629fSBarry Smith   PetscInt    i;
921447629fSBarry Smith 
931447629fSBarry Smith   /* It would be possible to use a single bisection search, which
941447629fSBarry Smith      recursively divided the indices to be converted, and searched
951447629fSBarry Smith      partitions which contained an index. This would result in
961447629fSBarry Smith      better running times if indices are clustered.
971447629fSBarry Smith   */
981447629fSBarry Smith   PetscFunctionBegin;
991447629fSBarry Smith   for (i = 0; i < n; i++) {
1001447629fSBarry Smith     idex = ia[i];
1011447629fSBarry Smith     if (idex < 0) continue;
1021447629fSBarry Smith     /* Use bisection since the array is sorted */
1031447629fSBarry Smith     low  = 0;
1041447629fSBarry Smith     high = N - 1;
1051447629fSBarry Smith     while (low <= high) {
1061447629fSBarry Smith       mid = (low + high) / 2;
1071447629fSBarry Smith       if (idex == app[mid]) break;
1081447629fSBarry Smith       else if (idex < app[mid]) high = mid - 1;
1091447629fSBarry Smith       else low = mid + 1;
1101447629fSBarry Smith     }
1111447629fSBarry Smith     if (low > high) ia[i] = -1;
1121447629fSBarry Smith     else ia[i] = petsc[perm[mid]];
1131447629fSBarry Smith   }
1141447629fSBarry Smith   PetscFunctionReturn(0);
1151447629fSBarry Smith }
1161447629fSBarry Smith 
117267267bdSJacob Faibussowitsch static struct _AOOps AOps = {
118267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(view, AOView_Mapping),
119267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(destroy, AODestroy_Mapping),
120267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(petsctoapplication, AOPetscToApplication_Mapping),
121267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(applicationtopetsc, AOApplicationToPetsc_Mapping),
122267267bdSJacob Faibussowitsch };
1231447629fSBarry Smith 
1241447629fSBarry Smith /*@C
1251447629fSBarry Smith   AOMappingHasApplicationIndex - Searches for the supplied application index.
1261447629fSBarry Smith 
1271447629fSBarry Smith   Input Parameters:
1281447629fSBarry Smith + ao       - The AOMapping
1291447629fSBarry Smith - index    - The application index
1301447629fSBarry Smith 
1311447629fSBarry Smith   Output Parameter:
1321447629fSBarry Smith . hasIndex - Flag is PETSC_TRUE if the index exists
1331447629fSBarry Smith 
1341447629fSBarry Smith   Level: intermediate
1351447629fSBarry Smith 
136db781477SPatrick Sanan .seealso: `AOMappingHasPetscIndex()`, `AOCreateMapping()`
1371447629fSBarry Smith @*/
138*d71ae5a4SJacob Faibussowitsch PetscErrorCode AOMappingHasApplicationIndex(AO ao, PetscInt idex, PetscBool *hasIndex)
139*d71ae5a4SJacob Faibussowitsch {
1401447629fSBarry Smith   AO_Mapping *aomap;
1411447629fSBarry Smith   PetscInt   *app;
1421447629fSBarry Smith   PetscInt    low, high, mid;
1431447629fSBarry Smith 
1441447629fSBarry Smith   PetscFunctionBegin;
1451447629fSBarry Smith   PetscValidHeaderSpecific(ao, AO_CLASSID, 1);
146dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(hasIndex, 3);
1471447629fSBarry Smith   aomap = (AO_Mapping *)ao->data;
1481447629fSBarry Smith   app   = aomap->app;
1491447629fSBarry Smith   /* Use bisection since the array is sorted */
1501447629fSBarry Smith   low  = 0;
1511447629fSBarry Smith   high = aomap->N - 1;
1521447629fSBarry Smith   while (low <= high) {
1531447629fSBarry Smith     mid = (low + high) / 2;
1541447629fSBarry Smith     if (idex == app[mid]) break;
1551447629fSBarry Smith     else if (idex < app[mid]) high = mid - 1;
1561447629fSBarry Smith     else low = mid + 1;
1571447629fSBarry Smith   }
1581447629fSBarry Smith   if (low > high) *hasIndex = PETSC_FALSE;
1591447629fSBarry Smith   else *hasIndex = PETSC_TRUE;
1601447629fSBarry Smith   PetscFunctionReturn(0);
1611447629fSBarry Smith }
1621447629fSBarry Smith 
1631447629fSBarry Smith /*@
1641447629fSBarry Smith   AOMappingHasPetscIndex - Searches for the supplied petsc index.
1651447629fSBarry Smith 
1661447629fSBarry Smith   Input Parameters:
1671447629fSBarry Smith + ao       - The AOMapping
1681447629fSBarry Smith - index    - The petsc index
1691447629fSBarry Smith 
1701447629fSBarry Smith   Output Parameter:
1711447629fSBarry Smith . hasIndex - Flag is PETSC_TRUE if the index exists
1721447629fSBarry Smith 
1731447629fSBarry Smith   Level: intermediate
1741447629fSBarry Smith 
175db781477SPatrick Sanan .seealso: `AOMappingHasApplicationIndex()`, `AOCreateMapping()`
1761447629fSBarry Smith @*/
177*d71ae5a4SJacob Faibussowitsch PetscErrorCode AOMappingHasPetscIndex(AO ao, PetscInt idex, PetscBool *hasIndex)
178*d71ae5a4SJacob Faibussowitsch {
1791447629fSBarry Smith   AO_Mapping *aomap;
1801447629fSBarry Smith   PetscInt   *petsc;
1811447629fSBarry Smith   PetscInt    low, high, mid;
1821447629fSBarry Smith 
1831447629fSBarry Smith   PetscFunctionBegin;
1841447629fSBarry Smith   PetscValidHeaderSpecific(ao, AO_CLASSID, 1);
185dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(hasIndex, 3);
1861447629fSBarry Smith   aomap = (AO_Mapping *)ao->data;
1871447629fSBarry Smith   petsc = aomap->petsc;
1881447629fSBarry Smith   /* Use bisection since the array is sorted */
1891447629fSBarry Smith   low  = 0;
1901447629fSBarry Smith   high = aomap->N - 1;
1911447629fSBarry Smith   while (low <= high) {
1921447629fSBarry Smith     mid = (low + high) / 2;
1931447629fSBarry Smith     if (idex == petsc[mid]) break;
1941447629fSBarry Smith     else if (idex < petsc[mid]) high = mid - 1;
1951447629fSBarry Smith     else low = mid + 1;
1961447629fSBarry Smith   }
1971447629fSBarry Smith   if (low > high) *hasIndex = PETSC_FALSE;
1981447629fSBarry Smith   else *hasIndex = PETSC_TRUE;
1991447629fSBarry Smith   PetscFunctionReturn(0);
2001447629fSBarry Smith }
2011447629fSBarry Smith 
2021447629fSBarry Smith /*@C
2031447629fSBarry Smith   AOCreateMapping - Creates a basic application mapping using two integer arrays.
2041447629fSBarry Smith 
2051447629fSBarry Smith   Input Parameters:
2061447629fSBarry Smith + comm    - MPI communicator that is to share AO
2071447629fSBarry Smith . napp    - size of integer arrays
2081447629fSBarry Smith . myapp   - integer array that defines an ordering
2091447629fSBarry Smith - mypetsc - integer array that defines another ordering (may be NULL to indicate the identity ordering)
2101447629fSBarry Smith 
2111447629fSBarry Smith   Output Parameter:
2121447629fSBarry Smith . aoout   - the new application mapping
2131447629fSBarry Smith 
2141447629fSBarry Smith   Options Database Key:
215147403d9SBarry Smith . -ao_view - call AOView() at the conclusion of AOCreateMapping()
2161447629fSBarry Smith 
2171447629fSBarry Smith   Level: beginner
2181447629fSBarry Smith 
21995452b02SPatrick Sanan     Notes:
22095452b02SPatrick 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.
2211447629fSBarry Smith        Use AOCreateBasic() or AOCreateBasicIS() if they do not have holes for better performance.
2221447629fSBarry Smith 
223db781477SPatrick Sanan .seealso: `AOCreateBasic()`, `AOCreateBasic()`, `AOCreateMappingIS()`, `AODestroy()`
2241447629fSBarry Smith @*/
225*d71ae5a4SJacob Faibussowitsch PetscErrorCode AOCreateMapping(MPI_Comm comm, PetscInt napp, const PetscInt myapp[], const PetscInt mypetsc[], AO *aoout)
226*d71ae5a4SJacob Faibussowitsch {
2271447629fSBarry Smith   AO          ao;
2281447629fSBarry Smith   AO_Mapping *aomap;
2291447629fSBarry Smith   PetscInt   *allpetsc, *allapp;
2301447629fSBarry Smith   PetscInt   *petscPerm, *appPerm;
2311447629fSBarry Smith   PetscInt   *petsc;
2321447629fSBarry Smith   PetscMPIInt size, rank, *lens, *disp, nnapp;
2331447629fSBarry Smith   PetscInt    N, start;
2341447629fSBarry Smith   PetscInt    i;
2351447629fSBarry Smith 
2361447629fSBarry Smith   PetscFunctionBegin;
2371447629fSBarry Smith   PetscValidPointer(aoout, 5);
2384c8fdceaSLisandro Dalcin   *aoout = NULL;
2399566063dSJacob Faibussowitsch   PetscCall(AOInitializePackage());
2401447629fSBarry Smith 
2419566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(ao, AO_CLASSID, "AO", "Application Ordering", "AO", comm, AODestroy, AOView));
2424dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&aomap));
2439566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(ao->ops, &AOps, sizeof(AOps)));
2441447629fSBarry Smith   ao->data = (void *)aomap;
2451447629fSBarry Smith 
2461447629fSBarry Smith   /* transmit all lengths to all processors */
2479566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
2489566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2499566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(size, &lens, size, &disp));
2501447629fSBarry Smith   nnapp = napp;
2519566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&nnapp, 1, MPI_INT, lens, 1, MPI_INT, comm));
2521447629fSBarry Smith   N = 0;
2531447629fSBarry Smith   for (i = 0; i < size; i++) {
2541447629fSBarry Smith     disp[i] = N;
2551447629fSBarry Smith     N += lens[i];
2561447629fSBarry Smith   }
2571447629fSBarry Smith   aomap->N = N;
2581447629fSBarry Smith   ao->N    = N;
2591447629fSBarry Smith   ao->n    = N;
2601447629fSBarry Smith 
2611447629fSBarry Smith   /* If mypetsc is 0 then use "natural" numbering */
2621447629fSBarry Smith   if (!mypetsc) {
2631447629fSBarry Smith     start = disp[rank];
2649566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(napp + 1, &petsc));
2651447629fSBarry Smith     for (i = 0; i < napp; i++) petsc[i] = start + i;
2661447629fSBarry Smith   } else {
2671447629fSBarry Smith     petsc = (PetscInt *)mypetsc;
2681447629fSBarry Smith   }
2691447629fSBarry Smith 
2701447629fSBarry Smith   /* get all indices on all processors */
2719566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(N, &allapp, N, &appPerm, N, &allpetsc, N, &petscPerm));
2729566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgatherv((void *)myapp, napp, MPIU_INT, allapp, lens, disp, MPIU_INT, comm));
2739566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgatherv((void *)petsc, napp, MPIU_INT, allpetsc, lens, disp, MPIU_INT, comm));
2749566063dSJacob Faibussowitsch   PetscCall(PetscFree2(lens, disp));
2751447629fSBarry Smith 
2761447629fSBarry Smith   /* generate a list of application and PETSc node numbers */
2779566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(N, &aomap->app, N, &aomap->appPerm, N, &aomap->petsc, N, &aomap->petscPerm));
2781447629fSBarry Smith   for (i = 0; i < N; i++) {
2791447629fSBarry Smith     appPerm[i]   = i;
2801447629fSBarry Smith     petscPerm[i] = i;
2811447629fSBarry Smith   }
2829566063dSJacob Faibussowitsch   PetscCall(PetscSortIntWithPermutation(N, allpetsc, petscPerm));
2839566063dSJacob Faibussowitsch   PetscCall(PetscSortIntWithPermutation(N, allapp, appPerm));
2841447629fSBarry Smith   /* Form sorted arrays of indices */
2851447629fSBarry Smith   for (i = 0; i < N; i++) {
2861447629fSBarry Smith     aomap->app[i]   = allapp[appPerm[i]];
2871447629fSBarry Smith     aomap->petsc[i] = allpetsc[petscPerm[i]];
2881447629fSBarry Smith   }
2891447629fSBarry Smith   /* Invert petscPerm[] into aomap->petscPerm[] */
2901447629fSBarry Smith   for (i = 0; i < N; i++) aomap->petscPerm[petscPerm[i]] = i;
2911447629fSBarry Smith 
2921447629fSBarry Smith   /* Form map between aomap->app[] and aomap->petsc[] */
2931447629fSBarry Smith   for (i = 0; i < N; i++) aomap->appPerm[i] = aomap->petscPerm[appPerm[i]];
2941447629fSBarry Smith 
2951447629fSBarry Smith   /* Invert appPerm[] into allapp[] */
2961447629fSBarry Smith   for (i = 0; i < N; i++) allapp[appPerm[i]] = i;
2971447629fSBarry Smith 
2981447629fSBarry Smith   /* Form map between aomap->petsc[] and aomap->app[] */
2991447629fSBarry Smith   for (i = 0; i < N; i++) aomap->petscPerm[i] = allapp[petscPerm[i]];
3001447629fSBarry Smith 
30176bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
3021447629fSBarry Smith     /* Check that the permutations are complementary */
303ad540459SPierre Jolivet     for (i = 0; i < N; i++) PetscCheck(i == aomap->appPerm[aomap->petscPerm[i]], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid ordering");
30476bd3646SJed Brown   }
3051447629fSBarry Smith   /* Cleanup */
30648a46eb9SPierre Jolivet   if (!mypetsc) PetscCall(PetscFree(petsc));
3079566063dSJacob Faibussowitsch   PetscCall(PetscFree4(allapp, appPerm, allpetsc, petscPerm));
3081447629fSBarry Smith 
3099566063dSJacob Faibussowitsch   PetscCall(AOViewFromOptions(ao, NULL, "-ao_view"));
3101447629fSBarry Smith 
3111447629fSBarry Smith   *aoout = ao;
3121447629fSBarry Smith   PetscFunctionReturn(0);
3131447629fSBarry Smith }
3141447629fSBarry Smith 
315487a658cSBarry Smith /*@
3161447629fSBarry Smith   AOCreateMappingIS - Creates a basic application ordering using two index sets.
3171447629fSBarry Smith 
3181447629fSBarry Smith   Input Parameters:
3191447629fSBarry Smith + comm    - MPI communicator that is to share AO
3201447629fSBarry Smith . isapp   - index set that defines an ordering
3211447629fSBarry Smith - ispetsc - index set that defines another ordering, maybe NULL for identity IS
3221447629fSBarry Smith 
3231447629fSBarry Smith   Output Parameter:
3241447629fSBarry Smith . aoout   - the new application ordering
3251447629fSBarry Smith 
3261447629fSBarry Smith   Options Database Key:
327147403d9SBarry Smith . -ao_view - call AOView() at the conclusion of AOCreateMappingIS()
3281447629fSBarry Smith 
3291447629fSBarry Smith   Level: beginner
3301447629fSBarry Smith 
33195452b02SPatrick Sanan     Notes:
33295452b02SPatrick 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.
3331447629fSBarry Smith        Use AOCreateBasic() or AOCreateBasicIS() if they do not have holes for better performance.
3341447629fSBarry Smith 
335db781477SPatrick Sanan .seealso: `AOCreateBasic()`, `AOCreateMapping()`, `AODestroy()`
3361447629fSBarry Smith @*/
337*d71ae5a4SJacob Faibussowitsch PetscErrorCode AOCreateMappingIS(IS isapp, IS ispetsc, AO *aoout)
338*d71ae5a4SJacob Faibussowitsch {
3391447629fSBarry Smith   MPI_Comm        comm;
3401447629fSBarry Smith   const PetscInt *mypetsc, *myapp;
3411447629fSBarry Smith   PetscInt        napp, npetsc;
3421447629fSBarry Smith 
3431447629fSBarry Smith   PetscFunctionBegin;
3449566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)isapp, &comm));
3459566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isapp, &napp));
3461447629fSBarry Smith   if (ispetsc) {
3479566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(ispetsc, &npetsc));
34808401ef6SPierre Jolivet     PetscCheck(napp == npetsc, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local IS lengths must match");
3499566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(ispetsc, &mypetsc));
3501447629fSBarry Smith   } else {
3511447629fSBarry Smith     mypetsc = NULL;
3521447629fSBarry Smith   }
3539566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isapp, &myapp));
3541447629fSBarry Smith 
3559566063dSJacob Faibussowitsch   PetscCall(AOCreateMapping(comm, napp, myapp, mypetsc, aoout));
3561447629fSBarry Smith 
3579566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isapp, &myapp));
35848a46eb9SPierre Jolivet   if (ispetsc) PetscCall(ISRestoreIndices(ispetsc, &mypetsc));
3591447629fSBarry Smith   PetscFunctionReturn(0);
3601447629fSBarry Smith }
361