xref: /petsc/src/vec/is/ao/impls/basic/aobasic.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
11447629fSBarry Smith 
21447629fSBarry Smith /*
31447629fSBarry Smith     The most basic AO application ordering routines. These store the
41447629fSBarry Smith   entire orderings 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 *app;   /* app[i] is the partner for the ith PETSc slot */
111447629fSBarry Smith   PetscInt *petsc; /* petsc[j] is the partner for the jth app slot */
121447629fSBarry Smith } AO_Basic;
131447629fSBarry Smith 
141447629fSBarry Smith /*
151447629fSBarry Smith        All processors have the same data so processor 1 prints it
161447629fSBarry Smith */
17d71ae5a4SJacob Faibussowitsch PetscErrorCode AOView_Basic(AO ao, PetscViewer viewer)
18d71ae5a4SJacob Faibussowitsch {
191447629fSBarry Smith   PetscMPIInt rank;
201447629fSBarry Smith   PetscInt    i;
211447629fSBarry Smith   AO_Basic   *aobasic = (AO_Basic *)ao->data;
221447629fSBarry Smith   PetscBool   iascii;
231447629fSBarry Smith 
241447629fSBarry Smith   PetscFunctionBegin;
259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)ao), &rank));
26dd400576SPatrick Sanan   if (rank == 0) {
279566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
281447629fSBarry Smith     if (iascii) {
299566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Number of elements in ordering %" PetscInt_FMT "\n", ao->N));
309566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "PETSc->App  App->PETSc\n"));
3148a46eb9SPierre Jolivet       for (i = 0; i < ao->N; i++) PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT "  %3" PetscInt_FMT "    %3" PetscInt_FMT "  %3" PetscInt_FMT "\n", i, aobasic->app[i], i, aobasic->petsc[i]));
321447629fSBarry Smith     }
331447629fSBarry Smith   }
349566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
35*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
361447629fSBarry Smith }
371447629fSBarry Smith 
38d71ae5a4SJacob Faibussowitsch PetscErrorCode AODestroy_Basic(AO ao)
39d71ae5a4SJacob Faibussowitsch {
401447629fSBarry Smith   AO_Basic *aobasic = (AO_Basic *)ao->data;
411447629fSBarry Smith 
421447629fSBarry Smith   PetscFunctionBegin;
439566063dSJacob Faibussowitsch   PetscCall(PetscFree2(aobasic->app, aobasic->petsc));
449566063dSJacob Faibussowitsch   PetscCall(PetscFree(aobasic));
45*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
461447629fSBarry Smith }
471447629fSBarry Smith 
48d71ae5a4SJacob Faibussowitsch PetscErrorCode AOBasicGetIndices_Private(AO ao, PetscInt **app, PetscInt **petsc)
49d71ae5a4SJacob Faibussowitsch {
501447629fSBarry Smith   AO_Basic *basic = (AO_Basic *)ao->data;
511447629fSBarry Smith 
521447629fSBarry Smith   PetscFunctionBegin;
531447629fSBarry Smith   if (app) *app = basic->app;
541447629fSBarry Smith   if (petsc) *petsc = basic->petsc;
55*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
561447629fSBarry Smith }
571447629fSBarry Smith 
58d71ae5a4SJacob Faibussowitsch PetscErrorCode AOPetscToApplication_Basic(AO ao, PetscInt n, PetscInt *ia)
59d71ae5a4SJacob Faibussowitsch {
601447629fSBarry Smith   PetscInt  i, N = ao->N;
611447629fSBarry Smith   AO_Basic *aobasic = (AO_Basic *)ao->data;
621447629fSBarry Smith 
631447629fSBarry Smith   PetscFunctionBegin;
641447629fSBarry Smith   for (i = 0; i < n; i++) {
651447629fSBarry Smith     if (ia[i] >= 0 && ia[i] < N) {
661447629fSBarry Smith       ia[i] = aobasic->app[ia[i]];
671447629fSBarry Smith     } else {
681447629fSBarry Smith       ia[i] = -1;
691447629fSBarry Smith     }
701447629fSBarry Smith   }
71*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
721447629fSBarry Smith }
731447629fSBarry Smith 
74d71ae5a4SJacob Faibussowitsch PetscErrorCode AOApplicationToPetsc_Basic(AO ao, PetscInt n, PetscInt *ia)
75d71ae5a4SJacob Faibussowitsch {
761447629fSBarry Smith   PetscInt  i, N = ao->N;
771447629fSBarry Smith   AO_Basic *aobasic = (AO_Basic *)ao->data;
781447629fSBarry Smith 
791447629fSBarry Smith   PetscFunctionBegin;
801447629fSBarry Smith   for (i = 0; i < n; i++) {
811447629fSBarry Smith     if (ia[i] >= 0 && ia[i] < N) {
821447629fSBarry Smith       ia[i] = aobasic->petsc[ia[i]];
831447629fSBarry Smith     } else {
841447629fSBarry Smith       ia[i] = -1;
851447629fSBarry Smith     }
861447629fSBarry Smith   }
87*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
881447629fSBarry Smith }
891447629fSBarry Smith 
90d71ae5a4SJacob Faibussowitsch PetscErrorCode AOPetscToApplicationPermuteInt_Basic(AO ao, PetscInt block, PetscInt *array)
91d71ae5a4SJacob Faibussowitsch {
921447629fSBarry Smith   AO_Basic *aobasic = (AO_Basic *)ao->data;
931447629fSBarry Smith   PetscInt *temp;
941447629fSBarry Smith   PetscInt  i, j;
951447629fSBarry Smith 
961447629fSBarry Smith   PetscFunctionBegin;
979566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ao->N * block, &temp));
981447629fSBarry Smith   for (i = 0; i < ao->N; i++) {
991447629fSBarry Smith     for (j = 0; j < block; j++) temp[i * block + j] = array[aobasic->petsc[i] * block + j];
1001447629fSBarry Smith   }
1019566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(array, temp, ao->N * block));
1029566063dSJacob Faibussowitsch   PetscCall(PetscFree(temp));
103*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1041447629fSBarry Smith }
1051447629fSBarry Smith 
106d71ae5a4SJacob Faibussowitsch PetscErrorCode AOApplicationToPetscPermuteInt_Basic(AO ao, PetscInt block, PetscInt *array)
107d71ae5a4SJacob Faibussowitsch {
1081447629fSBarry Smith   AO_Basic *aobasic = (AO_Basic *)ao->data;
1091447629fSBarry Smith   PetscInt *temp;
1101447629fSBarry Smith   PetscInt  i, j;
1111447629fSBarry Smith 
1121447629fSBarry Smith   PetscFunctionBegin;
1139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ao->N * block, &temp));
1141447629fSBarry Smith   for (i = 0; i < ao->N; i++) {
1151447629fSBarry Smith     for (j = 0; j < block; j++) temp[i * block + j] = array[aobasic->app[i] * block + j];
1161447629fSBarry Smith   }
1179566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(array, temp, ao->N * block));
1189566063dSJacob Faibussowitsch   PetscCall(PetscFree(temp));
119*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1201447629fSBarry Smith }
1211447629fSBarry Smith 
122d71ae5a4SJacob Faibussowitsch PetscErrorCode AOPetscToApplicationPermuteReal_Basic(AO ao, PetscInt block, PetscReal *array)
123d71ae5a4SJacob Faibussowitsch {
1241447629fSBarry Smith   AO_Basic  *aobasic = (AO_Basic *)ao->data;
1251447629fSBarry Smith   PetscReal *temp;
1261447629fSBarry Smith   PetscInt   i, j;
1271447629fSBarry Smith 
1281447629fSBarry Smith   PetscFunctionBegin;
1299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ao->N * block, &temp));
1301447629fSBarry Smith   for (i = 0; i < ao->N; i++) {
1311447629fSBarry Smith     for (j = 0; j < block; j++) temp[i * block + j] = array[aobasic->petsc[i] * block + j];
1321447629fSBarry Smith   }
1339566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(array, temp, ao->N * block));
1349566063dSJacob Faibussowitsch   PetscCall(PetscFree(temp));
135*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1361447629fSBarry Smith }
1371447629fSBarry Smith 
138d71ae5a4SJacob Faibussowitsch PetscErrorCode AOApplicationToPetscPermuteReal_Basic(AO ao, PetscInt block, PetscReal *array)
139d71ae5a4SJacob Faibussowitsch {
1401447629fSBarry Smith   AO_Basic  *aobasic = (AO_Basic *)ao->data;
1411447629fSBarry Smith   PetscReal *temp;
1421447629fSBarry Smith   PetscInt   i, j;
1431447629fSBarry Smith 
1441447629fSBarry Smith   PetscFunctionBegin;
1459566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ao->N * block, &temp));
1461447629fSBarry Smith   for (i = 0; i < ao->N; i++) {
1471447629fSBarry Smith     for (j = 0; j < block; j++) temp[i * block + j] = array[aobasic->app[i] * block + j];
1481447629fSBarry Smith   }
1499566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(array, temp, ao->N * block));
1509566063dSJacob Faibussowitsch   PetscCall(PetscFree(temp));
151*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1521447629fSBarry Smith }
1531447629fSBarry Smith 
1541447629fSBarry Smith static struct _AOOps AOOps_Basic = {
155267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(view, AOView_Basic),
156267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(destroy, AODestroy_Basic),
157267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(petsctoapplication, AOPetscToApplication_Basic),
158267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(applicationtopetsc, AOApplicationToPetsc_Basic),
159267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(petsctoapplicationpermuteint, AOPetscToApplicationPermuteInt_Basic),
160267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(applicationtopetscpermuteint, AOApplicationToPetscPermuteInt_Basic),
161267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(petsctoapplicationpermutereal, AOPetscToApplicationPermuteReal_Basic),
162267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(applicationtopetscpermutereal, AOApplicationToPetscPermuteReal_Basic),
1631447629fSBarry Smith };
1641447629fSBarry Smith 
165d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode AOCreate_Basic(AO ao)
166d71ae5a4SJacob Faibussowitsch {
1671447629fSBarry Smith   AO_Basic       *aobasic;
1681447629fSBarry Smith   PetscMPIInt     size, rank, count, *lens, *disp;
1691447629fSBarry Smith   PetscInt        napp, *allpetsc, *allapp, ip, ia, N, i, *petsc = NULL, start;
1701447629fSBarry Smith   IS              isapp = ao->isapp, ispetsc = ao->ispetsc;
1711447629fSBarry Smith   MPI_Comm        comm;
1721447629fSBarry Smith   const PetscInt *myapp, *mypetsc = NULL;
1731447629fSBarry Smith 
1741447629fSBarry Smith   PetscFunctionBegin;
1751447629fSBarry Smith   /* create special struct aobasic */
1764dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&aobasic));
1771447629fSBarry Smith   ao->data = (void *)aobasic;
1789566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(ao->ops, &AOOps_Basic, sizeof(struct _AOOps)));
1799566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)ao, AOBASIC));
1801447629fSBarry Smith 
1819566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isapp, &napp));
1829566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isapp, &myapp));
1831447629fSBarry Smith 
1849566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(napp, &count));
1851447629fSBarry Smith 
1861447629fSBarry Smith   /* transmit all lengths to all processors */
1879566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)isapp, &comm));
1889566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
1899566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1909566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(size, &lens, size, &disp));
1919566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&count, 1, MPI_INT, lens, 1, MPI_INT, comm));
1921447629fSBarry Smith   N = 0;
1931447629fSBarry Smith   for (i = 0; i < size; i++) {
1949566063dSJacob Faibussowitsch     PetscCall(PetscMPIIntCast(N, disp + i)); /* = sum(lens[j]), j< i */
1951447629fSBarry Smith     N += lens[i];
1961447629fSBarry Smith   }
1971447629fSBarry Smith   ao->N = N;
1981447629fSBarry Smith   ao->n = N;
1991447629fSBarry Smith 
2001447629fSBarry Smith   /* If mypetsc is 0 then use "natural" numbering */
2011447629fSBarry Smith   if (napp) {
2021447629fSBarry Smith     if (!ispetsc) {
2031447629fSBarry Smith       start = disp[rank];
2049566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(napp + 1, &petsc));
2051447629fSBarry Smith       for (i = 0; i < napp; i++) petsc[i] = start + i;
2061447629fSBarry Smith     } else {
2079566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(ispetsc, &mypetsc));
2081447629fSBarry Smith       petsc = (PetscInt *)mypetsc;
2091447629fSBarry Smith     }
2101447629fSBarry Smith   }
2111447629fSBarry Smith 
2121447629fSBarry Smith   /* get all indices on all processors */
2139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(N, &allpetsc, N, &allapp));
2149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgatherv(petsc, count, MPIU_INT, allpetsc, lens, disp, MPIU_INT, comm));
2159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgatherv((void *)myapp, count, MPIU_INT, allapp, lens, disp, MPIU_INT, comm));
2169566063dSJacob Faibussowitsch   PetscCall(PetscFree2(lens, disp));
2171447629fSBarry Smith 
2184e1ad211SJed Brown   if (PetscDefined(USE_DEBUG)) {
2191447629fSBarry Smith     PetscInt *sorted;
2209566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(N, &sorted));
2211447629fSBarry Smith 
2229566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(sorted, allpetsc, N));
2239566063dSJacob Faibussowitsch     PetscCall(PetscSortInt(N, sorted));
224ad540459SPierre Jolivet     for (i = 0; i < N; i++) PetscCheck(sorted[i] == i, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "PETSc ordering requires a permutation of numbers 0 to N-1\n it is missing %" PetscInt_FMT " has %" PetscInt_FMT, i, sorted[i]);
2251447629fSBarry Smith 
2269566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(sorted, allapp, N));
2279566063dSJacob Faibussowitsch     PetscCall(PetscSortInt(N, sorted));
228ad540459SPierre Jolivet     for (i = 0; i < N; i++) PetscCheck(sorted[i] == i, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Application ordering requires a permutation of numbers 0 to N-1\n it is missing %" PetscInt_FMT " has %" PetscInt_FMT, i, sorted[i]);
2291447629fSBarry Smith 
2309566063dSJacob Faibussowitsch     PetscCall(PetscFree(sorted));
2311447629fSBarry Smith   }
2321447629fSBarry Smith 
2331447629fSBarry Smith   /* generate a list of application and PETSc node numbers */
2349566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(N, &aobasic->app, N, &aobasic->petsc));
2351447629fSBarry Smith   for (i = 0; i < N; i++) {
2361447629fSBarry Smith     ip = allpetsc[i];
2371447629fSBarry Smith     ia = allapp[i];
2381447629fSBarry Smith     /* check there are no duplicates */
239c9cc58a2SBarry Smith     PetscCheck(!aobasic->app[ip], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Duplicate in PETSc ordering at position %" PetscInt_FMT ". Already mapped to %" PetscInt_FMT ", not %" PetscInt_FMT ".", i, aobasic->app[ip] - 1, ia);
2401447629fSBarry Smith     aobasic->app[ip] = ia + 1;
241c9cc58a2SBarry Smith     PetscCheck(!aobasic->petsc[ia], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Duplicate in Application ordering at position %" PetscInt_FMT ". Already mapped to %" PetscInt_FMT ", not %" PetscInt_FMT ".", i, aobasic->petsc[ia] - 1, ip);
2421447629fSBarry Smith     aobasic->petsc[ia] = ip + 1;
2431447629fSBarry Smith   }
24448a46eb9SPierre Jolivet   if (napp && !mypetsc) PetscCall(PetscFree(petsc));
2459566063dSJacob Faibussowitsch   PetscCall(PetscFree2(allpetsc, allapp));
2461447629fSBarry Smith   /* shift indices down by one */
2471447629fSBarry Smith   for (i = 0; i < N; i++) {
2481447629fSBarry Smith     aobasic->app[i]--;
2491447629fSBarry Smith     aobasic->petsc[i]--;
2501447629fSBarry Smith   }
2511447629fSBarry Smith 
2529566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isapp, &myapp));
2531447629fSBarry Smith   if (napp) {
2541447629fSBarry Smith     if (ispetsc) {
2559566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(ispetsc, &mypetsc));
2561447629fSBarry Smith     } else {
2579566063dSJacob Faibussowitsch       PetscCall(PetscFree(petsc));
2581447629fSBarry Smith     }
2591447629fSBarry Smith   }
260*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2611447629fSBarry Smith }
2621447629fSBarry Smith 
2631447629fSBarry Smith /*@C
2641447629fSBarry Smith    AOCreateBasic - Creates a basic application ordering using two integer arrays.
2651447629fSBarry Smith 
266d083f849SBarry Smith    Collective
2671447629fSBarry Smith 
2681447629fSBarry Smith    Input Parameters:
269cab54364SBarry Smith +  comm - MPI communicator that is to share `AO`
2701447629fSBarry Smith .  napp - size of integer arrays
2711447629fSBarry Smith .  myapp - integer array that defines an ordering
2721447629fSBarry Smith -  mypetsc - integer array that defines another ordering (may be NULL to
2731447629fSBarry Smith              indicate the natural ordering, that is 0,1,2,3,...)
2741447629fSBarry Smith 
2751447629fSBarry Smith    Output Parameter:
2761447629fSBarry Smith .  aoout - the new application ordering
2771447629fSBarry Smith 
2781447629fSBarry Smith    Level: beginner
2791447629fSBarry Smith 
280cab54364SBarry Smith    Note:
281cab54364SBarry Smith    The arrays myapp and mypetsc must contain the all the integers 0 to napp-1 with no duplicates; that is there cannot be any "holes"
282cab54364SBarry Smith    in the indices. Use `AOCreateMapping()` or `AOCreateMappingIS()` if you wish to have "holes" in the indices.
2831447629fSBarry Smith 
284cab54364SBarry Smith .seealso: [](sec_ao), [](sec_scatter), `AO`, `AOCreateBasicIS()`, `AODestroy()`, `AOPetscToApplication()`, `AOApplicationToPetsc()`
2851447629fSBarry Smith @*/
286d71ae5a4SJacob Faibussowitsch PetscErrorCode AOCreateBasic(MPI_Comm comm, PetscInt napp, const PetscInt myapp[], const PetscInt mypetsc[], AO *aoout)
287d71ae5a4SJacob Faibussowitsch {
2881447629fSBarry Smith   IS              isapp, ispetsc;
2891447629fSBarry Smith   const PetscInt *app = myapp, *petsc = mypetsc;
2901447629fSBarry Smith 
2911447629fSBarry Smith   PetscFunctionBegin;
2929566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm, napp, app, PETSC_USE_POINTER, &isapp));
2931447629fSBarry Smith   if (mypetsc) {
2949566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, napp, petsc, PETSC_USE_POINTER, &ispetsc));
2951447629fSBarry Smith   } else {
2961447629fSBarry Smith     ispetsc = NULL;
2971447629fSBarry Smith   }
2989566063dSJacob Faibussowitsch   PetscCall(AOCreateBasicIS(isapp, ispetsc, aoout));
2999566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isapp));
30048a46eb9SPierre Jolivet   if (mypetsc) PetscCall(ISDestroy(&ispetsc));
301*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3021447629fSBarry Smith }
3031447629fSBarry Smith 
3041447629fSBarry Smith /*@C
305cab54364SBarry Smith    AOCreateBasicIS - Creates a basic application ordering using two `IS` index sets.
3061447629fSBarry Smith 
307c3339decSBarry Smith    Collective
3081447629fSBarry Smith 
3091447629fSBarry Smith    Input Parameters:
3101447629fSBarry Smith +  isapp - index set that defines an ordering
3111447629fSBarry Smith -  ispetsc - index set that defines another ordering (may be NULL to use the
3121447629fSBarry Smith              natural ordering)
3131447629fSBarry Smith 
3141447629fSBarry Smith    Output Parameter:
3151447629fSBarry Smith .  aoout - the new application ordering
3161447629fSBarry Smith 
3171447629fSBarry Smith    Level: beginner
3181447629fSBarry Smith 
319cab54364SBarry Smith     Note:
320cab54364SBarry Smith     The index sets isapp and ispetsc must contain the all the integers 0 to napp-1 (where napp is the length of the index sets) with no duplicates;
3211447629fSBarry Smith     that is there cannot be any "holes"
3221447629fSBarry Smith 
323cab54364SBarry Smith .seealso: [](sec_ao), [](sec_scatter), `IS`, `AO`, `AOCreateBasic()`, `AODestroy()`
3241447629fSBarry Smith @*/
325d71ae5a4SJacob Faibussowitsch PetscErrorCode AOCreateBasicIS(IS isapp, IS ispetsc, AO *aoout)
326d71ae5a4SJacob Faibussowitsch {
3271447629fSBarry Smith   MPI_Comm comm;
3281447629fSBarry Smith   AO       ao;
3291447629fSBarry Smith 
3301447629fSBarry Smith   PetscFunctionBegin;
3319566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)isapp, &comm));
3329566063dSJacob Faibussowitsch   PetscCall(AOCreate(comm, &ao));
3339566063dSJacob Faibussowitsch   PetscCall(AOSetIS(ao, isapp, ispetsc));
3349566063dSJacob Faibussowitsch   PetscCall(AOSetType(ao, AOBASIC));
3359566063dSJacob Faibussowitsch   PetscCall(AOViewFromOptions(ao, NULL, "-ao_view"));
3361447629fSBarry Smith   *aoout = ao;
337*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3381447629fSBarry Smith }
339