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