xref: /petsc/src/vec/is/ao/impls/mapping/aomapping.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
189371c9d4SSatish Balay PetscErrorCode AODestroy_Mapping(AO ao) {
191447629fSBarry Smith   AO_Mapping *aomap = (AO_Mapping *)ao->data;
201447629fSBarry Smith 
211447629fSBarry Smith   PetscFunctionBegin;
229566063dSJacob Faibussowitsch   PetscCall(PetscFree4(aomap->app, aomap->appPerm, aomap->petsc, aomap->petscPerm));
239566063dSJacob Faibussowitsch   PetscCall(PetscFree(aomap));
241447629fSBarry Smith   PetscFunctionReturn(0);
251447629fSBarry Smith }
261447629fSBarry Smith 
279371c9d4SSatish Balay PetscErrorCode AOView_Mapping(AO ao, PetscViewer viewer) {
281447629fSBarry Smith   AO_Mapping *aomap = (AO_Mapping *)ao->data;
291447629fSBarry Smith   PetscMPIInt rank;
301447629fSBarry Smith   PetscInt    i;
311447629fSBarry Smith   PetscBool   iascii;
321447629fSBarry Smith 
331447629fSBarry Smith   PetscFunctionBegin;
349566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)ao), &rank));
351447629fSBarry Smith   if (rank) PetscFunctionReturn(0);
369566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
371447629fSBarry Smith   if (iascii) {
382abc8c78SJacob Faibussowitsch     PetscViewerASCIIPrintf(viewer, "Number of elements in ordering %" PetscInt_FMT "\n", aomap->N);
391447629fSBarry Smith     PetscViewerASCIIPrintf(viewer, "   App.   PETSc\n");
409371c9d4SSatish Balay     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]]); }
411447629fSBarry Smith   }
421447629fSBarry Smith   PetscFunctionReturn(0);
431447629fSBarry Smith }
441447629fSBarry Smith 
459371c9d4SSatish Balay PetscErrorCode AOPetscToApplication_Mapping(AO ao, PetscInt n, PetscInt *ia) {
461447629fSBarry Smith   AO_Mapping *aomap = (AO_Mapping *)ao->data;
471447629fSBarry Smith   PetscInt   *app   = aomap->app;
481447629fSBarry Smith   PetscInt   *petsc = aomap->petsc;
491447629fSBarry Smith   PetscInt   *perm  = aomap->petscPerm;
501447629fSBarry Smith   PetscInt    N     = aomap->N;
511447629fSBarry Smith   PetscInt    low, high, mid = 0;
521447629fSBarry Smith   PetscInt    idex;
531447629fSBarry Smith   PetscInt    i;
541447629fSBarry Smith 
551447629fSBarry Smith   /* It would be possible to use a single bisection search, which
561447629fSBarry Smith      recursively divided the indices to be converted, and searched
571447629fSBarry Smith      partitions which contained an index. This would result in
581447629fSBarry Smith      better running times if indices are clustered.
591447629fSBarry Smith   */
601447629fSBarry Smith   PetscFunctionBegin;
611447629fSBarry Smith   for (i = 0; i < n; i++) {
621447629fSBarry Smith     idex = ia[i];
631447629fSBarry Smith     if (idex < 0) continue;
641447629fSBarry Smith     /* Use bisection since the array is sorted */
651447629fSBarry Smith     low  = 0;
661447629fSBarry Smith     high = N - 1;
671447629fSBarry Smith     while (low <= high) {
681447629fSBarry Smith       mid = (low + high) / 2;
691447629fSBarry Smith       if (idex == petsc[mid]) break;
701447629fSBarry Smith       else if (idex < petsc[mid]) high = mid - 1;
711447629fSBarry Smith       else low = mid + 1;
721447629fSBarry Smith     }
731447629fSBarry Smith     if (low > high) ia[i] = -1;
741447629fSBarry Smith     else ia[i] = app[perm[mid]];
751447629fSBarry Smith   }
761447629fSBarry Smith   PetscFunctionReturn(0);
771447629fSBarry Smith }
781447629fSBarry Smith 
799371c9d4SSatish Balay PetscErrorCode AOApplicationToPetsc_Mapping(AO ao, PetscInt n, PetscInt *ia) {
801447629fSBarry Smith   AO_Mapping *aomap = (AO_Mapping *)ao->data;
811447629fSBarry Smith   PetscInt   *app   = aomap->app;
821447629fSBarry Smith   PetscInt   *petsc = aomap->petsc;
831447629fSBarry Smith   PetscInt   *perm  = aomap->appPerm;
841447629fSBarry Smith   PetscInt    N     = aomap->N;
851447629fSBarry Smith   PetscInt    low, high, mid = 0;
861447629fSBarry Smith   PetscInt    idex;
871447629fSBarry Smith   PetscInt    i;
881447629fSBarry Smith 
891447629fSBarry Smith   /* It would be possible to use a single bisection search, which
901447629fSBarry Smith      recursively divided the indices to be converted, and searched
911447629fSBarry Smith      partitions which contained an index. This would result in
921447629fSBarry Smith      better running times if indices are clustered.
931447629fSBarry Smith   */
941447629fSBarry Smith   PetscFunctionBegin;
951447629fSBarry Smith   for (i = 0; i < n; i++) {
961447629fSBarry Smith     idex = ia[i];
971447629fSBarry Smith     if (idex < 0) continue;
981447629fSBarry Smith     /* Use bisection since the array is sorted */
991447629fSBarry Smith     low  = 0;
1001447629fSBarry Smith     high = N - 1;
1011447629fSBarry Smith     while (low <= high) {
1021447629fSBarry Smith       mid = (low + high) / 2;
1031447629fSBarry Smith       if (idex == app[mid]) break;
1041447629fSBarry Smith       else if (idex < app[mid]) high = mid - 1;
1051447629fSBarry Smith       else low = mid + 1;
1061447629fSBarry Smith     }
1071447629fSBarry Smith     if (low > high) ia[i] = -1;
1081447629fSBarry Smith     else ia[i] = petsc[perm[mid]];
1091447629fSBarry Smith   }
1101447629fSBarry Smith   PetscFunctionReturn(0);
1111447629fSBarry Smith }
1121447629fSBarry Smith 
113267267bdSJacob Faibussowitsch static struct _AOOps AOps = {
114267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(view, AOView_Mapping),
115267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(destroy, AODestroy_Mapping),
116267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(petsctoapplication, AOPetscToApplication_Mapping),
117267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(applicationtopetsc, AOApplicationToPetsc_Mapping),
118267267bdSJacob Faibussowitsch };
1191447629fSBarry Smith 
1201447629fSBarry Smith /*@C
1211447629fSBarry Smith   AOMappingHasApplicationIndex - Searches for the supplied application index.
1221447629fSBarry Smith 
1231447629fSBarry Smith   Input Parameters:
1241447629fSBarry Smith + ao       - The AOMapping
1251447629fSBarry Smith - index    - The application index
1261447629fSBarry Smith 
1271447629fSBarry Smith   Output Parameter:
1281447629fSBarry Smith . hasIndex - Flag is PETSC_TRUE if the index exists
1291447629fSBarry Smith 
1301447629fSBarry Smith   Level: intermediate
1311447629fSBarry Smith 
132db781477SPatrick Sanan .seealso: `AOMappingHasPetscIndex()`, `AOCreateMapping()`
1331447629fSBarry Smith @*/
1349371c9d4SSatish Balay PetscErrorCode AOMappingHasApplicationIndex(AO ao, PetscInt idex, PetscBool *hasIndex) {
1351447629fSBarry Smith   AO_Mapping *aomap;
1361447629fSBarry Smith   PetscInt   *app;
1371447629fSBarry Smith   PetscInt    low, high, mid;
1381447629fSBarry Smith 
1391447629fSBarry Smith   PetscFunctionBegin;
1401447629fSBarry Smith   PetscValidHeaderSpecific(ao, AO_CLASSID, 1);
141dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(hasIndex, 3);
1421447629fSBarry Smith   aomap = (AO_Mapping *)ao->data;
1431447629fSBarry Smith   app   = aomap->app;
1441447629fSBarry Smith   /* Use bisection since the array is sorted */
1451447629fSBarry Smith   low   = 0;
1461447629fSBarry Smith   high  = aomap->N - 1;
1471447629fSBarry Smith   while (low <= high) {
1481447629fSBarry Smith     mid = (low + high) / 2;
1491447629fSBarry Smith     if (idex == app[mid]) break;
1501447629fSBarry Smith     else if (idex < app[mid]) high = mid - 1;
1511447629fSBarry Smith     else low = mid + 1;
1521447629fSBarry Smith   }
1531447629fSBarry Smith   if (low > high) *hasIndex = PETSC_FALSE;
1541447629fSBarry Smith   else *hasIndex = PETSC_TRUE;
1551447629fSBarry Smith   PetscFunctionReturn(0);
1561447629fSBarry Smith }
1571447629fSBarry Smith 
1581447629fSBarry Smith /*@
1591447629fSBarry Smith   AOMappingHasPetscIndex - Searches for the supplied petsc index.
1601447629fSBarry Smith 
1611447629fSBarry Smith   Input Parameters:
1621447629fSBarry Smith + ao       - The AOMapping
1631447629fSBarry Smith - index    - The petsc index
1641447629fSBarry Smith 
1651447629fSBarry Smith   Output Parameter:
1661447629fSBarry Smith . hasIndex - Flag is PETSC_TRUE if the index exists
1671447629fSBarry Smith 
1681447629fSBarry Smith   Level: intermediate
1691447629fSBarry Smith 
170db781477SPatrick Sanan .seealso: `AOMappingHasApplicationIndex()`, `AOCreateMapping()`
1711447629fSBarry Smith @*/
1729371c9d4SSatish Balay PetscErrorCode AOMappingHasPetscIndex(AO ao, PetscInt idex, PetscBool *hasIndex) {
1731447629fSBarry Smith   AO_Mapping *aomap;
1741447629fSBarry Smith   PetscInt   *petsc;
1751447629fSBarry Smith   PetscInt    low, high, mid;
1761447629fSBarry Smith 
1771447629fSBarry Smith   PetscFunctionBegin;
1781447629fSBarry Smith   PetscValidHeaderSpecific(ao, AO_CLASSID, 1);
179dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(hasIndex, 3);
1801447629fSBarry Smith   aomap = (AO_Mapping *)ao->data;
1811447629fSBarry Smith   petsc = aomap->petsc;
1821447629fSBarry Smith   /* Use bisection since the array is sorted */
1831447629fSBarry Smith   low   = 0;
1841447629fSBarry Smith   high  = aomap->N - 1;
1851447629fSBarry Smith   while (low <= high) {
1861447629fSBarry Smith     mid = (low + high) / 2;
1871447629fSBarry Smith     if (idex == petsc[mid]) break;
1881447629fSBarry Smith     else if (idex < petsc[mid]) high = mid - 1;
1891447629fSBarry Smith     else low = mid + 1;
1901447629fSBarry Smith   }
1911447629fSBarry Smith   if (low > high) *hasIndex = PETSC_FALSE;
1921447629fSBarry Smith   else *hasIndex = PETSC_TRUE;
1931447629fSBarry Smith   PetscFunctionReturn(0);
1941447629fSBarry Smith }
1951447629fSBarry Smith 
1961447629fSBarry Smith /*@C
1971447629fSBarry Smith   AOCreateMapping - Creates a basic application mapping using two integer arrays.
1981447629fSBarry Smith 
1991447629fSBarry Smith   Input Parameters:
2001447629fSBarry Smith + comm    - MPI communicator that is to share AO
2011447629fSBarry Smith . napp    - size of integer arrays
2021447629fSBarry Smith . myapp   - integer array that defines an ordering
2031447629fSBarry Smith - mypetsc - integer array that defines another ordering (may be NULL to indicate the identity ordering)
2041447629fSBarry Smith 
2051447629fSBarry Smith   Output Parameter:
2061447629fSBarry Smith . aoout   - the new application mapping
2071447629fSBarry Smith 
2081447629fSBarry Smith   Options Database Key:
209147403d9SBarry Smith . -ao_view - call AOView() at the conclusion of AOCreateMapping()
2101447629fSBarry Smith 
2111447629fSBarry Smith   Level: beginner
2121447629fSBarry Smith 
21395452b02SPatrick Sanan     Notes:
21495452b02SPatrick 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.
2151447629fSBarry Smith        Use AOCreateBasic() or AOCreateBasicIS() if they do not have holes for better performance.
2161447629fSBarry Smith 
217db781477SPatrick Sanan .seealso: `AOCreateBasic()`, `AOCreateBasic()`, `AOCreateMappingIS()`, `AODestroy()`
2181447629fSBarry Smith @*/
2199371c9d4SSatish Balay PetscErrorCode AOCreateMapping(MPI_Comm comm, PetscInt napp, const PetscInt myapp[], const PetscInt mypetsc[], AO *aoout) {
2201447629fSBarry Smith   AO          ao;
2211447629fSBarry Smith   AO_Mapping *aomap;
2221447629fSBarry Smith   PetscInt   *allpetsc, *allapp;
2231447629fSBarry Smith   PetscInt   *petscPerm, *appPerm;
2241447629fSBarry Smith   PetscInt   *petsc;
2251447629fSBarry Smith   PetscMPIInt size, rank, *lens, *disp, nnapp;
2261447629fSBarry Smith   PetscInt    N, start;
2271447629fSBarry Smith   PetscInt    i;
2281447629fSBarry Smith 
2291447629fSBarry Smith   PetscFunctionBegin;
2301447629fSBarry Smith   PetscValidPointer(aoout, 5);
2314c8fdceaSLisandro Dalcin   *aoout = NULL;
2329566063dSJacob Faibussowitsch   PetscCall(AOInitializePackage());
2331447629fSBarry Smith 
2349566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(ao, AO_CLASSID, "AO", "Application Ordering", "AO", comm, AODestroy, AOView));
2359566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(ao, &aomap));
2369566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(ao->ops, &AOps, sizeof(AOps)));
2371447629fSBarry Smith   ao->data = (void *)aomap;
2381447629fSBarry Smith 
2391447629fSBarry Smith   /* transmit all lengths to all processors */
2409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
2419566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2429566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(size, &lens, size, &disp));
2431447629fSBarry Smith   nnapp = napp;
2449566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&nnapp, 1, MPI_INT, lens, 1, MPI_INT, comm));
2451447629fSBarry Smith   N = 0;
2461447629fSBarry Smith   for (i = 0; i < size; i++) {
2471447629fSBarry Smith     disp[i] = N;
2481447629fSBarry Smith     N += lens[i];
2491447629fSBarry Smith   }
2501447629fSBarry Smith   aomap->N = N;
2511447629fSBarry Smith   ao->N    = N;
2521447629fSBarry Smith   ao->n    = N;
2531447629fSBarry Smith 
2541447629fSBarry Smith   /* If mypetsc is 0 then use "natural" numbering */
2551447629fSBarry Smith   if (!mypetsc) {
2561447629fSBarry Smith     start = disp[rank];
2579566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(napp + 1, &petsc));
2581447629fSBarry Smith     for (i = 0; i < napp; i++) petsc[i] = start + i;
2591447629fSBarry Smith   } else {
2601447629fSBarry Smith     petsc = (PetscInt *)mypetsc;
2611447629fSBarry Smith   }
2621447629fSBarry Smith 
2631447629fSBarry Smith   /* get all indices on all processors */
2649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(N, &allapp, N, &appPerm, N, &allpetsc, N, &petscPerm));
2659566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgatherv((void *)myapp, napp, MPIU_INT, allapp, lens, disp, MPIU_INT, comm));
2669566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgatherv((void *)petsc, napp, MPIU_INT, allpetsc, lens, disp, MPIU_INT, comm));
2679566063dSJacob Faibussowitsch   PetscCall(PetscFree2(lens, disp));
2681447629fSBarry Smith 
2691447629fSBarry Smith   /* generate a list of application and PETSc node numbers */
2709566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(N, &aomap->app, N, &aomap->appPerm, N, &aomap->petsc, N, &aomap->petscPerm));
2719566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)ao, 4 * N * sizeof(PetscInt)));
2721447629fSBarry Smith   for (i = 0; i < N; i++) {
2731447629fSBarry Smith     appPerm[i]   = i;
2741447629fSBarry Smith     petscPerm[i] = i;
2751447629fSBarry Smith   }
2769566063dSJacob Faibussowitsch   PetscCall(PetscSortIntWithPermutation(N, allpetsc, petscPerm));
2779566063dSJacob Faibussowitsch   PetscCall(PetscSortIntWithPermutation(N, allapp, appPerm));
2781447629fSBarry Smith   /* Form sorted arrays of indices */
2791447629fSBarry Smith   for (i = 0; i < N; i++) {
2801447629fSBarry Smith     aomap->app[i]   = allapp[appPerm[i]];
2811447629fSBarry Smith     aomap->petsc[i] = allpetsc[petscPerm[i]];
2821447629fSBarry Smith   }
2831447629fSBarry Smith   /* Invert petscPerm[] into aomap->petscPerm[] */
2841447629fSBarry Smith   for (i = 0; i < N; i++) aomap->petscPerm[petscPerm[i]] = i;
2851447629fSBarry Smith 
2861447629fSBarry Smith   /* Form map between aomap->app[] and aomap->petsc[] */
2871447629fSBarry Smith   for (i = 0; i < N; i++) aomap->appPerm[i] = aomap->petscPerm[appPerm[i]];
2881447629fSBarry Smith 
2891447629fSBarry Smith   /* Invert appPerm[] into allapp[] */
2901447629fSBarry Smith   for (i = 0; i < N; i++) allapp[appPerm[i]] = i;
2911447629fSBarry Smith 
2921447629fSBarry Smith   /* Form map between aomap->petsc[] and aomap->app[] */
2931447629fSBarry Smith   for (i = 0; i < N; i++) aomap->petscPerm[i] = allapp[petscPerm[i]];
2941447629fSBarry Smith 
29576bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
2961447629fSBarry Smith     /* Check that the permutations are complementary */
2979371c9d4SSatish Balay     for (i = 0; i < N; i++) { PetscCheck(i == aomap->appPerm[aomap->petscPerm[i]], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid ordering"); }
29876bd3646SJed Brown   }
2991447629fSBarry Smith   /* Cleanup */
300*48a46eb9SPierre Jolivet   if (!mypetsc) PetscCall(PetscFree(petsc));
3019566063dSJacob Faibussowitsch   PetscCall(PetscFree4(allapp, appPerm, allpetsc, petscPerm));
3021447629fSBarry Smith 
3039566063dSJacob Faibussowitsch   PetscCall(AOViewFromOptions(ao, NULL, "-ao_view"));
3041447629fSBarry Smith 
3051447629fSBarry Smith   *aoout = ao;
3061447629fSBarry Smith   PetscFunctionReturn(0);
3071447629fSBarry Smith }
3081447629fSBarry Smith 
309487a658cSBarry Smith /*@
3101447629fSBarry Smith   AOCreateMappingIS - Creates a basic application ordering using two index sets.
3111447629fSBarry Smith 
3121447629fSBarry Smith   Input Parameters:
3131447629fSBarry Smith + comm    - MPI communicator that is to share AO
3141447629fSBarry Smith . isapp   - index set that defines an ordering
3151447629fSBarry Smith - ispetsc - index set that defines another ordering, maybe NULL for identity IS
3161447629fSBarry Smith 
3171447629fSBarry Smith   Output Parameter:
3181447629fSBarry Smith . aoout   - the new application ordering
3191447629fSBarry Smith 
3201447629fSBarry Smith   Options Database Key:
321147403d9SBarry Smith . -ao_view - call AOView() at the conclusion of AOCreateMappingIS()
3221447629fSBarry Smith 
3231447629fSBarry Smith   Level: beginner
3241447629fSBarry Smith 
32595452b02SPatrick Sanan     Notes:
32695452b02SPatrick 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.
3271447629fSBarry Smith        Use AOCreateBasic() or AOCreateBasicIS() if they do not have holes for better performance.
3281447629fSBarry Smith 
329db781477SPatrick Sanan .seealso: `AOCreateBasic()`, `AOCreateMapping()`, `AODestroy()`
3301447629fSBarry Smith @*/
3319371c9d4SSatish Balay PetscErrorCode AOCreateMappingIS(IS isapp, IS ispetsc, AO *aoout) {
3321447629fSBarry Smith   MPI_Comm        comm;
3331447629fSBarry Smith   const PetscInt *mypetsc, *myapp;
3341447629fSBarry Smith   PetscInt        napp, npetsc;
3351447629fSBarry Smith 
3361447629fSBarry Smith   PetscFunctionBegin;
3379566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)isapp, &comm));
3389566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isapp, &napp));
3391447629fSBarry Smith   if (ispetsc) {
3409566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(ispetsc, &npetsc));
34108401ef6SPierre Jolivet     PetscCheck(napp == npetsc, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local IS lengths must match");
3429566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(ispetsc, &mypetsc));
3431447629fSBarry Smith   } else {
3441447629fSBarry Smith     mypetsc = NULL;
3451447629fSBarry Smith   }
3469566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isapp, &myapp));
3471447629fSBarry Smith 
3489566063dSJacob Faibussowitsch   PetscCall(AOCreateMapping(comm, napp, myapp, mypetsc, aoout));
3491447629fSBarry Smith 
3509566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isapp, &myapp));
351*48a46eb9SPierre Jolivet   if (ispetsc) PetscCall(ISRestoreIndices(ispetsc, &mypetsc));
3521447629fSBarry Smith   PetscFunctionReturn(0);
3531447629fSBarry Smith }
354