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 18d71ae5a4SJacob Faibussowitsch PetscErrorCode AODestroy_Mapping(AO ao) 19d71ae5a4SJacob 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)); 253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 261447629fSBarry Smith } 271447629fSBarry Smith 28d71ae5a4SJacob Faibussowitsch PetscErrorCode AOView_Mapping(AO ao, PetscViewer viewer) 29d71ae5a4SJacob 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)); 373ba16761SJacob Faibussowitsch if (rank) PetscFunctionReturn(PETSC_SUCCESS); 389566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 391447629fSBarry Smith if (iascii) { 403ba16761SJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Number of elements in ordering %" PetscInt_FMT "\n", aomap->N)); 413ba16761SJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " App. PETSc\n")); 423ba16761SJacob 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]])); 431447629fSBarry Smith } 443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 451447629fSBarry Smith } 461447629fSBarry Smith 47d71ae5a4SJacob Faibussowitsch PetscErrorCode AOPetscToApplication_Mapping(AO ao, PetscInt n, PetscInt *ia) 48d71ae5a4SJacob 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 } 793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 801447629fSBarry Smith } 811447629fSBarry Smith 82d71ae5a4SJacob Faibussowitsch PetscErrorCode AOApplicationToPetsc_Mapping(AO ao, PetscInt n, PetscInt *ia) 83d71ae5a4SJacob 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 } 1143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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 125cab54364SBarry Smith AOMappingHasApplicationIndex - Checks if an `AO` has a requested application index. 126cab54364SBarry Smith 127cab54364SBarry Smith Not Collective 1281447629fSBarry Smith 1291447629fSBarry Smith Input Parameters: 130cab54364SBarry Smith + ao - The `AO` 13138b5cf2dSJacob Faibussowitsch - idex - The application index 1321447629fSBarry Smith 1331447629fSBarry Smith Output Parameter: 134cab54364SBarry Smith . hasIndex - Flag is `PETSC_TRUE` if the index exists 1351447629fSBarry Smith 1361447629fSBarry Smith Level: intermediate 1371447629fSBarry Smith 13838b5cf2dSJacob Faibussowitsch Developer Notes: 139cab54364SBarry Smith The name of the function is wrong, it should be `AOHasApplicationIndex` 140cab54364SBarry Smith 141cab54364SBarry Smith .seealso: [](sec_ao), `AOMappingHasPetscIndex()`, `AOCreateMapping()`, `AO` 1421447629fSBarry Smith @*/ 143d71ae5a4SJacob Faibussowitsch PetscErrorCode AOMappingHasApplicationIndex(AO ao, PetscInt idex, PetscBool *hasIndex) 144d71ae5a4SJacob Faibussowitsch { 1451447629fSBarry Smith AO_Mapping *aomap; 1461447629fSBarry Smith PetscInt *app; 1471447629fSBarry Smith PetscInt low, high, mid; 1481447629fSBarry Smith 1491447629fSBarry Smith PetscFunctionBegin; 1501447629fSBarry Smith PetscValidHeaderSpecific(ao, AO_CLASSID, 1); 151dadcf809SJacob Faibussowitsch PetscValidBoolPointer(hasIndex, 3); 1521447629fSBarry Smith aomap = (AO_Mapping *)ao->data; 1531447629fSBarry Smith app = aomap->app; 1541447629fSBarry Smith /* Use bisection since the array is sorted */ 1551447629fSBarry Smith low = 0; 1561447629fSBarry Smith high = aomap->N - 1; 1571447629fSBarry Smith while (low <= high) { 1581447629fSBarry Smith mid = (low + high) / 2; 1591447629fSBarry Smith if (idex == app[mid]) break; 1601447629fSBarry Smith else if (idex < app[mid]) high = mid - 1; 1611447629fSBarry Smith else low = mid + 1; 1621447629fSBarry Smith } 1631447629fSBarry Smith if (low > high) *hasIndex = PETSC_FALSE; 1641447629fSBarry Smith else *hasIndex = PETSC_TRUE; 1653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1661447629fSBarry Smith } 1671447629fSBarry Smith 1681447629fSBarry Smith /*@ 169cab54364SBarry Smith AOMappingHasPetscIndex - checks if an `AO` has a requested petsc index. 170cab54364SBarry Smith 171cab54364SBarry Smith Not Collective 1721447629fSBarry Smith 1731447629fSBarry Smith Input Parameters: 174cab54364SBarry Smith + ao - The `AO` 17538b5cf2dSJacob Faibussowitsch - idex - The petsc index 1761447629fSBarry Smith 1771447629fSBarry Smith Output Parameter: 178cab54364SBarry Smith . hasIndex - Flag is `PETSC_TRUE` if the index exists 1791447629fSBarry Smith 1801447629fSBarry Smith Level: intermediate 1811447629fSBarry Smith 18238b5cf2dSJacob Faibussowitsch Developer Notes: 183cab54364SBarry Smith The name of the function is wrong, it should be `AOHasPetscIndex` 184cab54364SBarry Smith 185cab54364SBarry Smith .seealso: [](sec_ao), `AOMappingHasApplicationIndex()`, `AOCreateMapping()` 1861447629fSBarry Smith @*/ 187d71ae5a4SJacob Faibussowitsch PetscErrorCode AOMappingHasPetscIndex(AO ao, PetscInt idex, PetscBool *hasIndex) 188d71ae5a4SJacob Faibussowitsch { 1891447629fSBarry Smith AO_Mapping *aomap; 1901447629fSBarry Smith PetscInt *petsc; 1911447629fSBarry Smith PetscInt low, high, mid; 1921447629fSBarry Smith 1931447629fSBarry Smith PetscFunctionBegin; 1941447629fSBarry Smith PetscValidHeaderSpecific(ao, AO_CLASSID, 1); 195dadcf809SJacob Faibussowitsch PetscValidBoolPointer(hasIndex, 3); 1961447629fSBarry Smith aomap = (AO_Mapping *)ao->data; 1971447629fSBarry Smith petsc = aomap->petsc; 1981447629fSBarry Smith /* Use bisection since the array is sorted */ 1991447629fSBarry Smith low = 0; 2001447629fSBarry Smith high = aomap->N - 1; 2011447629fSBarry Smith while (low <= high) { 2021447629fSBarry Smith mid = (low + high) / 2; 2031447629fSBarry Smith if (idex == petsc[mid]) break; 2041447629fSBarry Smith else if (idex < petsc[mid]) high = mid - 1; 2051447629fSBarry Smith else low = mid + 1; 2061447629fSBarry Smith } 2071447629fSBarry Smith if (low > high) *hasIndex = PETSC_FALSE; 2081447629fSBarry Smith else *hasIndex = PETSC_TRUE; 2093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2101447629fSBarry Smith } 2111447629fSBarry Smith 2121447629fSBarry Smith /*@C 213cab54364SBarry Smith AOCreateMapping - Creates an application mapping using two integer arrays. 2141447629fSBarry Smith 2151447629fSBarry Smith Input Parameters: 216cab54364SBarry Smith + comm - MPI communicator that is to share the `AO` 2171447629fSBarry Smith . napp - size of integer arrays 2181447629fSBarry Smith . myapp - integer array that defines an ordering 2191447629fSBarry Smith - mypetsc - integer array that defines another ordering (may be NULL to indicate the identity ordering) 2201447629fSBarry Smith 2211447629fSBarry Smith Output Parameter: 2221447629fSBarry Smith . aoout - the new application mapping 2231447629fSBarry Smith 2241447629fSBarry Smith Options Database Key: 225cab54364SBarry Smith . -ao_view - call `AOView()` at the conclusion of `AOCreateMapping()` 2261447629fSBarry Smith 2271447629fSBarry Smith Level: beginner 2281447629fSBarry Smith 229cab54364SBarry Smith Note: 230cab54364SBarry 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. 231cab54364SBarry Smith Use `AOCreateBasic()` or `AOCreateBasicIS()` if they do not have holes for better performance. 2321447629fSBarry Smith 23338b5cf2dSJacob Faibussowitsch .seealso: [](sec_ao), `AOCreateBasic()`, `AOCreateMappingIS()`, `AODestroy()` 2341447629fSBarry Smith @*/ 235d71ae5a4SJacob Faibussowitsch PetscErrorCode AOCreateMapping(MPI_Comm comm, PetscInt napp, const PetscInt myapp[], const PetscInt mypetsc[], AO *aoout) 236d71ae5a4SJacob Faibussowitsch { 2371447629fSBarry Smith AO ao; 2381447629fSBarry Smith AO_Mapping *aomap; 2391447629fSBarry Smith PetscInt *allpetsc, *allapp; 2401447629fSBarry Smith PetscInt *petscPerm, *appPerm; 2411447629fSBarry Smith PetscInt *petsc; 2421447629fSBarry Smith PetscMPIInt size, rank, *lens, *disp, nnapp; 2431447629fSBarry Smith PetscInt N, start; 2441447629fSBarry Smith PetscInt i; 2451447629fSBarry Smith 2461447629fSBarry Smith PetscFunctionBegin; 2471447629fSBarry Smith PetscValidPointer(aoout, 5); 2484c8fdceaSLisandro Dalcin *aoout = NULL; 2499566063dSJacob Faibussowitsch PetscCall(AOInitializePackage()); 2501447629fSBarry Smith 2519566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(ao, AO_CLASSID, "AO", "Application Ordering", "AO", comm, AODestroy, AOView)); 2524dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&aomap)); 253aea10558SJacob Faibussowitsch ao->ops[0] = AOps; 2541447629fSBarry Smith ao->data = (void *)aomap; 2551447629fSBarry Smith 2561447629fSBarry Smith /* transmit all lengths to all processors */ 2579566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 2589566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2599566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size, &lens, size, &disp)); 2601447629fSBarry Smith nnapp = napp; 2619566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nnapp, 1, MPI_INT, lens, 1, MPI_INT, comm)); 2621447629fSBarry Smith N = 0; 2631447629fSBarry Smith for (i = 0; i < size; i++) { 2641447629fSBarry Smith disp[i] = N; 2651447629fSBarry Smith N += lens[i]; 2661447629fSBarry Smith } 2671447629fSBarry Smith aomap->N = N; 2681447629fSBarry Smith ao->N = N; 2691447629fSBarry Smith ao->n = N; 2701447629fSBarry Smith 2711447629fSBarry Smith /* If mypetsc is 0 then use "natural" numbering */ 2721447629fSBarry Smith if (!mypetsc) { 2731447629fSBarry Smith start = disp[rank]; 2749566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(napp + 1, &petsc)); 2751447629fSBarry Smith for (i = 0; i < napp; i++) petsc[i] = start + i; 2761447629fSBarry Smith } else { 2771447629fSBarry Smith petsc = (PetscInt *)mypetsc; 2781447629fSBarry Smith } 2791447629fSBarry Smith 2801447629fSBarry Smith /* get all indices on all processors */ 2819566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(N, &allapp, N, &appPerm, N, &allpetsc, N, &petscPerm)); 2829566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv((void *)myapp, napp, MPIU_INT, allapp, lens, disp, MPIU_INT, comm)); 2839566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv((void *)petsc, napp, MPIU_INT, allpetsc, lens, disp, MPIU_INT, comm)); 2849566063dSJacob Faibussowitsch PetscCall(PetscFree2(lens, disp)); 2851447629fSBarry Smith 2861447629fSBarry Smith /* generate a list of application and PETSc node numbers */ 2879566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(N, &aomap->app, N, &aomap->appPerm, N, &aomap->petsc, N, &aomap->petscPerm)); 2881447629fSBarry Smith for (i = 0; i < N; i++) { 2891447629fSBarry Smith appPerm[i] = i; 2901447629fSBarry Smith petscPerm[i] = i; 2911447629fSBarry Smith } 2929566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(N, allpetsc, petscPerm)); 2939566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(N, allapp, appPerm)); 2941447629fSBarry Smith /* Form sorted arrays of indices */ 2951447629fSBarry Smith for (i = 0; i < N; i++) { 2961447629fSBarry Smith aomap->app[i] = allapp[appPerm[i]]; 2971447629fSBarry Smith aomap->petsc[i] = allpetsc[petscPerm[i]]; 2981447629fSBarry Smith } 2991447629fSBarry Smith /* Invert petscPerm[] into aomap->petscPerm[] */ 3001447629fSBarry Smith for (i = 0; i < N; i++) aomap->petscPerm[petscPerm[i]] = i; 3011447629fSBarry Smith 3021447629fSBarry Smith /* Form map between aomap->app[] and aomap->petsc[] */ 3031447629fSBarry Smith for (i = 0; i < N; i++) aomap->appPerm[i] = aomap->petscPerm[appPerm[i]]; 3041447629fSBarry Smith 3051447629fSBarry Smith /* Invert appPerm[] into allapp[] */ 3061447629fSBarry Smith for (i = 0; i < N; i++) allapp[appPerm[i]] = i; 3071447629fSBarry Smith 3081447629fSBarry Smith /* Form map between aomap->petsc[] and aomap->app[] */ 3091447629fSBarry Smith for (i = 0; i < N; i++) aomap->petscPerm[i] = allapp[petscPerm[i]]; 3101447629fSBarry Smith 31176bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 3121447629fSBarry Smith /* Check that the permutations are complementary */ 313ad540459SPierre Jolivet for (i = 0; i < N; i++) PetscCheck(i == aomap->appPerm[aomap->petscPerm[i]], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid ordering"); 31476bd3646SJed Brown } 3151447629fSBarry Smith /* Cleanup */ 31648a46eb9SPierre Jolivet if (!mypetsc) PetscCall(PetscFree(petsc)); 3179566063dSJacob Faibussowitsch PetscCall(PetscFree4(allapp, appPerm, allpetsc, petscPerm)); 3181447629fSBarry Smith 3199566063dSJacob Faibussowitsch PetscCall(AOViewFromOptions(ao, NULL, "-ao_view")); 3201447629fSBarry Smith 3211447629fSBarry Smith *aoout = ao; 3223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3231447629fSBarry Smith } 3241447629fSBarry Smith 325487a658cSBarry Smith /*@ 326cab54364SBarry Smith AOCreateMappingIS - Creates an application mapping using two index sets. 3271447629fSBarry Smith 3281447629fSBarry Smith Input Parameters: 329*e33f79d8SJacob Faibussowitsch + isapp - index set that defines an ordering 330cab54364SBarry Smith - ispetsc - index set that defines another ordering, maybe NULL for identity `IS` 3311447629fSBarry Smith 3321447629fSBarry Smith Output Parameter: 3331447629fSBarry Smith . aoout - the new application ordering 3341447629fSBarry Smith 3351447629fSBarry Smith Options Database Key: 336cab54364SBarry Smith . -ao_view - call `AOView()` at the conclusion of `AOCreateMappingIS()` 3371447629fSBarry Smith 3381447629fSBarry Smith Level: beginner 3391447629fSBarry Smith 340cab54364SBarry Smith Note: 341*e33f79d8SJacob 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. 342cab54364SBarry Smith Use `AOCreateBasic()` or `AOCreateBasicIS()` if they do not have holes for better performance. 3431447629fSBarry Smith 344cab54364SBarry Smith .seealso: [](sec_ao), [](sec_scatter), `AOCreateBasic()`, `AOCreateMapping()`, `AODestroy()` 3451447629fSBarry Smith @*/ 346d71ae5a4SJacob Faibussowitsch PetscErrorCode AOCreateMappingIS(IS isapp, IS ispetsc, AO *aoout) 347d71ae5a4SJacob Faibussowitsch { 3481447629fSBarry Smith MPI_Comm comm; 3491447629fSBarry Smith const PetscInt *mypetsc, *myapp; 3501447629fSBarry Smith PetscInt napp, npetsc; 3511447629fSBarry Smith 3521447629fSBarry Smith PetscFunctionBegin; 3539566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)isapp, &comm)); 3549566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isapp, &napp)); 3551447629fSBarry Smith if (ispetsc) { 3569566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(ispetsc, &npetsc)); 35708401ef6SPierre Jolivet PetscCheck(napp == npetsc, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local IS lengths must match"); 3589566063dSJacob Faibussowitsch PetscCall(ISGetIndices(ispetsc, &mypetsc)); 3591447629fSBarry Smith } else { 3601447629fSBarry Smith mypetsc = NULL; 3611447629fSBarry Smith } 3629566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isapp, &myapp)); 3631447629fSBarry Smith 3649566063dSJacob Faibussowitsch PetscCall(AOCreateMapping(comm, napp, myapp, mypetsc, aoout)); 3651447629fSBarry Smith 3669566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isapp, &myapp)); 36748a46eb9SPierre Jolivet if (ispetsc) PetscCall(ISRestoreIndices(ispetsc, &mypetsc)); 3683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3691447629fSBarry Smith } 370