1 2 /* 3 The memory scalable AO application ordering routines. These store the 4 orderings on each processor for that processor's range of values 5 */ 6 7 #include <../src/vec/is/ao/aoimpl.h> /*I "petscao.h" I*/ 8 9 typedef struct { 10 PetscInt *app_loc; /* app_loc[i] is the partner for the ith local PETSc slot */ 11 PetscInt *petsc_loc; /* petsc_loc[j] is the partner for the jth local app slot */ 12 PetscLayout map; /* determines the local sizes of ao */ 13 } AO_MemoryScalable; 14 15 /* 16 All processors ship the data to process 0 to be printed; note that this is not scalable because 17 process 0 allocates space for all the orderings entry across all the processes 18 */ 19 PetscErrorCode AOView_MemoryScalable(AO ao, PetscViewer viewer) { 20 PetscMPIInt rank, size; 21 AO_MemoryScalable *aomems = (AO_MemoryScalable *)ao->data; 22 PetscBool iascii; 23 PetscMPIInt tag_app, tag_petsc; 24 PetscLayout map = aomems->map; 25 PetscInt *app, *app_loc, *petsc, *petsc_loc, len, i, j; 26 MPI_Status status; 27 28 PetscFunctionBegin; 29 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 30 PetscCheck(iascii, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Viewer type %s not supported for AO MemoryScalable", ((PetscObject)viewer)->type_name); 31 32 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)ao), &rank)); 33 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)ao), &size)); 34 35 PetscCall(PetscObjectGetNewTag((PetscObject)ao, &tag_app)); 36 PetscCall(PetscObjectGetNewTag((PetscObject)ao, &tag_petsc)); 37 38 if (rank == 0) { 39 PetscCall(PetscViewerASCIIPrintf(viewer, "Number of elements in ordering %" PetscInt_FMT "\n", ao->N)); 40 PetscCall(PetscViewerASCIIPrintf(viewer, "PETSc->App App->PETSc\n")); 41 42 PetscCall(PetscMalloc2(map->N, &app, map->N, &petsc)); 43 len = map->n; 44 /* print local AO */ 45 PetscCall(PetscViewerASCIIPrintf(viewer, "Process [%d]\n", rank)); 46 for (i = 0; i < len; i++) { PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " %3" PetscInt_FMT " %3" PetscInt_FMT " %3" PetscInt_FMT "\n", i, aomems->app_loc[i], i, aomems->petsc_loc[i])); } 47 48 /* recv and print off-processor's AO */ 49 for (i = 1; i < size; i++) { 50 len = map->range[i + 1] - map->range[i]; 51 app_loc = app + map->range[i]; 52 petsc_loc = petsc + map->range[i]; 53 PetscCallMPI(MPI_Recv(app_loc, (PetscMPIInt)len, MPIU_INT, i, tag_app, PetscObjectComm((PetscObject)ao), &status)); 54 PetscCallMPI(MPI_Recv(petsc_loc, (PetscMPIInt)len, MPIU_INT, i, tag_petsc, PetscObjectComm((PetscObject)ao), &status)); 55 PetscCall(PetscViewerASCIIPrintf(viewer, "Process [%" PetscInt_FMT "]\n", i)); 56 for (j = 0; j < len; j++) { PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " %3" PetscInt_FMT " %3" PetscInt_FMT " %3" PetscInt_FMT "\n", map->range[i] + j, app_loc[j], map->range[i] + j, petsc_loc[j])); } 57 } 58 PetscCall(PetscFree2(app, petsc)); 59 60 } else { 61 /* send values */ 62 PetscCallMPI(MPI_Send((void *)aomems->app_loc, map->n, MPIU_INT, 0, tag_app, PetscObjectComm((PetscObject)ao))); 63 PetscCallMPI(MPI_Send((void *)aomems->petsc_loc, map->n, MPIU_INT, 0, tag_petsc, PetscObjectComm((PetscObject)ao))); 64 } 65 PetscCall(PetscViewerFlush(viewer)); 66 PetscFunctionReturn(0); 67 } 68 69 PetscErrorCode AODestroy_MemoryScalable(AO ao) { 70 AO_MemoryScalable *aomems = (AO_MemoryScalable *)ao->data; 71 72 PetscFunctionBegin; 73 PetscCall(PetscFree2(aomems->app_loc, aomems->petsc_loc)); 74 PetscCall(PetscLayoutDestroy(&aomems->map)); 75 PetscCall(PetscFree(aomems)); 76 PetscFunctionReturn(0); 77 } 78 79 /* 80 Input Parameters: 81 + ao - the application ordering context 82 . n - the number of integers in ia[] 83 . ia - the integers; these are replaced with their mapped value 84 - maploc - app_loc or petsc_loc in struct "AO_MemoryScalable" 85 86 Output Parameter: 87 . ia - the mapped interges 88 */ 89 PetscErrorCode AOMap_MemoryScalable_private(AO ao, PetscInt n, PetscInt *ia, const PetscInt *maploc) { 90 AO_MemoryScalable *aomems = (AO_MemoryScalable *)ao->data; 91 MPI_Comm comm; 92 PetscMPIInt rank, size, tag1, tag2; 93 PetscInt *owner, *start, *sizes, nsends, nreceives; 94 PetscInt nmax, count, *sindices, *rindices, i, j, idx, lastidx, *sindices2, *rindices2; 95 const PetscInt *owners = aomems->map->range; 96 MPI_Request *send_waits, *recv_waits, *send_waits2, *recv_waits2; 97 MPI_Status recv_status; 98 PetscMPIInt nindices, source, widx; 99 PetscInt *rbuf, *sbuf; 100 MPI_Status *send_status, *send_status2; 101 102 PetscFunctionBegin; 103 PetscCall(PetscObjectGetComm((PetscObject)ao, &comm)); 104 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 105 PetscCallMPI(MPI_Comm_size(comm, &size)); 106 107 /* first count number of contributors to each processor */ 108 PetscCall(PetscMalloc1(size, &start)); 109 PetscCall(PetscCalloc2(2 * size, &sizes, n, &owner)); 110 111 j = 0; 112 lastidx = -1; 113 for (i = 0; i < n; i++) { 114 if (ia[i] < 0) owner[i] = -1; /* mark negative entries (which are not to be mapped) with a special negative value */ 115 if (ia[i] >= ao->N) owner[i] = -2; /* mark out of range entries with special negative value */ 116 else { 117 /* if indices are NOT locally sorted, need to start search at the beginning */ 118 if (lastidx > (idx = ia[i])) j = 0; 119 lastidx = idx; 120 for (; j < size; j++) { 121 if (idx >= owners[j] && idx < owners[j + 1]) { 122 sizes[2 * j]++; /* num of indices to be sent */ 123 sizes[2 * j + 1] = 1; /* send to proc[j] */ 124 owner[i] = j; 125 break; 126 } 127 } 128 } 129 } 130 sizes[2 * rank] = sizes[2 * rank + 1] = 0; /* do not receive from self! */ 131 nsends = 0; 132 for (i = 0; i < size; i++) nsends += sizes[2 * i + 1]; 133 134 /* inform other processors of number of messages and max length*/ 135 PetscCall(PetscMaxSum(comm, sizes, &nmax, &nreceives)); 136 137 /* allocate arrays */ 138 PetscCall(PetscObjectGetNewTag((PetscObject)ao, &tag1)); 139 PetscCall(PetscObjectGetNewTag((PetscObject)ao, &tag2)); 140 141 PetscCall(PetscMalloc2(nreceives * nmax, &rindices, nreceives, &recv_waits)); 142 PetscCall(PetscMalloc2(nsends * nmax, &rindices2, nsends, &recv_waits2)); 143 144 PetscCall(PetscMalloc3(n, &sindices, nsends, &send_waits, nsends, &send_status)); 145 PetscCall(PetscMalloc3(n, &sindices2, nreceives, &send_waits2, nreceives, &send_status2)); 146 147 /* post 1st receives: receive others requests 148 since we don't know how long each individual message is we 149 allocate the largest needed buffer for each receive. Potentially 150 this is a lot of wasted space. 151 */ 152 for (i = 0, count = 0; i < nreceives; i++) { PetscCallMPI(MPI_Irecv(rindices + nmax * i, nmax, MPIU_INT, MPI_ANY_SOURCE, tag1, comm, recv_waits + count++)); } 153 154 /* do 1st sends: 155 1) starts[i] gives the starting index in svalues for stuff going to 156 the ith processor 157 */ 158 start[0] = 0; 159 for (i = 1; i < size; i++) start[i] = start[i - 1] + sizes[2 * i - 2]; 160 for (i = 0; i < n; i++) { 161 j = owner[i]; 162 if (j == -1) continue; /* do not remap negative entries in ia[] */ 163 else if (j == -2) { /* out of range entries get mapped to -1 */ ia[i] = -1; 164 continue; 165 } else if (j != rank) { 166 sindices[start[j]++] = ia[i]; 167 } else { /* compute my own map */ 168 ia[i] = maploc[ia[i] - owners[rank]]; 169 } 170 } 171 172 start[0] = 0; 173 for (i = 1; i < size; i++) start[i] = start[i - 1] + sizes[2 * i - 2]; 174 for (i = 0, count = 0; i < size; i++) { 175 if (sizes[2 * i + 1]) { 176 /* send my request to others */ 177 PetscCallMPI(MPI_Isend(sindices + start[i], sizes[2 * i], MPIU_INT, i, tag1, comm, send_waits + count)); 178 /* post receive for the answer of my request */ 179 PetscCallMPI(MPI_Irecv(sindices2 + start[i], sizes[2 * i], MPIU_INT, i, tag2, comm, recv_waits2 + count)); 180 count++; 181 } 182 } 183 PetscCheck(nsends == count, comm, PETSC_ERR_SUP, "nsends %" PetscInt_FMT " != count %" PetscInt_FMT, nsends, count); 184 185 /* wait on 1st sends */ 186 if (nsends) { PetscCallMPI(MPI_Waitall(nsends, send_waits, send_status)); } 187 188 /* 1st recvs: other's requests */ 189 for (j = 0; j < nreceives; j++) { 190 PetscCallMPI(MPI_Waitany(nreceives, recv_waits, &widx, &recv_status)); /* idx: index of handle for operation that completed */ 191 PetscCallMPI(MPI_Get_count(&recv_status, MPIU_INT, &nindices)); 192 rbuf = rindices + nmax * widx; /* global index */ 193 source = recv_status.MPI_SOURCE; 194 195 /* compute mapping */ 196 sbuf = rbuf; 197 for (i = 0; i < nindices; i++) sbuf[i] = maploc[rbuf[i] - owners[rank]]; 198 199 /* send mapping back to the sender */ 200 PetscCallMPI(MPI_Isend(sbuf, nindices, MPIU_INT, source, tag2, comm, send_waits2 + widx)); 201 } 202 203 /* wait on 2nd sends */ 204 if (nreceives) { PetscCallMPI(MPI_Waitall(nreceives, send_waits2, send_status2)); } 205 206 /* 2nd recvs: for the answer of my request */ 207 for (j = 0; j < nsends; j++) { 208 PetscCallMPI(MPI_Waitany(nsends, recv_waits2, &widx, &recv_status)); 209 PetscCallMPI(MPI_Get_count(&recv_status, MPIU_INT, &nindices)); 210 source = recv_status.MPI_SOURCE; 211 /* pack output ia[] */ 212 rbuf = sindices2 + start[source]; 213 count = 0; 214 for (i = 0; i < n; i++) { 215 if (source == owner[i]) ia[i] = rbuf[count++]; 216 } 217 } 218 219 /* free arrays */ 220 PetscCall(PetscFree(start)); 221 PetscCall(PetscFree2(sizes, owner)); 222 PetscCall(PetscFree2(rindices, recv_waits)); 223 PetscCall(PetscFree2(rindices2, recv_waits2)); 224 PetscCall(PetscFree3(sindices, send_waits, send_status)); 225 PetscCall(PetscFree3(sindices2, send_waits2, send_status2)); 226 PetscFunctionReturn(0); 227 } 228 229 PetscErrorCode AOPetscToApplication_MemoryScalable(AO ao, PetscInt n, PetscInt *ia) { 230 AO_MemoryScalable *aomems = (AO_MemoryScalable *)ao->data; 231 PetscInt *app_loc = aomems->app_loc; 232 233 PetscFunctionBegin; 234 PetscCall(AOMap_MemoryScalable_private(ao, n, ia, app_loc)); 235 PetscFunctionReturn(0); 236 } 237 238 PetscErrorCode AOApplicationToPetsc_MemoryScalable(AO ao, PetscInt n, PetscInt *ia) { 239 AO_MemoryScalable *aomems = (AO_MemoryScalable *)ao->data; 240 PetscInt *petsc_loc = aomems->petsc_loc; 241 242 PetscFunctionBegin; 243 PetscCall(AOMap_MemoryScalable_private(ao, n, ia, petsc_loc)); 244 PetscFunctionReturn(0); 245 } 246 247 static struct _AOOps AOOps_MemoryScalable = { 248 PetscDesignatedInitializer(view, AOView_MemoryScalable), 249 PetscDesignatedInitializer(destroy, AODestroy_MemoryScalable), 250 PetscDesignatedInitializer(petsctoapplication, AOPetscToApplication_MemoryScalable), 251 PetscDesignatedInitializer(applicationtopetsc, AOApplicationToPetsc_MemoryScalable), 252 }; 253 254 PetscErrorCode AOCreateMemoryScalable_private(MPI_Comm comm, PetscInt napp, const PetscInt from_array[], const PetscInt to_array[], AO ao, PetscInt *aomap_loc) { 255 AO_MemoryScalable *aomems = (AO_MemoryScalable *)ao->data; 256 PetscLayout map = aomems->map; 257 PetscInt n_local = map->n, i, j; 258 PetscMPIInt rank, size, tag; 259 PetscInt *owner, *start, *sizes, nsends, nreceives; 260 PetscInt nmax, count, *sindices, *rindices, idx, lastidx; 261 PetscInt *owners = aomems->map->range; 262 MPI_Request *send_waits, *recv_waits; 263 MPI_Status recv_status; 264 PetscMPIInt nindices, widx; 265 PetscInt *rbuf; 266 PetscInt n = napp, ip, ia; 267 MPI_Status *send_status; 268 269 PetscFunctionBegin; 270 PetscCall(PetscArrayzero(aomap_loc, n_local)); 271 272 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 273 PetscCallMPI(MPI_Comm_size(comm, &size)); 274 275 /* first count number of contributors (of from_array[]) to each processor */ 276 PetscCall(PetscCalloc1(2 * size, &sizes)); 277 PetscCall(PetscMalloc1(n, &owner)); 278 279 j = 0; 280 lastidx = -1; 281 for (i = 0; i < n; i++) { 282 /* if indices are NOT locally sorted, need to start search at the beginning */ 283 if (lastidx > (idx = from_array[i])) j = 0; 284 lastidx = idx; 285 for (; j < size; j++) { 286 if (idx >= owners[j] && idx < owners[j + 1]) { 287 sizes[2 * j] += 2; /* num of indices to be sent - in pairs (ip,ia) */ 288 sizes[2 * j + 1] = 1; /* send to proc[j] */ 289 owner[i] = j; 290 break; 291 } 292 } 293 } 294 sizes[2 * rank] = sizes[2 * rank + 1] = 0; /* do not receive from self! */ 295 nsends = 0; 296 for (i = 0; i < size; i++) nsends += sizes[2 * i + 1]; 297 298 /* inform other processors of number of messages and max length*/ 299 PetscCall(PetscMaxSum(comm, sizes, &nmax, &nreceives)); 300 301 /* allocate arrays */ 302 PetscCall(PetscObjectGetNewTag((PetscObject)ao, &tag)); 303 PetscCall(PetscMalloc2(nreceives * nmax, &rindices, nreceives, &recv_waits)); 304 PetscCall(PetscMalloc3(2 * n, &sindices, nsends, &send_waits, nsends, &send_status)); 305 PetscCall(PetscMalloc1(size, &start)); 306 307 /* post receives: */ 308 for (i = 0; i < nreceives; i++) { PetscCallMPI(MPI_Irecv(rindices + nmax * i, nmax, MPIU_INT, MPI_ANY_SOURCE, tag, comm, recv_waits + i)); } 309 310 /* do sends: 311 1) starts[i] gives the starting index in svalues for stuff going to 312 the ith processor 313 */ 314 start[0] = 0; 315 for (i = 1; i < size; i++) start[i] = start[i - 1] + sizes[2 * i - 2]; 316 for (i = 0; i < n; i++) { 317 j = owner[i]; 318 if (j != rank) { 319 ip = from_array[i]; 320 ia = to_array[i]; 321 sindices[start[j]++] = ip; 322 sindices[start[j]++] = ia; 323 } else { /* compute my own map */ 324 ip = from_array[i] - owners[rank]; 325 ia = to_array[i]; 326 aomap_loc[ip] = ia; 327 } 328 } 329 330 start[0] = 0; 331 for (i = 1; i < size; i++) start[i] = start[i - 1] + sizes[2 * i - 2]; 332 for (i = 0, count = 0; i < size; i++) { 333 if (sizes[2 * i + 1]) { 334 PetscCallMPI(MPI_Isend(sindices + start[i], sizes[2 * i], MPIU_INT, i, tag, comm, send_waits + count)); 335 count++; 336 } 337 } 338 PetscCheck(nsends == count, comm, PETSC_ERR_SUP, "nsends %" PetscInt_FMT " != count %" PetscInt_FMT, nsends, count); 339 340 /* wait on sends */ 341 if (nsends) { PetscCallMPI(MPI_Waitall(nsends, send_waits, send_status)); } 342 343 /* recvs */ 344 count = 0; 345 for (j = nreceives; j > 0; j--) { 346 PetscCallMPI(MPI_Waitany(nreceives, recv_waits, &widx, &recv_status)); 347 PetscCallMPI(MPI_Get_count(&recv_status, MPIU_INT, &nindices)); 348 rbuf = rindices + nmax * widx; /* global index */ 349 350 /* compute local mapping */ 351 for (i = 0; i < nindices; i += 2) { /* pack aomap_loc */ 352 ip = rbuf[i] - owners[rank]; /* local index */ 353 ia = rbuf[i + 1]; 354 aomap_loc[ip] = ia; 355 } 356 count++; 357 } 358 359 PetscCall(PetscFree(start)); 360 PetscCall(PetscFree3(sindices, send_waits, send_status)); 361 PetscCall(PetscFree2(rindices, recv_waits)); 362 PetscCall(PetscFree(owner)); 363 PetscCall(PetscFree(sizes)); 364 PetscFunctionReturn(0); 365 } 366 367 PETSC_EXTERN PetscErrorCode AOCreate_MemoryScalable(AO ao) { 368 IS isapp = ao->isapp, ispetsc = ao->ispetsc; 369 const PetscInt *mypetsc, *myapp; 370 PetscInt napp, n_local, N, i, start, *petsc, *lens, *disp; 371 MPI_Comm comm; 372 AO_MemoryScalable *aomems; 373 PetscLayout map; 374 PetscMPIInt size, rank; 375 376 PetscFunctionBegin; 377 PetscCheck(isapp, PetscObjectComm((PetscObject)ao), PETSC_ERR_ARG_WRONGSTATE, "AOSetIS() must be called before AOSetType()"); 378 /* create special struct aomems */ 379 PetscCall(PetscNewLog(ao, &aomems)); 380 ao->data = (void *)aomems; 381 PetscCall(PetscMemcpy(ao->ops, &AOOps_MemoryScalable, sizeof(struct _AOOps))); 382 PetscCall(PetscObjectChangeTypeName((PetscObject)ao, AOMEMORYSCALABLE)); 383 384 /* transmit all local lengths of isapp to all processors */ 385 PetscCall(PetscObjectGetComm((PetscObject)isapp, &comm)); 386 PetscCallMPI(MPI_Comm_size(comm, &size)); 387 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 388 PetscCall(PetscMalloc2(size, &lens, size, &disp)); 389 PetscCall(ISGetLocalSize(isapp, &napp)); 390 PetscCallMPI(MPI_Allgather(&napp, 1, MPIU_INT, lens, 1, MPIU_INT, comm)); 391 392 N = 0; 393 for (i = 0; i < size; i++) { 394 disp[i] = N; 395 N += lens[i]; 396 } 397 398 /* If ispetsc is 0 then use "natural" numbering */ 399 if (napp) { 400 if (!ispetsc) { 401 start = disp[rank]; 402 PetscCall(PetscMalloc1(napp + 1, &petsc)); 403 for (i = 0; i < napp; i++) petsc[i] = start + i; 404 } else { 405 PetscCall(ISGetIndices(ispetsc, &mypetsc)); 406 petsc = (PetscInt *)mypetsc; 407 } 408 } else { 409 petsc = NULL; 410 } 411 412 /* create a map with global size N - used to determine the local sizes of ao - shall we use local napp instead of N? */ 413 PetscCall(PetscLayoutCreate(comm, &map)); 414 map->bs = 1; 415 map->N = N; 416 PetscCall(PetscLayoutSetUp(map)); 417 418 ao->N = N; 419 ao->n = map->n; 420 aomems->map = map; 421 422 /* create distributed indices app_loc: petsc->app and petsc_loc: app->petsc */ 423 n_local = map->n; 424 PetscCall(PetscCalloc2(n_local, &aomems->app_loc, n_local, &aomems->petsc_loc)); 425 PetscCall(PetscLogObjectMemory((PetscObject)ao, 2 * n_local * sizeof(PetscInt))); 426 PetscCall(ISGetIndices(isapp, &myapp)); 427 428 PetscCall(AOCreateMemoryScalable_private(comm, napp, petsc, myapp, ao, aomems->app_loc)); 429 PetscCall(AOCreateMemoryScalable_private(comm, napp, myapp, petsc, ao, aomems->petsc_loc)); 430 431 PetscCall(ISRestoreIndices(isapp, &myapp)); 432 if (napp) { 433 if (ispetsc) { 434 PetscCall(ISRestoreIndices(ispetsc, &mypetsc)); 435 } else { 436 PetscCall(PetscFree(petsc)); 437 } 438 } 439 PetscCall(PetscFree2(lens, disp)); 440 PetscFunctionReturn(0); 441 } 442 443 /*@C 444 AOCreateMemoryScalable - Creates a memory scalable application ordering using two integer arrays. 445 446 Collective 447 448 Input Parameters: 449 + comm - MPI communicator that is to share AO 450 . napp - size of integer arrays 451 . myapp - integer array that defines an ordering 452 - mypetsc - integer array that defines another ordering (may be NULL to 453 indicate the natural ordering, that is 0,1,2,3,...) 454 455 Output Parameter: 456 . aoout - the new application ordering 457 458 Level: beginner 459 460 Notes: 461 The arrays myapp and mypetsc must contain the all the integers 0 to napp-1 with no duplicates; that is there cannot be any "holes" 462 in the indices. Use AOCreateMapping() or AOCreateMappingIS() if you wish to have "holes" in the indices. 463 Comparing with AOCreateBasic(), this routine trades memory with message communication. 464 465 .seealso: `AOCreateMemoryScalableIS()`, `AODestroy()`, `AOPetscToApplication()`, `AOApplicationToPetsc()` 466 @*/ 467 PetscErrorCode AOCreateMemoryScalable(MPI_Comm comm, PetscInt napp, const PetscInt myapp[], const PetscInt mypetsc[], AO *aoout) { 468 IS isapp, ispetsc; 469 const PetscInt *app = myapp, *petsc = mypetsc; 470 471 PetscFunctionBegin; 472 PetscCall(ISCreateGeneral(comm, napp, app, PETSC_USE_POINTER, &isapp)); 473 if (mypetsc) { 474 PetscCall(ISCreateGeneral(comm, napp, petsc, PETSC_USE_POINTER, &ispetsc)); 475 } else { 476 ispetsc = NULL; 477 } 478 PetscCall(AOCreateMemoryScalableIS(isapp, ispetsc, aoout)); 479 PetscCall(ISDestroy(&isapp)); 480 if (mypetsc) { PetscCall(ISDestroy(&ispetsc)); } 481 PetscFunctionReturn(0); 482 } 483 484 /*@C 485 AOCreateMemoryScalableIS - Creates a memory scalable application ordering using two index sets. 486 487 Collective on IS 488 489 Input Parameters: 490 + isapp - index set that defines an ordering 491 - ispetsc - index set that defines another ordering (may be NULL to use the 492 natural ordering) 493 494 Output Parameter: 495 . aoout - the new application ordering 496 497 Level: beginner 498 499 Notes: 500 The index sets isapp and ispetsc must contain the all the integers 0 to napp-1 (where napp is the length of the index sets) with no duplicates; 501 that is there cannot be any "holes". 502 Comparing with AOCreateBasicIS(), this routine trades memory with message communication. 503 .seealso: `AOCreateMemoryScalable()`, `AODestroy()` 504 @*/ 505 PetscErrorCode AOCreateMemoryScalableIS(IS isapp, IS ispetsc, AO *aoout) { 506 MPI_Comm comm; 507 AO ao; 508 509 PetscFunctionBegin; 510 PetscCall(PetscObjectGetComm((PetscObject)isapp, &comm)); 511 PetscCall(AOCreate(comm, &ao)); 512 PetscCall(AOSetIS(ao, isapp, ispetsc)); 513 PetscCall(AOSetType(ao, AOMEMORYSCALABLE)); 514 PetscCall(AOViewFromOptions(ao, NULL, "-ao_view")); 515 *aoout = ao; 516 PetscFunctionReturn(0); 517 } 518