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