1*1447629fSBarry Smith 2*1447629fSBarry Smith /* 3*1447629fSBarry Smith These AO application ordering routines do not require that the input 4*1447629fSBarry Smith be a permutation, but merely a 1-1 mapping. This implementation still 5*1447629fSBarry Smith keeps the entire ordering on each processor. 6*1447629fSBarry Smith */ 7*1447629fSBarry Smith 8*1447629fSBarry Smith #include <../src/vec/is/ao/aoimpl.h> /*I "petscao.h" I*/ 9*1447629fSBarry Smith 10*1447629fSBarry Smith typedef struct { 11*1447629fSBarry Smith PetscInt N; 12*1447629fSBarry Smith PetscInt *app; /* app[i] is the partner for petsc[appPerm[i]] */ 13*1447629fSBarry Smith PetscInt *appPerm; 14*1447629fSBarry Smith PetscInt *petsc; /* petsc[j] is the partner for app[petscPerm[j]] */ 15*1447629fSBarry Smith PetscInt *petscPerm; 16*1447629fSBarry Smith } AO_Mapping; 17*1447629fSBarry Smith 18*1447629fSBarry Smith #undef __FUNCT__ 19*1447629fSBarry Smith #define __FUNCT__ "AODestroy_Mapping" 20*1447629fSBarry Smith PetscErrorCode AODestroy_Mapping(AO ao) 21*1447629fSBarry Smith { 22*1447629fSBarry Smith AO_Mapping *aomap = (AO_Mapping*) ao->data; 23*1447629fSBarry Smith PetscErrorCode ierr; 24*1447629fSBarry Smith 25*1447629fSBarry Smith PetscFunctionBegin; 26*1447629fSBarry Smith ierr = PetscFree4(aomap->app,aomap->appPerm,aomap->petsc,aomap->petscPerm);CHKERRQ(ierr); 27*1447629fSBarry Smith ierr = PetscFree(aomap);CHKERRQ(ierr); 28*1447629fSBarry Smith PetscFunctionReturn(0); 29*1447629fSBarry Smith } 30*1447629fSBarry Smith 31*1447629fSBarry Smith #undef __FUNCT__ 32*1447629fSBarry Smith #define __FUNCT__ "AOView_Mapping" 33*1447629fSBarry Smith PetscErrorCode AOView_Mapping(AO ao, PetscViewer viewer) 34*1447629fSBarry Smith { 35*1447629fSBarry Smith AO_Mapping *aomap = (AO_Mapping*) ao->data; 36*1447629fSBarry Smith PetscMPIInt rank; 37*1447629fSBarry Smith PetscInt i; 38*1447629fSBarry Smith PetscBool iascii; 39*1447629fSBarry Smith PetscErrorCode ierr; 40*1447629fSBarry Smith 41*1447629fSBarry Smith PetscFunctionBegin; 42*1447629fSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)ao), &rank);CHKERRQ(ierr); 43*1447629fSBarry Smith if (rank) PetscFunctionReturn(0); 44*1447629fSBarry Smith 45*1447629fSBarry Smith if (!viewer) viewer = PETSC_VIEWER_STDOUT_SELF; 46*1447629fSBarry Smith 47*1447629fSBarry Smith ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 48*1447629fSBarry Smith if (iascii) { 49*1447629fSBarry Smith PetscViewerASCIIPrintf(viewer, "Number of elements in ordering %D\n", aomap->N); 50*1447629fSBarry Smith PetscViewerASCIIPrintf(viewer, " App. PETSc\n"); 51*1447629fSBarry Smith for (i = 0; i < aomap->N; i++) { 52*1447629fSBarry Smith PetscViewerASCIIPrintf(viewer, "%D %D %D\n", i, aomap->app[i], aomap->petsc[aomap->appPerm[i]]); 53*1447629fSBarry Smith } 54*1447629fSBarry Smith } 55*1447629fSBarry Smith PetscFunctionReturn(0); 56*1447629fSBarry Smith } 57*1447629fSBarry Smith 58*1447629fSBarry Smith #undef __FUNCT__ 59*1447629fSBarry Smith #define __FUNCT__ "AOPetscToApplication_Mapping" 60*1447629fSBarry Smith PetscErrorCode AOPetscToApplication_Mapping(AO ao, PetscInt n, PetscInt *ia) 61*1447629fSBarry Smith { 62*1447629fSBarry Smith AO_Mapping *aomap = (AO_Mapping*) ao->data; 63*1447629fSBarry Smith PetscInt *app = aomap->app; 64*1447629fSBarry Smith PetscInt *petsc = aomap->petsc; 65*1447629fSBarry Smith PetscInt *perm = aomap->petscPerm; 66*1447629fSBarry Smith PetscInt N = aomap->N; 67*1447629fSBarry Smith PetscInt low, high, mid=0; 68*1447629fSBarry Smith PetscInt idex; 69*1447629fSBarry Smith PetscInt i; 70*1447629fSBarry Smith 71*1447629fSBarry Smith /* It would be possible to use a single bisection search, which 72*1447629fSBarry Smith recursively divided the indices to be converted, and searched 73*1447629fSBarry Smith partitions which contained an index. This would result in 74*1447629fSBarry Smith better running times if indices are clustered. 75*1447629fSBarry Smith */ 76*1447629fSBarry Smith PetscFunctionBegin; 77*1447629fSBarry Smith for (i = 0; i < n; i++) { 78*1447629fSBarry Smith idex = ia[i]; 79*1447629fSBarry Smith if (idex < 0) continue; 80*1447629fSBarry Smith /* Use bisection since the array is sorted */ 81*1447629fSBarry Smith low = 0; 82*1447629fSBarry Smith high = N - 1; 83*1447629fSBarry Smith while (low <= high) { 84*1447629fSBarry Smith mid = (low + high)/2; 85*1447629fSBarry Smith if (idex == petsc[mid]) break; 86*1447629fSBarry Smith else if (idex < petsc[mid]) high = mid - 1; 87*1447629fSBarry Smith else low = mid + 1; 88*1447629fSBarry Smith } 89*1447629fSBarry Smith if (low > high) ia[i] = -1; 90*1447629fSBarry Smith else ia[i] = app[perm[mid]]; 91*1447629fSBarry Smith } 92*1447629fSBarry Smith PetscFunctionReturn(0); 93*1447629fSBarry Smith } 94*1447629fSBarry Smith 95*1447629fSBarry Smith #undef __FUNCT__ 96*1447629fSBarry Smith #define __FUNCT__ "AOApplicationToPetsc_Mapping" 97*1447629fSBarry Smith PetscErrorCode AOApplicationToPetsc_Mapping(AO ao, PetscInt n, PetscInt *ia) 98*1447629fSBarry Smith { 99*1447629fSBarry Smith AO_Mapping *aomap = (AO_Mapping*) ao->data; 100*1447629fSBarry Smith PetscInt *app = aomap->app; 101*1447629fSBarry Smith PetscInt *petsc = aomap->petsc; 102*1447629fSBarry Smith PetscInt *perm = aomap->appPerm; 103*1447629fSBarry Smith PetscInt N = aomap->N; 104*1447629fSBarry Smith PetscInt low, high, mid=0; 105*1447629fSBarry Smith PetscInt idex; 106*1447629fSBarry Smith PetscInt i; 107*1447629fSBarry Smith 108*1447629fSBarry Smith /* It would be possible to use a single bisection search, which 109*1447629fSBarry Smith recursively divided the indices to be converted, and searched 110*1447629fSBarry Smith partitions which contained an index. This would result in 111*1447629fSBarry Smith better running times if indices are clustered. 112*1447629fSBarry Smith */ 113*1447629fSBarry Smith PetscFunctionBegin; 114*1447629fSBarry Smith for (i = 0; i < n; i++) { 115*1447629fSBarry Smith idex = ia[i]; 116*1447629fSBarry Smith if (idex < 0) continue; 117*1447629fSBarry Smith /* Use bisection since the array is sorted */ 118*1447629fSBarry Smith low = 0; 119*1447629fSBarry Smith high = N - 1; 120*1447629fSBarry Smith while (low <= high) { 121*1447629fSBarry Smith mid = (low + high)/2; 122*1447629fSBarry Smith if (idex == app[mid]) break; 123*1447629fSBarry Smith else if (idex < app[mid]) high = mid - 1; 124*1447629fSBarry Smith else low = mid + 1; 125*1447629fSBarry Smith } 126*1447629fSBarry Smith if (low > high) ia[i] = -1; 127*1447629fSBarry Smith else ia[i] = petsc[perm[mid]]; 128*1447629fSBarry Smith } 129*1447629fSBarry Smith PetscFunctionReturn(0); 130*1447629fSBarry Smith } 131*1447629fSBarry Smith 132*1447629fSBarry Smith static struct _AOOps AOps = {AOView_Mapping, 133*1447629fSBarry Smith AODestroy_Mapping, 134*1447629fSBarry Smith AOPetscToApplication_Mapping, 135*1447629fSBarry Smith AOApplicationToPetsc_Mapping, 136*1447629fSBarry Smith NULL, 137*1447629fSBarry Smith NULL, 138*1447629fSBarry Smith NULL, 139*1447629fSBarry Smith NULL}; 140*1447629fSBarry Smith 141*1447629fSBarry Smith #undef __FUNCT__ 142*1447629fSBarry Smith #define __FUNCT__ "AOMappingHasApplicationIndex" 143*1447629fSBarry Smith /*@C 144*1447629fSBarry Smith AOMappingHasApplicationIndex - Searches for the supplied application index. 145*1447629fSBarry Smith 146*1447629fSBarry Smith Input Parameters: 147*1447629fSBarry Smith + ao - The AOMapping 148*1447629fSBarry Smith - index - The application index 149*1447629fSBarry Smith 150*1447629fSBarry Smith Output Parameter: 151*1447629fSBarry Smith . hasIndex - Flag is PETSC_TRUE if the index exists 152*1447629fSBarry Smith 153*1447629fSBarry Smith Level: intermediate 154*1447629fSBarry Smith 155*1447629fSBarry Smith .keywords: AO, index 156*1447629fSBarry Smith .seealso: AOMappingHasPetscIndex(), AOCreateMapping() 157*1447629fSBarry Smith @*/ 158*1447629fSBarry Smith PetscErrorCode AOMappingHasApplicationIndex(AO ao, PetscInt idex, PetscBool *hasIndex) 159*1447629fSBarry Smith { 160*1447629fSBarry Smith AO_Mapping *aomap; 161*1447629fSBarry Smith PetscInt *app; 162*1447629fSBarry Smith PetscInt low, high, mid; 163*1447629fSBarry Smith 164*1447629fSBarry Smith PetscFunctionBegin; 165*1447629fSBarry Smith PetscValidHeaderSpecific(ao, AO_CLASSID,1); 166*1447629fSBarry Smith PetscValidPointer(hasIndex,3); 167*1447629fSBarry Smith aomap = (AO_Mapping*) ao->data; 168*1447629fSBarry Smith app = aomap->app; 169*1447629fSBarry Smith /* Use bisection since the array is sorted */ 170*1447629fSBarry Smith low = 0; 171*1447629fSBarry Smith high = aomap->N - 1; 172*1447629fSBarry Smith while (low <= high) { 173*1447629fSBarry Smith mid = (low + high)/2; 174*1447629fSBarry Smith if (idex == app[mid]) break; 175*1447629fSBarry Smith else if (idex < app[mid]) high = mid - 1; 176*1447629fSBarry Smith else low = mid + 1; 177*1447629fSBarry Smith } 178*1447629fSBarry Smith if (low > high) *hasIndex = PETSC_FALSE; 179*1447629fSBarry Smith else *hasIndex = PETSC_TRUE; 180*1447629fSBarry Smith PetscFunctionReturn(0); 181*1447629fSBarry Smith } 182*1447629fSBarry Smith 183*1447629fSBarry Smith #undef __FUNCT__ 184*1447629fSBarry Smith #define __FUNCT__ "AOMappingHasPetscIndex" 185*1447629fSBarry Smith /*@ 186*1447629fSBarry Smith AOMappingHasPetscIndex - Searches for the supplied petsc index. 187*1447629fSBarry Smith 188*1447629fSBarry Smith Input Parameters: 189*1447629fSBarry Smith + ao - The AOMapping 190*1447629fSBarry Smith - index - The petsc index 191*1447629fSBarry Smith 192*1447629fSBarry Smith Output Parameter: 193*1447629fSBarry Smith . hasIndex - Flag is PETSC_TRUE if the index exists 194*1447629fSBarry Smith 195*1447629fSBarry Smith Level: intermediate 196*1447629fSBarry Smith 197*1447629fSBarry Smith .keywords: AO, index 198*1447629fSBarry Smith .seealso: AOMappingHasApplicationIndex(), AOCreateMapping() 199*1447629fSBarry Smith @*/ 200*1447629fSBarry Smith PetscErrorCode AOMappingHasPetscIndex(AO ao, PetscInt idex, PetscBool *hasIndex) 201*1447629fSBarry Smith { 202*1447629fSBarry Smith AO_Mapping *aomap; 203*1447629fSBarry Smith PetscInt *petsc; 204*1447629fSBarry Smith PetscInt low, high, mid; 205*1447629fSBarry Smith 206*1447629fSBarry Smith PetscFunctionBegin; 207*1447629fSBarry Smith PetscValidHeaderSpecific(ao, AO_CLASSID,1); 208*1447629fSBarry Smith PetscValidPointer(hasIndex,3); 209*1447629fSBarry Smith aomap = (AO_Mapping*) ao->data; 210*1447629fSBarry Smith petsc = aomap->petsc; 211*1447629fSBarry Smith /* Use bisection since the array is sorted */ 212*1447629fSBarry Smith low = 0; 213*1447629fSBarry Smith high = aomap->N - 1; 214*1447629fSBarry Smith while (low <= high) { 215*1447629fSBarry Smith mid = (low + high)/2; 216*1447629fSBarry Smith if (idex == petsc[mid]) break; 217*1447629fSBarry Smith else if (idex < petsc[mid]) high = mid - 1; 218*1447629fSBarry Smith else low = mid + 1; 219*1447629fSBarry Smith } 220*1447629fSBarry Smith if (low > high) *hasIndex = PETSC_FALSE; 221*1447629fSBarry Smith else *hasIndex = PETSC_TRUE; 222*1447629fSBarry Smith PetscFunctionReturn(0); 223*1447629fSBarry Smith } 224*1447629fSBarry Smith 225*1447629fSBarry Smith #undef __FUNCT__ 226*1447629fSBarry Smith #define __FUNCT__ "AOCreateMapping" 227*1447629fSBarry Smith /*@C 228*1447629fSBarry Smith AOCreateMapping - Creates a basic application mapping using two integer arrays. 229*1447629fSBarry Smith 230*1447629fSBarry Smith Input Parameters: 231*1447629fSBarry Smith + comm - MPI communicator that is to share AO 232*1447629fSBarry Smith . napp - size of integer arrays 233*1447629fSBarry Smith . myapp - integer array that defines an ordering 234*1447629fSBarry Smith - mypetsc - integer array that defines another ordering (may be NULL to indicate the identity ordering) 235*1447629fSBarry Smith 236*1447629fSBarry Smith Output Parameter: 237*1447629fSBarry Smith . aoout - the new application mapping 238*1447629fSBarry Smith 239*1447629fSBarry Smith Options Database Key: 240*1447629fSBarry Smith $ -ao_view : call AOView() at the conclusion of AOCreateMapping() 241*1447629fSBarry Smith 242*1447629fSBarry Smith Level: beginner 243*1447629fSBarry Smith 244*1447629fSBarry 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. 245*1447629fSBarry Smith Use AOCreateBasic() or AOCreateBasicIS() if they do not have holes for better performance. 246*1447629fSBarry Smith 247*1447629fSBarry Smith .keywords: AO, create 248*1447629fSBarry Smith .seealso: AOCreateBasic(), AOCreateBasic(), AOCreateMappingIS(), AODestroy() 249*1447629fSBarry Smith @*/ 250*1447629fSBarry Smith PetscErrorCode AOCreateMapping(MPI_Comm comm,PetscInt napp,const PetscInt myapp[],const PetscInt mypetsc[],AO *aoout) 251*1447629fSBarry Smith { 252*1447629fSBarry Smith AO ao; 253*1447629fSBarry Smith AO_Mapping *aomap; 254*1447629fSBarry Smith PetscInt *allpetsc, *allapp; 255*1447629fSBarry Smith PetscInt *petscPerm, *appPerm; 256*1447629fSBarry Smith PetscInt *petsc; 257*1447629fSBarry Smith PetscMPIInt size, rank,*lens, *disp,nnapp; 258*1447629fSBarry Smith PetscInt N, start; 259*1447629fSBarry Smith PetscInt i; 260*1447629fSBarry Smith PetscBool opt; 261*1447629fSBarry Smith PetscErrorCode ierr; 262*1447629fSBarry Smith 263*1447629fSBarry Smith PetscFunctionBegin; 264*1447629fSBarry Smith PetscValidPointer(aoout,5); 265*1447629fSBarry Smith *aoout = 0; 266*1447629fSBarry Smith #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 267*1447629fSBarry Smith ierr = AOInitializePackage(NULL);CHKERRQ(ierr); 268*1447629fSBarry Smith #endif 269*1447629fSBarry Smith 270*1447629fSBarry Smith ierr = PetscHeaderCreate(ao, _p_AO, struct _AOOps, AO_CLASSID, "AO", "Application Ordering", "AO", comm, AODestroy, AOView);CHKERRQ(ierr); 271*1447629fSBarry Smith ierr = PetscNewLog(ao, AO_Mapping, &aomap);CHKERRQ(ierr); 272*1447629fSBarry Smith ierr = PetscMemcpy(ao->ops, &AOps, sizeof(AOps));CHKERRQ(ierr); 273*1447629fSBarry Smith ao->data = (void*) aomap; 274*1447629fSBarry Smith 275*1447629fSBarry Smith /* transmit all lengths to all processors */ 276*1447629fSBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 277*1447629fSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 278*1447629fSBarry Smith ierr = PetscMalloc2(size,PetscMPIInt, &lens,size,PetscMPIInt,&disp);CHKERRQ(ierr); 279*1447629fSBarry Smith nnapp = napp; 280*1447629fSBarry Smith ierr = MPI_Allgather(&nnapp, 1, MPI_INT, lens, 1, MPI_INT, comm);CHKERRQ(ierr); 281*1447629fSBarry Smith N = 0; 282*1447629fSBarry Smith for (i = 0; i < size; i++) { 283*1447629fSBarry Smith disp[i] = N; 284*1447629fSBarry Smith N += lens[i]; 285*1447629fSBarry Smith } 286*1447629fSBarry Smith aomap->N = N; 287*1447629fSBarry Smith ao->N = N; 288*1447629fSBarry Smith ao->n = N; 289*1447629fSBarry Smith 290*1447629fSBarry Smith /* If mypetsc is 0 then use "natural" numbering */ 291*1447629fSBarry Smith if (!mypetsc) { 292*1447629fSBarry Smith start = disp[rank]; 293*1447629fSBarry Smith ierr = PetscMalloc((napp+1) * sizeof(PetscInt), &petsc);CHKERRQ(ierr); 294*1447629fSBarry Smith for (i = 0; i < napp; i++) petsc[i] = start + i; 295*1447629fSBarry Smith } else { 296*1447629fSBarry Smith petsc = (PetscInt*)mypetsc; 297*1447629fSBarry Smith } 298*1447629fSBarry Smith 299*1447629fSBarry Smith /* get all indices on all processors */ 300*1447629fSBarry Smith ierr = PetscMalloc4(N,PetscInt, &allapp,N,PetscInt,&appPerm,N,PetscInt,&allpetsc,N,PetscInt,&petscPerm);CHKERRQ(ierr); 301*1447629fSBarry Smith ierr = MPI_Allgatherv((void*)myapp, napp, MPIU_INT, allapp, lens, disp, MPIU_INT, comm);CHKERRQ(ierr); 302*1447629fSBarry Smith ierr = MPI_Allgatherv((void*)petsc, napp, MPIU_INT, allpetsc, lens, disp, MPIU_INT, comm);CHKERRQ(ierr); 303*1447629fSBarry Smith ierr = PetscFree2(lens,disp);CHKERRQ(ierr); 304*1447629fSBarry Smith 305*1447629fSBarry Smith /* generate a list of application and PETSc node numbers */ 306*1447629fSBarry Smith ierr = PetscMalloc4(N,PetscInt, &aomap->app,N,PetscInt,&aomap->appPerm,N,PetscInt,&aomap->petsc,N,PetscInt,&aomap->petscPerm);CHKERRQ(ierr); 307*1447629fSBarry Smith ierr = PetscLogObjectMemory(ao, 4*N * sizeof(PetscInt));CHKERRQ(ierr); 308*1447629fSBarry Smith for (i = 0; i < N; i++) { 309*1447629fSBarry Smith appPerm[i] = i; 310*1447629fSBarry Smith petscPerm[i] = i; 311*1447629fSBarry Smith } 312*1447629fSBarry Smith ierr = PetscSortIntWithPermutation(N, allpetsc, petscPerm);CHKERRQ(ierr); 313*1447629fSBarry Smith ierr = PetscSortIntWithPermutation(N, allapp, appPerm);CHKERRQ(ierr); 314*1447629fSBarry Smith /* Form sorted arrays of indices */ 315*1447629fSBarry Smith for (i = 0; i < N; i++) { 316*1447629fSBarry Smith aomap->app[i] = allapp[appPerm[i]]; 317*1447629fSBarry Smith aomap->petsc[i] = allpetsc[petscPerm[i]]; 318*1447629fSBarry Smith } 319*1447629fSBarry Smith /* Invert petscPerm[] into aomap->petscPerm[] */ 320*1447629fSBarry Smith for (i = 0; i < N; i++) aomap->petscPerm[petscPerm[i]] = i; 321*1447629fSBarry Smith 322*1447629fSBarry Smith /* Form map between aomap->app[] and aomap->petsc[] */ 323*1447629fSBarry Smith for (i = 0; i < N; i++) aomap->appPerm[i] = aomap->petscPerm[appPerm[i]]; 324*1447629fSBarry Smith 325*1447629fSBarry Smith /* Invert appPerm[] into allapp[] */ 326*1447629fSBarry Smith for (i = 0; i < N; i++) allapp[appPerm[i]] = i; 327*1447629fSBarry Smith 328*1447629fSBarry Smith /* Form map between aomap->petsc[] and aomap->app[] */ 329*1447629fSBarry Smith for (i = 0; i < N; i++) aomap->petscPerm[i] = allapp[petscPerm[i]]; 330*1447629fSBarry Smith 331*1447629fSBarry Smith #if defined(PETSC_USE_DEBUG) 332*1447629fSBarry Smith /* Check that the permutations are complementary */ 333*1447629fSBarry Smith for (i = 0; i < N; i++) { 334*1447629fSBarry Smith if (i != aomap->appPerm[aomap->petscPerm[i]]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Invalid ordering"); 335*1447629fSBarry Smith } 336*1447629fSBarry Smith #endif 337*1447629fSBarry Smith /* Cleanup */ 338*1447629fSBarry Smith if (!mypetsc) { 339*1447629fSBarry Smith ierr = PetscFree(petsc);CHKERRQ(ierr); 340*1447629fSBarry Smith } 341*1447629fSBarry Smith ierr = PetscFree4(allapp,appPerm,allpetsc,petscPerm);CHKERRQ(ierr); 342*1447629fSBarry Smith 343*1447629fSBarry Smith opt = PETSC_FALSE; 344*1447629fSBarry Smith ierr = PetscOptionsGetBool(NULL, "-ao_view", &opt,NULL);CHKERRQ(ierr); 345*1447629fSBarry Smith if (opt) { 346*1447629fSBarry Smith ierr = AOView(ao, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 347*1447629fSBarry Smith } 348*1447629fSBarry Smith 349*1447629fSBarry Smith *aoout = ao; 350*1447629fSBarry Smith PetscFunctionReturn(0); 351*1447629fSBarry Smith } 352*1447629fSBarry Smith 353*1447629fSBarry Smith #undef __FUNCT__ 354*1447629fSBarry Smith #define __FUNCT__ "AOCreateMappingIS" 355*1447629fSBarry Smith /*@C 356*1447629fSBarry Smith AOCreateMappingIS - Creates a basic application ordering using two index sets. 357*1447629fSBarry Smith 358*1447629fSBarry Smith Input Parameters: 359*1447629fSBarry Smith + comm - MPI communicator that is to share AO 360*1447629fSBarry Smith . isapp - index set that defines an ordering 361*1447629fSBarry Smith - ispetsc - index set that defines another ordering, maybe NULL for identity IS 362*1447629fSBarry Smith 363*1447629fSBarry Smith Output Parameter: 364*1447629fSBarry Smith . aoout - the new application ordering 365*1447629fSBarry Smith 366*1447629fSBarry Smith Options Database Key: 367*1447629fSBarry Smith $ -ao_view : call AOView() at the conclusion of AOCreateMappingIS() 368*1447629fSBarry Smith 369*1447629fSBarry Smith Level: beginner 370*1447629fSBarry Smith 371*1447629fSBarry 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. 372*1447629fSBarry Smith Use AOCreateBasic() or AOCreateBasicIS() if they do not have holes for better performance. 373*1447629fSBarry Smith 374*1447629fSBarry Smith .keywords: AO, create 375*1447629fSBarry Smith .seealso: AOCreateBasic(), AOCreateMapping(), AODestroy() 376*1447629fSBarry Smith @*/ 377*1447629fSBarry Smith PetscErrorCode AOCreateMappingIS(IS isapp, IS ispetsc, AO *aoout) 378*1447629fSBarry Smith { 379*1447629fSBarry Smith MPI_Comm comm; 380*1447629fSBarry Smith const PetscInt *mypetsc, *myapp; 381*1447629fSBarry Smith PetscInt napp, npetsc; 382*1447629fSBarry Smith PetscErrorCode ierr; 383*1447629fSBarry Smith 384*1447629fSBarry Smith PetscFunctionBegin; 385*1447629fSBarry Smith ierr = PetscObjectGetComm((PetscObject) isapp, &comm);CHKERRQ(ierr); 386*1447629fSBarry Smith ierr = ISGetLocalSize(isapp, &napp);CHKERRQ(ierr); 387*1447629fSBarry Smith if (ispetsc) { 388*1447629fSBarry Smith ierr = ISGetLocalSize(ispetsc, &npetsc);CHKERRQ(ierr); 389*1447629fSBarry Smith if (napp != npetsc) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "Local IS lengths must match"); 390*1447629fSBarry Smith ierr = ISGetIndices(ispetsc, &mypetsc);CHKERRQ(ierr); 391*1447629fSBarry Smith } else { 392*1447629fSBarry Smith mypetsc = NULL; 393*1447629fSBarry Smith } 394*1447629fSBarry Smith ierr = ISGetIndices(isapp, &myapp);CHKERRQ(ierr); 395*1447629fSBarry Smith 396*1447629fSBarry Smith ierr = AOCreateMapping(comm, napp, myapp, mypetsc, aoout);CHKERRQ(ierr); 397*1447629fSBarry Smith 398*1447629fSBarry Smith ierr = ISRestoreIndices(isapp, &myapp);CHKERRQ(ierr); 399*1447629fSBarry Smith if (ispetsc) { 400*1447629fSBarry Smith ierr = ISRestoreIndices(ispetsc, &mypetsc);CHKERRQ(ierr); 401*1447629fSBarry Smith } 402*1447629fSBarry Smith PetscFunctionReturn(0); 403*1447629fSBarry Smith } 404