xref: /petsc/src/vec/is/ao/impls/mapping/aomapping.c (revision dcca6d9d80ebd869fe6029bd05a3aa9faafef49e)
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 
181447629fSBarry Smith #undef __FUNCT__
191447629fSBarry Smith #define __FUNCT__ "AODestroy_Mapping"
201447629fSBarry Smith PetscErrorCode AODestroy_Mapping(AO ao)
211447629fSBarry Smith {
221447629fSBarry Smith   AO_Mapping     *aomap = (AO_Mapping*) ao->data;
231447629fSBarry Smith   PetscErrorCode ierr;
241447629fSBarry Smith 
251447629fSBarry Smith   PetscFunctionBegin;
261447629fSBarry Smith   ierr = PetscFree4(aomap->app,aomap->appPerm,aomap->petsc,aomap->petscPerm);CHKERRQ(ierr);
271447629fSBarry Smith   ierr = PetscFree(aomap);CHKERRQ(ierr);
281447629fSBarry Smith   PetscFunctionReturn(0);
291447629fSBarry Smith }
301447629fSBarry Smith 
311447629fSBarry Smith #undef __FUNCT__
321447629fSBarry Smith #define __FUNCT__ "AOView_Mapping"
331447629fSBarry Smith PetscErrorCode AOView_Mapping(AO ao, PetscViewer viewer)
341447629fSBarry Smith {
351447629fSBarry Smith   AO_Mapping     *aomap = (AO_Mapping*) ao->data;
361447629fSBarry Smith   PetscMPIInt    rank;
371447629fSBarry Smith   PetscInt       i;
381447629fSBarry Smith   PetscBool      iascii;
391447629fSBarry Smith   PetscErrorCode ierr;
401447629fSBarry Smith 
411447629fSBarry Smith   PetscFunctionBegin;
421447629fSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)ao), &rank);CHKERRQ(ierr);
431447629fSBarry Smith   if (rank) PetscFunctionReturn(0);
441447629fSBarry Smith 
451447629fSBarry Smith   if (!viewer) viewer = PETSC_VIEWER_STDOUT_SELF;
461447629fSBarry Smith 
471447629fSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
481447629fSBarry Smith   if (iascii) {
491447629fSBarry Smith     PetscViewerASCIIPrintf(viewer, "Number of elements in ordering %D\n", aomap->N);
501447629fSBarry Smith     PetscViewerASCIIPrintf(viewer, "   App.   PETSc\n");
511447629fSBarry Smith     for (i = 0; i < aomap->N; i++) {
521447629fSBarry Smith       PetscViewerASCIIPrintf(viewer, "%D   %D    %D\n", i, aomap->app[i], aomap->petsc[aomap->appPerm[i]]);
531447629fSBarry Smith     }
541447629fSBarry Smith   }
551447629fSBarry Smith   PetscFunctionReturn(0);
561447629fSBarry Smith }
571447629fSBarry Smith 
581447629fSBarry Smith #undef __FUNCT__
591447629fSBarry Smith #define __FUNCT__ "AOPetscToApplication_Mapping"
601447629fSBarry Smith PetscErrorCode AOPetscToApplication_Mapping(AO ao, PetscInt n, PetscInt *ia)
611447629fSBarry Smith {
621447629fSBarry Smith   AO_Mapping *aomap = (AO_Mapping*) ao->data;
631447629fSBarry Smith   PetscInt   *app   = aomap->app;
641447629fSBarry Smith   PetscInt   *petsc = aomap->petsc;
651447629fSBarry Smith   PetscInt   *perm  = aomap->petscPerm;
661447629fSBarry Smith   PetscInt   N      = aomap->N;
671447629fSBarry Smith   PetscInt   low, high, mid=0;
681447629fSBarry Smith   PetscInt   idex;
691447629fSBarry Smith   PetscInt   i;
701447629fSBarry Smith 
711447629fSBarry Smith   /* It would be possible to use a single bisection search, which
721447629fSBarry Smith      recursively divided the indices to be converted, and searched
731447629fSBarry Smith      partitions which contained an index. This would result in
741447629fSBarry Smith      better running times if indices are clustered.
751447629fSBarry Smith   */
761447629fSBarry Smith   PetscFunctionBegin;
771447629fSBarry Smith   for (i = 0; i < n; i++) {
781447629fSBarry Smith     idex = ia[i];
791447629fSBarry Smith     if (idex < 0) continue;
801447629fSBarry Smith     /* Use bisection since the array is sorted */
811447629fSBarry Smith     low  = 0;
821447629fSBarry Smith     high = N - 1;
831447629fSBarry Smith     while (low <= high) {
841447629fSBarry Smith       mid = (low + high)/2;
851447629fSBarry Smith       if (idex == petsc[mid]) break;
861447629fSBarry Smith       else if (idex < petsc[mid]) high = mid - 1;
871447629fSBarry Smith       else low = mid + 1;
881447629fSBarry Smith     }
891447629fSBarry Smith     if (low > high) ia[i] = -1;
901447629fSBarry Smith     else            ia[i] = app[perm[mid]];
911447629fSBarry Smith   }
921447629fSBarry Smith   PetscFunctionReturn(0);
931447629fSBarry Smith }
941447629fSBarry Smith 
951447629fSBarry Smith #undef __FUNCT__
961447629fSBarry Smith #define __FUNCT__ "AOApplicationToPetsc_Mapping"
971447629fSBarry Smith PetscErrorCode AOApplicationToPetsc_Mapping(AO ao, PetscInt n, PetscInt *ia)
981447629fSBarry Smith {
991447629fSBarry Smith   AO_Mapping *aomap = (AO_Mapping*) ao->data;
1001447629fSBarry Smith   PetscInt   *app   = aomap->app;
1011447629fSBarry Smith   PetscInt   *petsc = aomap->petsc;
1021447629fSBarry Smith   PetscInt   *perm  = aomap->appPerm;
1031447629fSBarry Smith   PetscInt   N      = aomap->N;
1041447629fSBarry Smith   PetscInt   low, high, mid=0;
1051447629fSBarry Smith   PetscInt   idex;
1061447629fSBarry Smith   PetscInt   i;
1071447629fSBarry Smith 
1081447629fSBarry Smith   /* It would be possible to use a single bisection search, which
1091447629fSBarry Smith      recursively divided the indices to be converted, and searched
1101447629fSBarry Smith      partitions which contained an index. This would result in
1111447629fSBarry Smith      better running times if indices are clustered.
1121447629fSBarry Smith   */
1131447629fSBarry Smith   PetscFunctionBegin;
1141447629fSBarry Smith   for (i = 0; i < n; i++) {
1151447629fSBarry Smith     idex = ia[i];
1161447629fSBarry Smith     if (idex < 0) continue;
1171447629fSBarry Smith     /* Use bisection since the array is sorted */
1181447629fSBarry Smith     low  = 0;
1191447629fSBarry Smith     high = N - 1;
1201447629fSBarry Smith     while (low <= high) {
1211447629fSBarry Smith       mid = (low + high)/2;
1221447629fSBarry Smith       if (idex == app[mid]) break;
1231447629fSBarry Smith       else if (idex < app[mid]) high = mid - 1;
1241447629fSBarry Smith       else low = mid + 1;
1251447629fSBarry Smith     }
1261447629fSBarry Smith     if (low > high) ia[i] = -1;
1271447629fSBarry Smith     else            ia[i] = petsc[perm[mid]];
1281447629fSBarry Smith   }
1291447629fSBarry Smith   PetscFunctionReturn(0);
1301447629fSBarry Smith }
1311447629fSBarry Smith 
1321447629fSBarry Smith static struct _AOOps AOps = {AOView_Mapping,
1331447629fSBarry Smith                              AODestroy_Mapping,
1341447629fSBarry Smith                              AOPetscToApplication_Mapping,
1351447629fSBarry Smith                              AOApplicationToPetsc_Mapping,
1361447629fSBarry Smith                              NULL,
1371447629fSBarry Smith                              NULL,
1381447629fSBarry Smith                              NULL,
1391447629fSBarry Smith                              NULL};
1401447629fSBarry Smith 
1411447629fSBarry Smith #undef __FUNCT__
1421447629fSBarry Smith #define __FUNCT__ "AOMappingHasApplicationIndex"
1431447629fSBarry Smith /*@C
1441447629fSBarry Smith   AOMappingHasApplicationIndex - Searches for the supplied application index.
1451447629fSBarry Smith 
1461447629fSBarry Smith   Input Parameters:
1471447629fSBarry Smith + ao       - The AOMapping
1481447629fSBarry Smith - index    - The application index
1491447629fSBarry Smith 
1501447629fSBarry Smith   Output Parameter:
1511447629fSBarry Smith . hasIndex - Flag is PETSC_TRUE if the index exists
1521447629fSBarry Smith 
1531447629fSBarry Smith   Level: intermediate
1541447629fSBarry Smith 
1551447629fSBarry Smith .keywords: AO, index
1561447629fSBarry Smith .seealso: AOMappingHasPetscIndex(), AOCreateMapping()
1571447629fSBarry Smith @*/
1581447629fSBarry Smith PetscErrorCode  AOMappingHasApplicationIndex(AO ao, PetscInt idex, PetscBool  *hasIndex)
1591447629fSBarry Smith {
1601447629fSBarry Smith   AO_Mapping *aomap;
1611447629fSBarry Smith   PetscInt   *app;
1621447629fSBarry Smith   PetscInt   low, high, mid;
1631447629fSBarry Smith 
1641447629fSBarry Smith   PetscFunctionBegin;
1651447629fSBarry Smith   PetscValidHeaderSpecific(ao, AO_CLASSID,1);
1661447629fSBarry Smith   PetscValidPointer(hasIndex,3);
1671447629fSBarry Smith   aomap = (AO_Mapping*) ao->data;
1681447629fSBarry Smith   app   = aomap->app;
1691447629fSBarry Smith   /* Use bisection since the array is sorted */
1701447629fSBarry Smith   low  = 0;
1711447629fSBarry Smith   high = aomap->N - 1;
1721447629fSBarry Smith   while (low <= high) {
1731447629fSBarry Smith     mid = (low + high)/2;
1741447629fSBarry Smith     if (idex == app[mid]) break;
1751447629fSBarry Smith     else if (idex < app[mid]) high = mid - 1;
1761447629fSBarry Smith     else low = mid + 1;
1771447629fSBarry Smith   }
1781447629fSBarry Smith   if (low > high) *hasIndex = PETSC_FALSE;
1791447629fSBarry Smith   else *hasIndex = PETSC_TRUE;
1801447629fSBarry Smith   PetscFunctionReturn(0);
1811447629fSBarry Smith }
1821447629fSBarry Smith 
1831447629fSBarry Smith #undef __FUNCT__
1841447629fSBarry Smith #define __FUNCT__ "AOMappingHasPetscIndex"
1851447629fSBarry Smith /*@
1861447629fSBarry Smith   AOMappingHasPetscIndex - Searches for the supplied petsc index.
1871447629fSBarry Smith 
1881447629fSBarry Smith   Input Parameters:
1891447629fSBarry Smith + ao       - The AOMapping
1901447629fSBarry Smith - index    - The petsc index
1911447629fSBarry Smith 
1921447629fSBarry Smith   Output Parameter:
1931447629fSBarry Smith . hasIndex - Flag is PETSC_TRUE if the index exists
1941447629fSBarry Smith 
1951447629fSBarry Smith   Level: intermediate
1961447629fSBarry Smith 
1971447629fSBarry Smith .keywords: AO, index
1981447629fSBarry Smith .seealso: AOMappingHasApplicationIndex(), AOCreateMapping()
1991447629fSBarry Smith @*/
2001447629fSBarry Smith PetscErrorCode  AOMappingHasPetscIndex(AO ao, PetscInt idex, PetscBool  *hasIndex)
2011447629fSBarry Smith {
2021447629fSBarry Smith   AO_Mapping *aomap;
2031447629fSBarry Smith   PetscInt   *petsc;
2041447629fSBarry Smith   PetscInt   low, high, mid;
2051447629fSBarry Smith 
2061447629fSBarry Smith   PetscFunctionBegin;
2071447629fSBarry Smith   PetscValidHeaderSpecific(ao, AO_CLASSID,1);
2081447629fSBarry Smith   PetscValidPointer(hasIndex,3);
2091447629fSBarry Smith   aomap = (AO_Mapping*) ao->data;
2101447629fSBarry Smith   petsc = aomap->petsc;
2111447629fSBarry Smith   /* Use bisection since the array is sorted */
2121447629fSBarry Smith   low  = 0;
2131447629fSBarry Smith   high = aomap->N - 1;
2141447629fSBarry Smith   while (low <= high) {
2151447629fSBarry Smith     mid = (low + high)/2;
2161447629fSBarry Smith     if (idex == petsc[mid]) break;
2171447629fSBarry Smith     else if (idex < petsc[mid]) high = mid - 1;
2181447629fSBarry Smith     else low = mid + 1;
2191447629fSBarry Smith   }
2201447629fSBarry Smith   if (low > high) *hasIndex = PETSC_FALSE;
2211447629fSBarry Smith   else *hasIndex = PETSC_TRUE;
2221447629fSBarry Smith   PetscFunctionReturn(0);
2231447629fSBarry Smith }
2241447629fSBarry Smith 
2251447629fSBarry Smith #undef __FUNCT__
2261447629fSBarry Smith #define __FUNCT__ "AOCreateMapping"
2271447629fSBarry Smith /*@C
2281447629fSBarry Smith   AOCreateMapping - Creates a basic application mapping using two integer arrays.
2291447629fSBarry Smith 
2301447629fSBarry Smith   Input Parameters:
2311447629fSBarry Smith + comm    - MPI communicator that is to share AO
2321447629fSBarry Smith . napp    - size of integer arrays
2331447629fSBarry Smith . myapp   - integer array that defines an ordering
2341447629fSBarry Smith - mypetsc - integer array that defines another ordering (may be NULL to indicate the identity ordering)
2351447629fSBarry Smith 
2361447629fSBarry Smith   Output Parameter:
2371447629fSBarry Smith . aoout   - the new application mapping
2381447629fSBarry Smith 
2391447629fSBarry Smith   Options Database Key:
2401447629fSBarry Smith $ -ao_view : call AOView() at the conclusion of AOCreateMapping()
2411447629fSBarry Smith 
2421447629fSBarry Smith   Level: beginner
2431447629fSBarry Smith 
2441447629fSBarry Smith     Notes: 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.
2451447629fSBarry Smith        Use AOCreateBasic() or AOCreateBasicIS() if they do not have holes for better performance.
2461447629fSBarry Smith 
2471447629fSBarry Smith .keywords: AO, create
2481447629fSBarry Smith .seealso: AOCreateBasic(), AOCreateBasic(), AOCreateMappingIS(), AODestroy()
2491447629fSBarry Smith @*/
2501447629fSBarry Smith PetscErrorCode  AOCreateMapping(MPI_Comm comm,PetscInt napp,const PetscInt myapp[],const PetscInt mypetsc[],AO *aoout)
2511447629fSBarry Smith {
2521447629fSBarry Smith   AO             ao;
2531447629fSBarry Smith   AO_Mapping     *aomap;
2541447629fSBarry Smith   PetscInt       *allpetsc,  *allapp;
2551447629fSBarry Smith   PetscInt       *petscPerm, *appPerm;
2561447629fSBarry Smith   PetscInt       *petsc;
2571447629fSBarry Smith   PetscMPIInt    size, rank,*lens, *disp,nnapp;
2581447629fSBarry Smith   PetscInt       N, start;
2591447629fSBarry Smith   PetscInt       i;
2601447629fSBarry Smith   PetscErrorCode ierr;
2611447629fSBarry Smith 
2621447629fSBarry Smith   PetscFunctionBegin;
2631447629fSBarry Smith   PetscValidPointer(aoout,5);
2641447629fSBarry Smith   *aoout = 0;
265607a6623SBarry Smith   ierr = AOInitializePackage();CHKERRQ(ierr);
2661447629fSBarry Smith 
2671447629fSBarry Smith   ierr     = PetscHeaderCreate(ao, _p_AO, struct _AOOps, AO_CLASSID, "AO", "Application Ordering", "AO", comm, AODestroy, AOView);CHKERRQ(ierr);
2681447629fSBarry Smith   ierr     = PetscNewLog(ao, AO_Mapping, &aomap);CHKERRQ(ierr);
2691447629fSBarry Smith   ierr     = PetscMemcpy(ao->ops, &AOps, sizeof(AOps));CHKERRQ(ierr);
2701447629fSBarry Smith   ao->data = (void*) aomap;
2711447629fSBarry Smith 
2721447629fSBarry Smith   /* transmit all lengths to all processors */
2731447629fSBarry Smith   ierr  = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
2741447629fSBarry Smith   ierr  = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
275*dcca6d9dSJed Brown   ierr  = PetscMalloc2(size, &lens,size,&disp);CHKERRQ(ierr);
2761447629fSBarry Smith   nnapp = napp;
2771447629fSBarry Smith   ierr  = MPI_Allgather(&nnapp, 1, MPI_INT, lens, 1, MPI_INT, comm);CHKERRQ(ierr);
2781447629fSBarry Smith   N     = 0;
2791447629fSBarry Smith   for (i = 0; i < size; i++) {
2801447629fSBarry Smith     disp[i] = N;
2811447629fSBarry Smith     N      += lens[i];
2821447629fSBarry Smith   }
2831447629fSBarry Smith   aomap->N = N;
2841447629fSBarry Smith   ao->N    = N;
2851447629fSBarry Smith   ao->n    = N;
2861447629fSBarry Smith 
2871447629fSBarry Smith   /* If mypetsc is 0 then use "natural" numbering */
2881447629fSBarry Smith   if (!mypetsc) {
2891447629fSBarry Smith     start = disp[rank];
2901447629fSBarry Smith     ierr  = PetscMalloc((napp+1) * sizeof(PetscInt), &petsc);CHKERRQ(ierr);
2911447629fSBarry Smith     for (i = 0; i < napp; i++) petsc[i] = start + i;
2921447629fSBarry Smith   } else {
2931447629fSBarry Smith     petsc = (PetscInt*)mypetsc;
2941447629fSBarry Smith   }
2951447629fSBarry Smith 
2961447629fSBarry Smith   /* get all indices on all processors */
297*dcca6d9dSJed Brown   ierr = PetscMalloc4(N, &allapp,N,&appPerm,N,&allpetsc,N,&petscPerm);CHKERRQ(ierr);
2981447629fSBarry Smith   ierr = MPI_Allgatherv((void*)myapp, napp, MPIU_INT, allapp,   lens, disp, MPIU_INT, comm);CHKERRQ(ierr);
2991447629fSBarry Smith   ierr = MPI_Allgatherv((void*)petsc, napp, MPIU_INT, allpetsc, lens, disp, MPIU_INT, comm);CHKERRQ(ierr);
3001447629fSBarry Smith   ierr = PetscFree2(lens,disp);CHKERRQ(ierr);
3011447629fSBarry Smith 
3021447629fSBarry Smith   /* generate a list of application and PETSc node numbers */
303*dcca6d9dSJed Brown   ierr = PetscMalloc4(N, &aomap->app,N,&aomap->appPerm,N,&aomap->petsc,N,&aomap->petscPerm);CHKERRQ(ierr);
3043bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)ao, 4*N * sizeof(PetscInt));CHKERRQ(ierr);
3051447629fSBarry Smith   for (i = 0; i < N; i++) {
3061447629fSBarry Smith     appPerm[i]   = i;
3071447629fSBarry Smith     petscPerm[i] = i;
3081447629fSBarry Smith   }
3091447629fSBarry Smith   ierr = PetscSortIntWithPermutation(N, allpetsc, petscPerm);CHKERRQ(ierr);
3101447629fSBarry Smith   ierr = PetscSortIntWithPermutation(N, allapp,   appPerm);CHKERRQ(ierr);
3111447629fSBarry Smith   /* Form sorted arrays of indices */
3121447629fSBarry Smith   for (i = 0; i < N; i++) {
3131447629fSBarry Smith     aomap->app[i]   = allapp[appPerm[i]];
3141447629fSBarry Smith     aomap->petsc[i] = allpetsc[petscPerm[i]];
3151447629fSBarry Smith   }
3161447629fSBarry Smith   /* Invert petscPerm[] into aomap->petscPerm[] */
3171447629fSBarry Smith   for (i = 0; i < N; i++) aomap->petscPerm[petscPerm[i]] = i;
3181447629fSBarry Smith 
3191447629fSBarry Smith   /* Form map between aomap->app[] and aomap->petsc[] */
3201447629fSBarry Smith   for (i = 0; i < N; i++) aomap->appPerm[i] = aomap->petscPerm[appPerm[i]];
3211447629fSBarry Smith 
3221447629fSBarry Smith   /* Invert appPerm[] into allapp[] */
3231447629fSBarry Smith   for (i = 0; i < N; i++) allapp[appPerm[i]] = i;
3241447629fSBarry Smith 
3251447629fSBarry Smith   /* Form map between aomap->petsc[] and aomap->app[] */
3261447629fSBarry Smith   for (i = 0; i < N; i++) aomap->petscPerm[i] = allapp[petscPerm[i]];
3271447629fSBarry Smith 
3281447629fSBarry Smith #if defined(PETSC_USE_DEBUG)
3291447629fSBarry Smith   /* Check that the permutations are complementary */
3301447629fSBarry Smith   for (i = 0; i < N; i++) {
3311447629fSBarry Smith     if (i != aomap->appPerm[aomap->petscPerm[i]]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Invalid ordering");
3321447629fSBarry Smith   }
3331447629fSBarry Smith #endif
3341447629fSBarry Smith   /* Cleanup */
3351447629fSBarry Smith   if (!mypetsc) {
3361447629fSBarry Smith     ierr = PetscFree(petsc);CHKERRQ(ierr);
3371447629fSBarry Smith   }
3381447629fSBarry Smith   ierr = PetscFree4(allapp,appPerm,allpetsc,petscPerm);CHKERRQ(ierr);
3391447629fSBarry Smith 
340817ea411SJed Brown   ierr = AOViewFromOptions(ao,NULL,"-ao_view");CHKERRQ(ierr);
3411447629fSBarry Smith 
3421447629fSBarry Smith   *aoout = ao;
3431447629fSBarry Smith   PetscFunctionReturn(0);
3441447629fSBarry Smith }
3451447629fSBarry Smith 
3461447629fSBarry Smith #undef __FUNCT__
3471447629fSBarry Smith #define __FUNCT__ "AOCreateMappingIS"
3481447629fSBarry Smith /*@C
3491447629fSBarry Smith   AOCreateMappingIS - Creates a basic application ordering using two index sets.
3501447629fSBarry Smith 
3511447629fSBarry Smith   Input Parameters:
3521447629fSBarry Smith + comm    - MPI communicator that is to share AO
3531447629fSBarry Smith . isapp   - index set that defines an ordering
3541447629fSBarry Smith - ispetsc - index set that defines another ordering, maybe NULL for identity IS
3551447629fSBarry Smith 
3561447629fSBarry Smith   Output Parameter:
3571447629fSBarry Smith . aoout   - the new application ordering
3581447629fSBarry Smith 
3591447629fSBarry Smith   Options Database Key:
3601447629fSBarry Smith $ -ao_view : call AOView() at the conclusion of AOCreateMappingIS()
3611447629fSBarry Smith 
3621447629fSBarry Smith   Level: beginner
3631447629fSBarry Smith 
3641447629fSBarry Smith     Notes: 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.
3651447629fSBarry Smith        Use AOCreateBasic() or AOCreateBasicIS() if they do not have holes for better performance.
3661447629fSBarry Smith 
3671447629fSBarry Smith .keywords: AO, create
3681447629fSBarry Smith .seealso: AOCreateBasic(), AOCreateMapping(), AODestroy()
3691447629fSBarry Smith @*/
3701447629fSBarry Smith PetscErrorCode  AOCreateMappingIS(IS isapp, IS ispetsc, AO *aoout)
3711447629fSBarry Smith {
3721447629fSBarry Smith   MPI_Comm       comm;
3731447629fSBarry Smith   const PetscInt *mypetsc, *myapp;
3741447629fSBarry Smith   PetscInt       napp, npetsc;
3751447629fSBarry Smith   PetscErrorCode ierr;
3761447629fSBarry Smith 
3771447629fSBarry Smith   PetscFunctionBegin;
3781447629fSBarry Smith   ierr = PetscObjectGetComm((PetscObject) isapp, &comm);CHKERRQ(ierr);
3791447629fSBarry Smith   ierr = ISGetLocalSize(isapp, &napp);CHKERRQ(ierr);
3801447629fSBarry Smith   if (ispetsc) {
3811447629fSBarry Smith     ierr = ISGetLocalSize(ispetsc, &npetsc);CHKERRQ(ierr);
3821447629fSBarry Smith     if (napp != npetsc) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "Local IS lengths must match");
3831447629fSBarry Smith     ierr = ISGetIndices(ispetsc, &mypetsc);CHKERRQ(ierr);
3841447629fSBarry Smith   } else {
3851447629fSBarry Smith     mypetsc = NULL;
3861447629fSBarry Smith   }
3871447629fSBarry Smith   ierr = ISGetIndices(isapp, &myapp);CHKERRQ(ierr);
3881447629fSBarry Smith 
3891447629fSBarry Smith   ierr = AOCreateMapping(comm, napp, myapp, mypetsc, aoout);CHKERRQ(ierr);
3901447629fSBarry Smith 
3911447629fSBarry Smith   ierr = ISRestoreIndices(isapp, &myapp);CHKERRQ(ierr);
3921447629fSBarry Smith   if (ispetsc) {
3931447629fSBarry Smith     ierr = ISRestoreIndices(ispetsc, &mypetsc);CHKERRQ(ierr);
3941447629fSBarry Smith   }
3951447629fSBarry Smith   PetscFunctionReturn(0);
3961447629fSBarry Smith }
397