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 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 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; 321447629fSBarry Smith PetscBool iascii; 331447629fSBarry Smith 341447629fSBarry Smith PetscFunctionBegin; 359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)ao), &rank)); 363ba16761SJacob Faibussowitsch if (rank) PetscFunctionReturn(PETSC_SUCCESS); 379566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 381447629fSBarry Smith if (iascii) { 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 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 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), 121267267bdSJacob Faibussowitsch }; 1221447629fSBarry Smith 1231447629fSBarry Smith /*@C 124cab54364SBarry Smith AOMappingHasApplicationIndex - Checks if an `AO` has a requested application index. 125cab54364SBarry Smith 126cab54364SBarry Smith Not Collective 1271447629fSBarry Smith 1281447629fSBarry Smith Input Parameters: 129cab54364SBarry Smith + ao - The `AO` 13038b5cf2dSJacob Faibussowitsch - idex - The application index 1311447629fSBarry Smith 1321447629fSBarry Smith Output Parameter: 133cab54364SBarry Smith . hasIndex - Flag is `PETSC_TRUE` if the index exists 1341447629fSBarry Smith 1351447629fSBarry Smith Level: intermediate 1361447629fSBarry Smith 13738b5cf2dSJacob Faibussowitsch Developer Notes: 138cab54364SBarry Smith The name of the function is wrong, it should be `AOHasApplicationIndex` 139cab54364SBarry Smith 140cab54364SBarry Smith .seealso: [](sec_ao), `AOMappingHasPetscIndex()`, `AOCreateMapping()`, `AO` 1411447629fSBarry Smith @*/ 142d71ae5a4SJacob Faibussowitsch PetscErrorCode AOMappingHasApplicationIndex(AO ao, PetscInt idex, PetscBool *hasIndex) 143d71ae5a4SJacob Faibussowitsch { 1441447629fSBarry Smith AO_Mapping *aomap; 1451447629fSBarry Smith PetscInt *app; 1461447629fSBarry Smith PetscInt low, high, mid; 1471447629fSBarry Smith 1481447629fSBarry Smith PetscFunctionBegin; 1491447629fSBarry Smith PetscValidHeaderSpecific(ao, AO_CLASSID, 1); 1504f572ea9SToby Isaac PetscAssertPointer(hasIndex, 3); 1511447629fSBarry Smith aomap = (AO_Mapping *)ao->data; 1521447629fSBarry Smith app = aomap->app; 1531447629fSBarry Smith /* Use bisection since the array is sorted */ 1541447629fSBarry Smith low = 0; 1551447629fSBarry Smith high = aomap->N - 1; 1561447629fSBarry Smith while (low <= high) { 1571447629fSBarry Smith mid = (low + high) / 2; 1581447629fSBarry Smith if (idex == app[mid]) break; 1591447629fSBarry Smith else if (idex < app[mid]) high = mid - 1; 1601447629fSBarry Smith else low = mid + 1; 1611447629fSBarry Smith } 1621447629fSBarry Smith if (low > high) *hasIndex = PETSC_FALSE; 1631447629fSBarry Smith else *hasIndex = PETSC_TRUE; 1643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1651447629fSBarry Smith } 1661447629fSBarry Smith 1671447629fSBarry Smith /*@ 168cab54364SBarry Smith AOMappingHasPetscIndex - checks if an `AO` has a requested petsc index. 169cab54364SBarry Smith 170cab54364SBarry Smith Not Collective 1711447629fSBarry Smith 1721447629fSBarry Smith Input Parameters: 173cab54364SBarry Smith + ao - The `AO` 17438b5cf2dSJacob Faibussowitsch - idex - The petsc index 1751447629fSBarry Smith 1761447629fSBarry Smith Output Parameter: 177cab54364SBarry Smith . hasIndex - Flag is `PETSC_TRUE` if the index exists 1781447629fSBarry Smith 1791447629fSBarry Smith Level: intermediate 1801447629fSBarry Smith 18138b5cf2dSJacob Faibussowitsch Developer Notes: 182cab54364SBarry Smith The name of the function is wrong, it should be `AOHasPetscIndex` 183cab54364SBarry Smith 184cab54364SBarry Smith .seealso: [](sec_ao), `AOMappingHasApplicationIndex()`, `AOCreateMapping()` 1851447629fSBarry Smith @*/ 186d71ae5a4SJacob Faibussowitsch PetscErrorCode AOMappingHasPetscIndex(AO ao, PetscInt idex, PetscBool *hasIndex) 187d71ae5a4SJacob Faibussowitsch { 1881447629fSBarry Smith AO_Mapping *aomap; 1891447629fSBarry Smith PetscInt *petsc; 1901447629fSBarry Smith PetscInt low, high, mid; 1911447629fSBarry Smith 1921447629fSBarry Smith PetscFunctionBegin; 1931447629fSBarry Smith PetscValidHeaderSpecific(ao, AO_CLASSID, 1); 1944f572ea9SToby Isaac PetscAssertPointer(hasIndex, 3); 1951447629fSBarry Smith aomap = (AO_Mapping *)ao->data; 1961447629fSBarry Smith petsc = aomap->petsc; 1971447629fSBarry Smith /* Use bisection since the array is sorted */ 1981447629fSBarry Smith low = 0; 1991447629fSBarry Smith high = aomap->N - 1; 2001447629fSBarry Smith while (low <= high) { 2011447629fSBarry Smith mid = (low + high) / 2; 2021447629fSBarry Smith if (idex == petsc[mid]) break; 2031447629fSBarry Smith else if (idex < petsc[mid]) high = mid - 1; 2041447629fSBarry Smith else low = mid + 1; 2051447629fSBarry Smith } 2061447629fSBarry Smith if (low > high) *hasIndex = PETSC_FALSE; 2071447629fSBarry Smith else *hasIndex = PETSC_TRUE; 2083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2091447629fSBarry Smith } 2101447629fSBarry Smith 2111447629fSBarry Smith /*@C 212cab54364SBarry Smith AOCreateMapping - Creates an application mapping using two integer arrays. 2131447629fSBarry Smith 2141447629fSBarry Smith Input Parameters: 215cab54364SBarry Smith + comm - MPI communicator that is to share the `AO` 2161447629fSBarry Smith . napp - size of integer arrays 2171447629fSBarry Smith . myapp - integer array that defines an ordering 218*dc9a610eSPierre Jolivet - mypetsc - integer array that defines another ordering (may be `NULL` to indicate the identity ordering) 2191447629fSBarry Smith 2201447629fSBarry Smith Output Parameter: 2211447629fSBarry Smith . aoout - the new application mapping 2221447629fSBarry Smith 2231447629fSBarry Smith Options Database Key: 224cab54364SBarry Smith . -ao_view - call `AOView()` at the conclusion of `AOCreateMapping()` 2251447629fSBarry Smith 2261447629fSBarry Smith Level: beginner 2271447629fSBarry Smith 228cab54364SBarry Smith Note: 229cab54364SBarry 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. 230cab54364SBarry Smith Use `AOCreateBasic()` or `AOCreateBasicIS()` if they do not have holes for better performance. 2311447629fSBarry Smith 23238b5cf2dSJacob Faibussowitsch .seealso: [](sec_ao), `AOCreateBasic()`, `AOCreateMappingIS()`, `AODestroy()` 2331447629fSBarry Smith @*/ 234d71ae5a4SJacob Faibussowitsch PetscErrorCode AOCreateMapping(MPI_Comm comm, PetscInt napp, const PetscInt myapp[], const PetscInt mypetsc[], AO *aoout) 235d71ae5a4SJacob Faibussowitsch { 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 2451447629fSBarry Smith PetscFunctionBegin; 2464f572ea9SToby Isaac PetscAssertPointer(aoout, 5); 2474c8fdceaSLisandro Dalcin *aoout = NULL; 2489566063dSJacob Faibussowitsch PetscCall(AOInitializePackage()); 2491447629fSBarry Smith 2509566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(ao, AO_CLASSID, "AO", "Application Ordering", "AO", comm, AODestroy, AOView)); 2514dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&aomap)); 252aea10558SJacob Faibussowitsch ao->ops[0] = AOps; 2531447629fSBarry Smith ao->data = (void *)aomap; 2541447629fSBarry Smith 2551447629fSBarry Smith /* transmit all lengths to all processors */ 2569566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 2579566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2589566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size, &lens, size, &disp)); 2591447629fSBarry Smith nnapp = napp; 2609566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nnapp, 1, MPI_INT, lens, 1, MPI_INT, comm)); 2611447629fSBarry Smith N = 0; 2621447629fSBarry Smith for (i = 0; i < size; i++) { 2631447629fSBarry Smith disp[i] = N; 2641447629fSBarry Smith N += lens[i]; 2651447629fSBarry Smith } 2661447629fSBarry Smith aomap->N = N; 2671447629fSBarry Smith ao->N = N; 2681447629fSBarry Smith ao->n = N; 2691447629fSBarry Smith 2701447629fSBarry Smith /* If mypetsc is 0 then use "natural" numbering */ 2711447629fSBarry Smith if (!mypetsc) { 2721447629fSBarry Smith start = disp[rank]; 2739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(napp + 1, &petsc)); 2741447629fSBarry Smith for (i = 0; i < napp; i++) petsc[i] = start + i; 2751447629fSBarry Smith } else { 2761447629fSBarry Smith petsc = (PetscInt *)mypetsc; 2771447629fSBarry Smith } 2781447629fSBarry Smith 2791447629fSBarry Smith /* get all indices on all processors */ 2809566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(N, &allapp, N, &appPerm, N, &allpetsc, N, &petscPerm)); 2819566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv((void *)myapp, napp, MPIU_INT, allapp, lens, disp, MPIU_INT, comm)); 2829566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv((void *)petsc, napp, MPIU_INT, allpetsc, lens, disp, MPIU_INT, comm)); 2839566063dSJacob Faibussowitsch PetscCall(PetscFree2(lens, disp)); 2841447629fSBarry Smith 2851447629fSBarry Smith /* generate a list of application and PETSc node numbers */ 2869566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(N, &aomap->app, N, &aomap->appPerm, N, &aomap->petsc, N, &aomap->petscPerm)); 2871447629fSBarry Smith for (i = 0; i < N; i++) { 2881447629fSBarry Smith appPerm[i] = i; 2891447629fSBarry Smith petscPerm[i] = i; 2901447629fSBarry Smith } 2919566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(N, allpetsc, petscPerm)); 2929566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(N, allapp, appPerm)); 2931447629fSBarry Smith /* Form sorted arrays of indices */ 2941447629fSBarry Smith for (i = 0; i < N; i++) { 2951447629fSBarry Smith aomap->app[i] = allapp[appPerm[i]]; 2961447629fSBarry Smith aomap->petsc[i] = allpetsc[petscPerm[i]]; 2971447629fSBarry Smith } 2981447629fSBarry Smith /* Invert petscPerm[] into aomap->petscPerm[] */ 2991447629fSBarry Smith for (i = 0; i < N; i++) aomap->petscPerm[petscPerm[i]] = i; 3001447629fSBarry Smith 3011447629fSBarry Smith /* Form map between aomap->app[] and aomap->petsc[] */ 3021447629fSBarry Smith for (i = 0; i < N; i++) aomap->appPerm[i] = aomap->petscPerm[appPerm[i]]; 3031447629fSBarry Smith 3041447629fSBarry Smith /* Invert appPerm[] into allapp[] */ 3051447629fSBarry Smith for (i = 0; i < N; i++) allapp[appPerm[i]] = i; 3061447629fSBarry Smith 3071447629fSBarry Smith /* Form map between aomap->petsc[] and aomap->app[] */ 3081447629fSBarry Smith for (i = 0; i < N; i++) aomap->petscPerm[i] = allapp[petscPerm[i]]; 3091447629fSBarry Smith 31076bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 3111447629fSBarry Smith /* Check that the permutations are complementary */ 312ad540459SPierre Jolivet for (i = 0; i < N; i++) PetscCheck(i == aomap->appPerm[aomap->petscPerm[i]], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid ordering"); 31376bd3646SJed Brown } 3141447629fSBarry Smith /* Cleanup */ 31548a46eb9SPierre Jolivet if (!mypetsc) PetscCall(PetscFree(petsc)); 3169566063dSJacob Faibussowitsch PetscCall(PetscFree4(allapp, appPerm, allpetsc, petscPerm)); 3171447629fSBarry Smith 3189566063dSJacob Faibussowitsch PetscCall(AOViewFromOptions(ao, NULL, "-ao_view")); 3191447629fSBarry Smith 3201447629fSBarry Smith *aoout = ao; 3213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3221447629fSBarry Smith } 3231447629fSBarry Smith 324487a658cSBarry Smith /*@ 325cab54364SBarry Smith AOCreateMappingIS - Creates an application mapping using two index sets. 3261447629fSBarry Smith 3271447629fSBarry Smith Input Parameters: 328e33f79d8SJacob Faibussowitsch + isapp - index set that defines an ordering 329*dc9a610eSPierre Jolivet - ispetsc - index set that defines another ordering, maybe `NULL` for identity `IS` 3301447629fSBarry Smith 3311447629fSBarry Smith Output Parameter: 3321447629fSBarry Smith . aoout - the new application ordering 3331447629fSBarry Smith 3341447629fSBarry Smith Options Database Key: 335cab54364SBarry Smith . -ao_view - call `AOView()` at the conclusion of `AOCreateMappingIS()` 3361447629fSBarry Smith 3371447629fSBarry Smith Level: beginner 3381447629fSBarry Smith 339cab54364SBarry Smith Note: 340e33f79d8SJacob 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. 341cab54364SBarry Smith Use `AOCreateBasic()` or `AOCreateBasicIS()` if they do not have holes for better performance. 3421447629fSBarry Smith 343cab54364SBarry Smith .seealso: [](sec_ao), [](sec_scatter), `AOCreateBasic()`, `AOCreateMapping()`, `AODestroy()` 3441447629fSBarry Smith @*/ 345d71ae5a4SJacob Faibussowitsch PetscErrorCode AOCreateMappingIS(IS isapp, IS ispetsc, AO *aoout) 346d71ae5a4SJacob Faibussowitsch { 3471447629fSBarry Smith MPI_Comm comm; 3481447629fSBarry Smith const PetscInt *mypetsc, *myapp; 3491447629fSBarry Smith PetscInt napp, npetsc; 3501447629fSBarry Smith 3511447629fSBarry Smith PetscFunctionBegin; 3529566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)isapp, &comm)); 3539566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isapp, &napp)); 3541447629fSBarry Smith if (ispetsc) { 3559566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(ispetsc, &npetsc)); 35608401ef6SPierre Jolivet PetscCheck(napp == npetsc, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local IS lengths must match"); 3579566063dSJacob Faibussowitsch PetscCall(ISGetIndices(ispetsc, &mypetsc)); 3581447629fSBarry Smith } else { 3591447629fSBarry Smith mypetsc = NULL; 3601447629fSBarry Smith } 3619566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isapp, &myapp)); 3621447629fSBarry Smith 3639566063dSJacob Faibussowitsch PetscCall(AOCreateMapping(comm, napp, myapp, mypetsc, aoout)); 3641447629fSBarry Smith 3659566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isapp, &myapp)); 36648a46eb9SPierre Jolivet if (ispetsc) PetscCall(ISRestoreIndices(ispetsc, &mypetsc)); 3673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3681447629fSBarry Smith } 369