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