xref: /petsc/src/vec/is/ao/impls/mapping/aomapping.c (revision bcda9346efad4e5ba2d553af84eb238771ba1e25)
11447629fSBarry Smith /*
21447629fSBarry Smith   These AO application ordering routines do not require that the input
31447629fSBarry Smith   be a permutation, but merely a 1-1 mapping. This implementation still
41447629fSBarry Smith   keeps the entire ordering on each processor.
51447629fSBarry Smith */
61447629fSBarry Smith 
71447629fSBarry Smith #include <../src/vec/is/ao/aoimpl.h> /*I  "petscao.h" I*/
81447629fSBarry Smith 
91447629fSBarry Smith typedef struct {
101447629fSBarry Smith   PetscInt  N;
111447629fSBarry Smith   PetscInt *app; /* app[i] is the partner for petsc[appPerm[i]] */
121447629fSBarry Smith   PetscInt *appPerm;
131447629fSBarry Smith   PetscInt *petsc; /* petsc[j] is the partner for app[petscPerm[j]] */
141447629fSBarry Smith   PetscInt *petscPerm;
151447629fSBarry Smith } AO_Mapping;
161447629fSBarry Smith 
AODestroy_Mapping(AO ao)17da8c939bSJacob Faibussowitsch static PetscErrorCode AODestroy_Mapping(AO ao)
18d71ae5a4SJacob Faibussowitsch {
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));
243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
251447629fSBarry Smith }
261447629fSBarry Smith 
AOView_Mapping(AO ao,PetscViewer viewer)27da8c939bSJacob Faibussowitsch static PetscErrorCode AOView_Mapping(AO ao, PetscViewer viewer)
28d71ae5a4SJacob Faibussowitsch {
291447629fSBarry Smith   AO_Mapping *aomap = (AO_Mapping *)ao->data;
301447629fSBarry Smith   PetscMPIInt rank;
311447629fSBarry Smith   PetscInt    i;
32*9f196a02SMartin Diehl   PetscBool   isascii;
331447629fSBarry Smith 
341447629fSBarry Smith   PetscFunctionBegin;
359566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)ao), &rank));
363ba16761SJacob Faibussowitsch   if (rank) PetscFunctionReturn(PETSC_SUCCESS);
37*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
38*9f196a02SMartin Diehl   if (isascii) {
393ba16761SJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Number of elements in ordering %" PetscInt_FMT "\n", aomap->N));
403ba16761SJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "   App.   PETSc\n"));
413ba16761SJacob Faibussowitsch     for (i = 0; i < aomap->N; i++) PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "   %" PetscInt_FMT "    %" PetscInt_FMT "\n", i, aomap->app[i], aomap->petsc[aomap->appPerm[i]]));
421447629fSBarry Smith   }
433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
441447629fSBarry Smith }
451447629fSBarry Smith 
AOPetscToApplication_Mapping(AO ao,PetscInt n,PetscInt * ia)46da8c939bSJacob Faibussowitsch static PetscErrorCode AOPetscToApplication_Mapping(AO ao, PetscInt n, PetscInt *ia)
47d71ae5a4SJacob Faibussowitsch {
481447629fSBarry Smith   AO_Mapping *aomap = (AO_Mapping *)ao->data;
491447629fSBarry Smith   PetscInt   *app   = aomap->app;
501447629fSBarry Smith   PetscInt   *petsc = aomap->petsc;
511447629fSBarry Smith   PetscInt   *perm  = aomap->petscPerm;
521447629fSBarry Smith   PetscInt    N     = aomap->N;
531447629fSBarry Smith   PetscInt    low, high, mid = 0;
541447629fSBarry Smith   PetscInt    idex;
551447629fSBarry Smith   PetscInt    i;
561447629fSBarry Smith 
571447629fSBarry Smith   /* It would be possible to use a single bisection search, which
581447629fSBarry Smith      recursively divided the indices to be converted, and searched
591447629fSBarry Smith      partitions which contained an index. This would result in
601447629fSBarry Smith      better running times if indices are clustered.
611447629fSBarry Smith   */
621447629fSBarry Smith   PetscFunctionBegin;
631447629fSBarry Smith   for (i = 0; i < n; i++) {
641447629fSBarry Smith     idex = ia[i];
651447629fSBarry Smith     if (idex < 0) continue;
661447629fSBarry Smith     /* Use bisection since the array is sorted */
671447629fSBarry Smith     low  = 0;
681447629fSBarry Smith     high = N - 1;
691447629fSBarry Smith     while (low <= high) {
701447629fSBarry Smith       mid = (low + high) / 2;
711447629fSBarry Smith       if (idex == petsc[mid]) break;
721447629fSBarry Smith       else if (idex < petsc[mid]) high = mid - 1;
731447629fSBarry Smith       else low = mid + 1;
741447629fSBarry Smith     }
751447629fSBarry Smith     if (low > high) ia[i] = -1;
761447629fSBarry Smith     else ia[i] = app[perm[mid]];
771447629fSBarry Smith   }
783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
791447629fSBarry Smith }
801447629fSBarry Smith 
AOApplicationToPetsc_Mapping(AO ao,PetscInt n,PetscInt * ia)81da8c939bSJacob Faibussowitsch static PetscErrorCode AOApplicationToPetsc_Mapping(AO ao, PetscInt n, PetscInt *ia)
82d71ae5a4SJacob Faibussowitsch {
831447629fSBarry Smith   AO_Mapping *aomap = (AO_Mapping *)ao->data;
841447629fSBarry Smith   PetscInt   *app   = aomap->app;
851447629fSBarry Smith   PetscInt   *petsc = aomap->petsc;
861447629fSBarry Smith   PetscInt   *perm  = aomap->appPerm;
871447629fSBarry Smith   PetscInt    N     = aomap->N;
881447629fSBarry Smith   PetscInt    low, high, mid = 0;
891447629fSBarry Smith   PetscInt    idex;
901447629fSBarry Smith   PetscInt    i;
911447629fSBarry Smith 
921447629fSBarry Smith   /* It would be possible to use a single bisection search, which
931447629fSBarry Smith      recursively divided the indices to be converted, and searched
941447629fSBarry Smith      partitions which contained an index. This would result in
951447629fSBarry Smith      better running times if indices are clustered.
961447629fSBarry Smith   */
971447629fSBarry Smith   PetscFunctionBegin;
981447629fSBarry Smith   for (i = 0; i < n; i++) {
991447629fSBarry Smith     idex = ia[i];
1001447629fSBarry Smith     if (idex < 0) continue;
1011447629fSBarry Smith     /* Use bisection since the array is sorted */
1021447629fSBarry Smith     low  = 0;
1031447629fSBarry Smith     high = N - 1;
1041447629fSBarry Smith     while (low <= high) {
1051447629fSBarry Smith       mid = (low + high) / 2;
1061447629fSBarry Smith       if (idex == app[mid]) break;
1071447629fSBarry Smith       else if (idex < app[mid]) high = mid - 1;
1081447629fSBarry Smith       else low = mid + 1;
1091447629fSBarry Smith     }
1101447629fSBarry Smith     if (low > high) ia[i] = -1;
1111447629fSBarry Smith     else ia[i] = petsc[perm[mid]];
1121447629fSBarry Smith   }
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1141447629fSBarry Smith }
1151447629fSBarry Smith 
116da8c939bSJacob Faibussowitsch static const struct _AOOps AOps = {
117267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(view, AOView_Mapping),
118267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(destroy, AODestroy_Mapping),
119267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(petsctoapplication, AOPetscToApplication_Mapping),
120267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(applicationtopetsc, AOApplicationToPetsc_Mapping),
1218a403008SStefano Zampini   PetscDesignatedInitializer(petsctoapplicationpermuteint, NULL),
1228a403008SStefano Zampini   PetscDesignatedInitializer(applicationtopetscpermuteint, NULL),
1238a403008SStefano Zampini   PetscDesignatedInitializer(petsctoapplicationpermutereal, NULL),
1248a403008SStefano Zampini   PetscDesignatedInitializer(applicationtopetscpermutereal, NULL),
125267267bdSJacob Faibussowitsch };
1261447629fSBarry Smith 
127cc4c1da9SBarry Smith /*@
128cab54364SBarry Smith   AOMappingHasApplicationIndex - Checks if an `AO` has a requested application index.
129cab54364SBarry Smith 
130cab54364SBarry Smith   Not Collective
1311447629fSBarry Smith 
1321447629fSBarry Smith   Input Parameters:
133cab54364SBarry Smith + ao   - The `AO`
13438b5cf2dSJacob Faibussowitsch - idex - The application index
1351447629fSBarry Smith 
1361447629fSBarry Smith   Output Parameter:
137cab54364SBarry Smith . hasIndex - Flag is `PETSC_TRUE` if the index exists
1381447629fSBarry Smith 
1391447629fSBarry Smith   Level: intermediate
1401447629fSBarry Smith 
14138b5cf2dSJacob Faibussowitsch   Developer Notes:
142cab54364SBarry Smith   The name of the function is wrong, it should be `AOHasApplicationIndex`
143cab54364SBarry Smith 
144cab54364SBarry Smith .seealso: [](sec_ao), `AOMappingHasPetscIndex()`, `AOCreateMapping()`, `AO`
1451447629fSBarry Smith @*/
AOMappingHasApplicationIndex(AO ao,PetscInt idex,PetscBool * hasIndex)146d71ae5a4SJacob Faibussowitsch PetscErrorCode AOMappingHasApplicationIndex(AO ao, PetscInt idex, PetscBool *hasIndex)
147d71ae5a4SJacob Faibussowitsch {
1481447629fSBarry Smith   AO_Mapping *aomap;
1491447629fSBarry Smith   PetscInt   *app;
1501447629fSBarry Smith   PetscInt    low, high, mid;
1511447629fSBarry Smith 
1521447629fSBarry Smith   PetscFunctionBegin;
1531447629fSBarry Smith   PetscValidHeaderSpecific(ao, AO_CLASSID, 1);
1544f572ea9SToby Isaac   PetscAssertPointer(hasIndex, 3);
1551447629fSBarry Smith   aomap = (AO_Mapping *)ao->data;
1561447629fSBarry Smith   app   = aomap->app;
1571447629fSBarry Smith   /* Use bisection since the array is sorted */
1581447629fSBarry Smith   low  = 0;
1591447629fSBarry Smith   high = aomap->N - 1;
1601447629fSBarry Smith   while (low <= high) {
1611447629fSBarry Smith     mid = (low + high) / 2;
1621447629fSBarry Smith     if (idex == app[mid]) break;
1631447629fSBarry Smith     else if (idex < app[mid]) high = mid - 1;
1641447629fSBarry Smith     else low = mid + 1;
1651447629fSBarry Smith   }
1661447629fSBarry Smith   if (low > high) *hasIndex = PETSC_FALSE;
1671447629fSBarry Smith   else *hasIndex = PETSC_TRUE;
1683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1691447629fSBarry Smith }
1701447629fSBarry Smith 
1711447629fSBarry Smith /*@
172f0b74427SPierre Jolivet   AOMappingHasPetscIndex - checks if an `AO` has a requested PETSc index.
173cab54364SBarry Smith 
174cab54364SBarry Smith   Not Collective
1751447629fSBarry Smith 
1761447629fSBarry Smith   Input Parameters:
177cab54364SBarry Smith + ao   - The `AO`
178f0b74427SPierre Jolivet - idex - The PETSc index
1791447629fSBarry Smith 
1801447629fSBarry Smith   Output Parameter:
181cab54364SBarry Smith . hasIndex - Flag is `PETSC_TRUE` if the index exists
1821447629fSBarry Smith 
1831447629fSBarry Smith   Level: intermediate
1841447629fSBarry Smith 
18538b5cf2dSJacob Faibussowitsch   Developer Notes:
186cab54364SBarry Smith   The name of the function is wrong, it should be `AOHasPetscIndex`
187cab54364SBarry Smith 
188cab54364SBarry Smith .seealso: [](sec_ao), `AOMappingHasApplicationIndex()`, `AOCreateMapping()`
1891447629fSBarry Smith @*/
AOMappingHasPetscIndex(AO ao,PetscInt idex,PetscBool * hasIndex)190d71ae5a4SJacob Faibussowitsch PetscErrorCode AOMappingHasPetscIndex(AO ao, PetscInt idex, PetscBool *hasIndex)
191d71ae5a4SJacob Faibussowitsch {
1921447629fSBarry Smith   AO_Mapping *aomap;
1931447629fSBarry Smith   PetscInt   *petsc;
1941447629fSBarry Smith   PetscInt    low, high, mid;
1951447629fSBarry Smith 
1961447629fSBarry Smith   PetscFunctionBegin;
1971447629fSBarry Smith   PetscValidHeaderSpecific(ao, AO_CLASSID, 1);
1984f572ea9SToby Isaac   PetscAssertPointer(hasIndex, 3);
1991447629fSBarry Smith   aomap = (AO_Mapping *)ao->data;
2001447629fSBarry Smith   petsc = aomap->petsc;
2011447629fSBarry Smith   /* Use bisection since the array is sorted */
2021447629fSBarry Smith   low  = 0;
2031447629fSBarry Smith   high = aomap->N - 1;
2041447629fSBarry Smith   while (low <= high) {
2051447629fSBarry Smith     mid = (low + high) / 2;
2061447629fSBarry Smith     if (idex == petsc[mid]) break;
2071447629fSBarry Smith     else if (idex < petsc[mid]) high = mid - 1;
2081447629fSBarry Smith     else low = mid + 1;
2091447629fSBarry Smith   }
2101447629fSBarry Smith   if (low > high) *hasIndex = PETSC_FALSE;
2111447629fSBarry Smith   else *hasIndex = PETSC_TRUE;
2123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2131447629fSBarry Smith }
2141447629fSBarry Smith 
215cc4c1da9SBarry Smith /*@
216cab54364SBarry Smith   AOCreateMapping - Creates an application mapping using two integer arrays.
2171447629fSBarry Smith 
2181447629fSBarry Smith   Input Parameters:
219cab54364SBarry Smith + comm    - MPI communicator that is to share the `AO`
2201447629fSBarry Smith . napp    - size of integer arrays
2211447629fSBarry Smith . myapp   - integer array that defines an ordering
222dc9a610eSPierre Jolivet - mypetsc - integer array that defines another ordering (may be `NULL` to indicate the identity ordering)
2231447629fSBarry Smith 
2241447629fSBarry Smith   Output Parameter:
2251447629fSBarry Smith . aoout - the new application mapping
2261447629fSBarry Smith 
2271447629fSBarry Smith   Options Database Key:
228cab54364SBarry Smith . -ao_view - call `AOView()` at the conclusion of `AOCreateMapping()`
2291447629fSBarry Smith 
2301447629fSBarry Smith   Level: beginner
2311447629fSBarry Smith 
232cab54364SBarry Smith   Note:
233b6971eaeSBarry Smith   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.
234cab54364SBarry Smith   Use `AOCreateBasic()` or `AOCreateBasicIS()` if they do not have holes for better performance.
2351447629fSBarry Smith 
23638b5cf2dSJacob Faibussowitsch .seealso: [](sec_ao), `AOCreateBasic()`, `AOCreateMappingIS()`, `AODestroy()`
2371447629fSBarry Smith @*/
AOCreateMapping(MPI_Comm comm,PetscInt napp,const PetscInt myapp[],const PetscInt mypetsc[],AO * aoout)238d71ae5a4SJacob Faibussowitsch PetscErrorCode AOCreateMapping(MPI_Comm comm, PetscInt napp, const PetscInt myapp[], const PetscInt mypetsc[], AO *aoout)
239d71ae5a4SJacob Faibussowitsch {
2401447629fSBarry Smith   AO          ao;
2411447629fSBarry Smith   AO_Mapping *aomap;
2421447629fSBarry Smith   PetscInt   *allpetsc, *allapp;
2431447629fSBarry Smith   PetscInt   *petscPerm, *appPerm;
2441447629fSBarry Smith   PetscInt   *petsc;
2451447629fSBarry Smith   PetscMPIInt size, rank, *lens, *disp, nnapp;
2466497c311SBarry Smith   PetscCount  N;
2476497c311SBarry Smith   PetscInt    start;
2481447629fSBarry Smith 
2491447629fSBarry Smith   PetscFunctionBegin;
2504f572ea9SToby Isaac   PetscAssertPointer(aoout, 5);
2514c8fdceaSLisandro Dalcin   *aoout = NULL;
2529566063dSJacob Faibussowitsch   PetscCall(AOInitializePackage());
2531447629fSBarry Smith 
2549566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(ao, AO_CLASSID, "AO", "Application Ordering", "AO", comm, AODestroy, AOView));
2554dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&aomap));
256aea10558SJacob Faibussowitsch   ao->ops[0] = AOps;
2571447629fSBarry Smith   ao->data   = (void *)aomap;
2581447629fSBarry Smith 
2591447629fSBarry Smith   /* transmit all lengths to all processors */
2609566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
2619566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2629566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(size, &lens, size, &disp));
2636497c311SBarry Smith   PetscCall(PetscMPIIntCast(napp, &nnapp));
2649566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&nnapp, 1, MPI_INT, lens, 1, MPI_INT, comm));
2651447629fSBarry Smith   N = 0;
2666497c311SBarry Smith   for (PetscMPIInt i = 0; i < size; i++) {
2676497c311SBarry Smith     PetscCall(PetscMPIIntCast(N, &disp[i]));
2681447629fSBarry Smith     N += lens[i];
2691447629fSBarry Smith   }
2706497c311SBarry Smith   PetscCall(PetscIntCast(N, &aomap->N));
271835f2295SStefano Zampini   ao->N = aomap->N;
272835f2295SStefano Zampini   ao->n = aomap->N;
2731447629fSBarry Smith 
2741447629fSBarry Smith   /* If mypetsc is 0 then use "natural" numbering */
2751447629fSBarry Smith   if (!mypetsc) {
2761447629fSBarry Smith     start = disp[rank];
2779566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(napp + 1, &petsc));
2786497c311SBarry Smith     for (PetscInt i = 0; i < napp; i++) petsc[i] = start + i;
2791447629fSBarry Smith   } else {
2801447629fSBarry Smith     petsc = (PetscInt *)mypetsc;
2811447629fSBarry Smith   }
2821447629fSBarry Smith 
2831447629fSBarry Smith   /* get all indices on all processors */
2849566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(N, &allapp, N, &appPerm, N, &allpetsc, N, &petscPerm));
2856497c311SBarry Smith   PetscCallMPI(MPI_Allgatherv((void *)myapp, nnapp, MPIU_INT, allapp, lens, disp, MPIU_INT, comm));
2866497c311SBarry Smith   PetscCallMPI(MPI_Allgatherv((void *)petsc, nnapp, MPIU_INT, allpetsc, lens, disp, MPIU_INT, comm));
2879566063dSJacob Faibussowitsch   PetscCall(PetscFree2(lens, disp));
2881447629fSBarry Smith 
2891447629fSBarry Smith   /* generate a list of application and PETSc node numbers */
2909566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(N, &aomap->app, N, &aomap->appPerm, N, &aomap->petsc, N, &aomap->petscPerm));
2916497c311SBarry Smith   for (PetscInt i = 0; i < N; i++) {
2921447629fSBarry Smith     appPerm[i]   = i;
2931447629fSBarry Smith     petscPerm[i] = i;
2941447629fSBarry Smith   }
295835f2295SStefano Zampini   PetscCall(PetscSortIntWithPermutation(aomap->N, allpetsc, petscPerm));
296835f2295SStefano Zampini   PetscCall(PetscSortIntWithPermutation(aomap->N, allapp, appPerm));
2971447629fSBarry Smith   /* Form sorted arrays of indices */
2986497c311SBarry Smith   for (PetscInt i = 0; i < N; i++) {
2991447629fSBarry Smith     aomap->app[i]   = allapp[appPerm[i]];
3001447629fSBarry Smith     aomap->petsc[i] = allpetsc[petscPerm[i]];
3011447629fSBarry Smith   }
3021447629fSBarry Smith   /* Invert petscPerm[] into aomap->petscPerm[] */
3036497c311SBarry Smith   for (PetscInt i = 0; i < N; i++) aomap->petscPerm[petscPerm[i]] = i;
3041447629fSBarry Smith 
3051447629fSBarry Smith   /* Form map between aomap->app[] and aomap->petsc[] */
3066497c311SBarry Smith   for (PetscInt i = 0; i < N; i++) aomap->appPerm[i] = aomap->petscPerm[appPerm[i]];
3071447629fSBarry Smith 
3081447629fSBarry Smith   /* Invert appPerm[] into allapp[] */
3096497c311SBarry Smith   for (PetscInt i = 0; i < N; i++) allapp[appPerm[i]] = i;
3101447629fSBarry Smith 
3111447629fSBarry Smith   /* Form map between aomap->petsc[] and aomap->app[] */
3126497c311SBarry Smith   for (PetscInt i = 0; i < N; i++) aomap->petscPerm[i] = allapp[petscPerm[i]];
3131447629fSBarry Smith 
31476bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
3151447629fSBarry Smith     /* Check that the permutations are complementary */
3166497c311SBarry Smith     for (PetscInt i = 0; i < N; i++) PetscCheck(i == aomap->appPerm[aomap->petscPerm[i]], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid ordering");
31776bd3646SJed Brown   }
3181447629fSBarry Smith   /* Cleanup */
31948a46eb9SPierre Jolivet   if (!mypetsc) PetscCall(PetscFree(petsc));
3209566063dSJacob Faibussowitsch   PetscCall(PetscFree4(allapp, appPerm, allpetsc, petscPerm));
3219566063dSJacob Faibussowitsch   PetscCall(AOViewFromOptions(ao, NULL, "-ao_view"));
3221447629fSBarry Smith   *aoout = ao;
3233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3241447629fSBarry Smith }
3251447629fSBarry Smith 
3265d83a8b1SBarry Smith /*@
327cab54364SBarry Smith   AOCreateMappingIS - Creates an application mapping using two index sets.
3281447629fSBarry Smith 
3291447629fSBarry Smith   Input Parameters:
330e33f79d8SJacob Faibussowitsch + isapp   - index set that defines an ordering
331dc9a610eSPierre Jolivet - ispetsc - index set that defines another ordering, maybe `NULL` for identity `IS`
3321447629fSBarry Smith 
3331447629fSBarry Smith   Output Parameter:
3341447629fSBarry Smith . aoout - the new application ordering
3351447629fSBarry Smith 
3361447629fSBarry Smith   Options Database Key:
337cab54364SBarry Smith . -ao_view - call `AOView()` at the conclusion of `AOCreateMappingIS()`
3381447629fSBarry Smith 
3391447629fSBarry Smith   Level: beginner
3401447629fSBarry Smith 
341cab54364SBarry Smith   Note:
342e33f79d8SJacob Faibussowitsch   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.
343cab54364SBarry Smith   Use `AOCreateBasic()` or `AOCreateBasicIS()` if they do not have holes for better performance.
3441447629fSBarry Smith 
345cab54364SBarry Smith .seealso: [](sec_ao), [](sec_scatter), `AOCreateBasic()`, `AOCreateMapping()`, `AODestroy()`
3461447629fSBarry Smith @*/
AOCreateMappingIS(IS isapp,IS ispetsc,AO * aoout)347d71ae5a4SJacob Faibussowitsch PetscErrorCode AOCreateMappingIS(IS isapp, IS ispetsc, AO *aoout)
348d71ae5a4SJacob Faibussowitsch {
3491447629fSBarry Smith   MPI_Comm        comm;
3501447629fSBarry Smith   const PetscInt *mypetsc, *myapp;
3511447629fSBarry Smith   PetscInt        napp, npetsc;
3521447629fSBarry Smith 
3531447629fSBarry Smith   PetscFunctionBegin;
3549566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)isapp, &comm));
3559566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isapp, &napp));
3561447629fSBarry Smith   if (ispetsc) {
3579566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(ispetsc, &npetsc));
35808401ef6SPierre Jolivet     PetscCheck(napp == npetsc, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local IS lengths must match");
3599566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(ispetsc, &mypetsc));
3601447629fSBarry Smith   } else {
3611447629fSBarry Smith     mypetsc = NULL;
3621447629fSBarry Smith   }
3639566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isapp, &myapp));
3641447629fSBarry Smith 
3659566063dSJacob Faibussowitsch   PetscCall(AOCreateMapping(comm, napp, myapp, mypetsc, aoout));
3661447629fSBarry Smith 
3679566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isapp, &myapp));
36848a46eb9SPierre Jolivet   if (ispetsc) PetscCall(ISRestoreIndices(ispetsc, &mypetsc));
3693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3701447629fSBarry Smith }
371