Actual source code: aomemscalable.c

  1: /*
  2:     The memory scalable AO application ordering routines. These store the
  3:   orderings on each process for that process' range of values, this is more memory-efficient than `AOBASIC`
  4: */

  6: #include <../src/vec/is/ao/aoimpl.h>

  8: typedef struct {
  9:   PetscInt   *app_loc;   /* app_loc[i] is the partner for the ith local PETSc slot */
 10:   PetscInt   *petsc_loc; /* petsc_loc[j] is the partner for the jth local app slot */
 11:   PetscLayout map;       /* determines the local sizes of ao */
 12: } AO_MemoryScalable;

 14: /*
 15:        All processors ship the data to process 0 to be printed; note that this is not scalable because
 16:        process 0 allocates space for all the orderings entry across all the processes
 17: */
 18: static PetscErrorCode AOView_MemoryScalable(AO ao, PetscViewer viewer)
 19: {
 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, j;
 26:   MPI_Status         status;

 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);

 32:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)ao), &rank));
 33:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)ao), &size));

 35:   PetscCall(PetscObjectGetNewTag((PetscObject)ao, &tag_app));
 36:   PetscCall(PetscObjectGetNewTag((PetscObject)ao, &tag_petsc));

 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"));

 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 (PetscInt 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]));

 48:     /* recv and print off-processor's AO */
 49:     for (PetscMPIInt 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(MPIU_Recv(app_loc, len, MPIU_INT, i, tag_app, PetscObjectComm((PetscObject)ao), &status));
 54:       PetscCallMPI(MPIU_Recv(petsc_loc, len, MPIU_INT, i, tag_petsc, PetscObjectComm((PetscObject)ao), &status));
 55:       PetscCall(PetscViewerASCIIPrintf(viewer, "Process [%d]\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));

 60:   } else {
 61:     /* send values */
 62:     PetscCallMPI(MPIU_Send((void *)aomems->app_loc, map->n, MPIU_INT, 0, tag_app, PetscObjectComm((PetscObject)ao)));
 63:     PetscCallMPI(MPIU_Send((void *)aomems->petsc_loc, map->n, MPIU_INT, 0, tag_petsc, PetscObjectComm((PetscObject)ao)));
 64:   }
 65:   PetscCall(PetscViewerFlush(viewer));
 66:   PetscFunctionReturn(PETSC_SUCCESS);
 67: }

 69: static PetscErrorCode AODestroy_MemoryScalable(AO ao)
 70: {
 71:   AO_MemoryScalable *aomems = (AO_MemoryScalable *)ao->data;

 73:   PetscFunctionBegin;
 74:   PetscCall(PetscFree2(aomems->app_loc, aomems->petsc_loc));
 75:   PetscCall(PetscLayoutDestroy(&aomems->map));
 76:   PetscCall(PetscFree(aomems));
 77:   PetscFunctionReturn(PETSC_SUCCESS);
 78: }

 80: /*
 81:    Input Parameters:
 82: +   ao - the application ordering context
 83: .   n  - the number of integers in ia[]
 84: .   ia - the integers; these are replaced with their mapped value
 85: -   maploc - app_loc or petsc_loc in struct "AO_MemoryScalable"

 87:    Output Parameter:
 88: .   ia - the mapped interges
 89:  */
 90: static PetscErrorCode AOMap_MemoryScalable_private(AO ao, PetscInt n, PetscInt *ia, const PetscInt *maploc)
 91: {
 92:   AO_MemoryScalable *aomems = (AO_MemoryScalable *)ao->data;
 93:   MPI_Comm           comm;
 94:   PetscMPIInt        rank, size, tag1, tag2, count;
 95:   PetscMPIInt       *owner, j;
 96:   PetscInt           nmax, *sindices, *rindices, idx, lastidx, *sindices2, *rindices2, *sizes, *start;
 97:   PetscMPIInt        nreceives, nsends;
 98:   const PetscInt    *owners = aomems->map->range;
 99:   MPI_Request       *send_waits, *recv_waits, *send_waits2, *recv_waits2;
100:   MPI_Status         recv_status;
101:   PetscMPIInt        source, widx;
102:   PetscCount         nindices;
103:   PetscInt          *rbuf, *sbuf, inreceives;
104:   MPI_Status        *send_status, *send_status2;

106:   PetscFunctionBegin;
107:   PetscCall(PetscObjectGetComm((PetscObject)ao, &comm));
108:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
109:   PetscCallMPI(MPI_Comm_size(comm, &size));

111:   /*  first count number of contributors to each processor */
112:   PetscCall(PetscMalloc1(size, &start));
113:   PetscCall(PetscCalloc2(2 * size, &sizes, n, &owner));

115:   j       = 0;
116:   lastidx = -1;
117:   for (PetscInt i = 0; i < n; i++) {
118:     if (ia[i] < 0) owner[i] = -1;      /* mark negative entries (which are not to be mapped) with a special negative value */
119:     if (ia[i] >= ao->N) owner[i] = -2; /* mark out of range entries with special negative value */
120:     else {
121:       /* if indices are NOT locally sorted, need to start search at the beginning */
122:       if (lastidx > (idx = ia[i])) j = 0;
123:       lastidx = idx;
124:       for (; j < size; j++) {
125:         if (idx >= owners[j] && idx < owners[j + 1]) {
126:           sizes[2 * j]++;       /* num of indices to be sent */
127:           sizes[2 * j + 1] = 1; /* send to proc[j] */
128:           owner[i]         = j;
129:           break;
130:         }
131:       }
132:     }
133:   }
134:   sizes[2 * rank] = sizes[2 * rank + 1] = 0; /* do not receive from self! */
135:   nsends                                = 0;
136:   for (PetscMPIInt i = 0; i < size; i++) nsends += sizes[2 * i + 1];

138:   /* inform other processors of number of messages and max length*/
139:   PetscCall(PetscMaxSum(comm, sizes, &nmax, &inreceives));
140:   PetscCall(PetscMPIIntCast(inreceives, &nreceives));

142:   /* allocate arrays */
143:   PetscCall(PetscObjectGetNewTag((PetscObject)ao, &tag1));
144:   PetscCall(PetscObjectGetNewTag((PetscObject)ao, &tag2));

146:   PetscCall(PetscMalloc2(nreceives * nmax, &rindices, nreceives, &recv_waits));
147:   PetscCall(PetscMalloc2(nsends * nmax, &rindices2, nsends, &recv_waits2));

149:   PetscCall(PetscMalloc3(n, &sindices, nsends, &send_waits, nsends, &send_status));
150:   PetscCall(PetscMalloc3(n, &sindices2, nreceives, &send_waits2, nreceives, &send_status2));

152:   /* post 1st receives: receive others requests
153:      since we don't know how long each individual message is we
154:      allocate the largest needed buffer for each receive. Potentially
155:      this is a lot of wasted space.
156:   */
157:   for (PetscMPIInt i = 0; i < nreceives; i++) PetscCallMPI(MPIU_Irecv(rindices + nmax * i, nmax, MPIU_INT, MPI_ANY_SOURCE, tag1, comm, recv_waits + i));

159:   /* do 1st sends:
160:       1) starts[i] gives the starting index in svalues for stuff going to
161:          the ith processor
162:   */
163:   start[0] = 0;
164:   for (PetscMPIInt i = 1; i < size; i++) start[i] = start[i - 1] + sizes[2 * i - 2];
165:   for (PetscInt i = 0; i < n; i++) {
166:     j = owner[i];
167:     if (j == -1) continue; /* do not remap negative entries in ia[] */
168:     else if (j == -2) { /* out of range entries get mapped to -1 */ ia[i] = -1;
169:       continue;
170:     } else if (j != rank) {
171:       sindices[start[j]++] = ia[i];
172:     } else { /* compute my own map */
173:       ia[i] = maploc[ia[i] - owners[rank]];
174:     }
175:   }

177:   start[0] = 0;
178:   for (PetscMPIInt i = 1; i < size; i++) start[i] = start[i - 1] + sizes[2 * i - 2];
179:   count = 0;
180:   for (PetscMPIInt i = 0; i < size; i++) {
181:     if (sizes[2 * i + 1]) {
182:       /* send my request to others */
183:       PetscCallMPI(MPIU_Isend(sindices + start[i], sizes[2 * i], MPIU_INT, i, tag1, comm, send_waits + count));
184:       /* post receive for the answer of my request */
185:       PetscCallMPI(MPIU_Irecv(sindices2 + start[i], sizes[2 * i], MPIU_INT, i, tag2, comm, recv_waits2 + count));
186:       count++;
187:     }
188:   }
189:   PetscCheck(nsends == count, comm, PETSC_ERR_SUP, "nsends %d != count %d", nsends, count);

191:   /* wait on 1st sends */
192:   if (nsends) PetscCallMPI(MPI_Waitall((PetscMPIInt)nsends, send_waits, send_status));

194:   /* 1st recvs: other's requests */
195:   for (j = 0; j < nreceives; j++) {
196:     PetscCallMPI(MPI_Waitany((PetscMPIInt)nreceives, recv_waits, &widx, &recv_status)); /* idx: index of handle for operation that completed */
197:     PetscCallMPI(MPIU_Get_count(&recv_status, MPIU_INT, &nindices));
198:     rbuf   = rindices + nmax * widx; /* global index */
199:     source = recv_status.MPI_SOURCE;

201:     /* compute mapping */
202:     sbuf = rbuf;
203:     for (PetscCount i = 0; i < nindices; i++) sbuf[i] = maploc[rbuf[i] - owners[rank]];

205:     /* send mapping back to the sender */
206:     PetscCallMPI(MPIU_Isend(sbuf, nindices, MPIU_INT, source, tag2, comm, send_waits2 + widx));
207:   }

209:   /* wait on 2nd sends */
210:   if (nreceives) PetscCallMPI(MPI_Waitall((PetscMPIInt)nreceives, send_waits2, send_status2));

212:   /* 2nd recvs: for the answer of my request */
213:   for (j = 0; j < nsends; j++) {
214:     PetscCallMPI(MPI_Waitany((PetscMPIInt)nsends, recv_waits2, &widx, &recv_status));
215:     PetscCallMPI(MPIU_Get_count(&recv_status, MPIU_INT, &nindices));
216:     source = recv_status.MPI_SOURCE;
217:     /* pack output ia[] */
218:     rbuf  = sindices2 + start[source];
219:     count = 0;
220:     for (PetscCount i = 0; i < n; i++) {
221:       if (source == owner[i]) ia[i] = rbuf[count++];
222:     }
223:   }

225:   /* free arrays */
226:   PetscCall(PetscFree(start));
227:   PetscCall(PetscFree2(sizes, owner));
228:   PetscCall(PetscFree2(rindices, recv_waits));
229:   PetscCall(PetscFree2(rindices2, recv_waits2));
230:   PetscCall(PetscFree3(sindices, send_waits, send_status));
231:   PetscCall(PetscFree3(sindices2, send_waits2, send_status2));
232:   PetscFunctionReturn(PETSC_SUCCESS);
233: }

235: static PetscErrorCode AOPetscToApplication_MemoryScalable(AO ao, PetscInt n, PetscInt *ia)
236: {
237:   AO_MemoryScalable *aomems  = (AO_MemoryScalable *)ao->data;
238:   PetscInt          *app_loc = aomems->app_loc;

240:   PetscFunctionBegin;
241:   PetscCall(AOMap_MemoryScalable_private(ao, n, ia, app_loc));
242:   PetscFunctionReturn(PETSC_SUCCESS);
243: }

245: static PetscErrorCode AOApplicationToPetsc_MemoryScalable(AO ao, PetscInt n, PetscInt *ia)
246: {
247:   AO_MemoryScalable *aomems    = (AO_MemoryScalable *)ao->data;
248:   PetscInt          *petsc_loc = aomems->petsc_loc;

250:   PetscFunctionBegin;
251:   PetscCall(AOMap_MemoryScalable_private(ao, n, ia, petsc_loc));
252:   PetscFunctionReturn(PETSC_SUCCESS);
253: }

255: static const struct _AOOps AOOps_MemoryScalable = {
256:   PetscDesignatedInitializer(view, AOView_MemoryScalable),
257:   PetscDesignatedInitializer(destroy, AODestroy_MemoryScalable),
258:   PetscDesignatedInitializer(petsctoapplication, AOPetscToApplication_MemoryScalable),
259:   PetscDesignatedInitializer(applicationtopetsc, AOApplicationToPetsc_MemoryScalable),
260: };

262: static PetscErrorCode AOCreateMemoryScalable_private(MPI_Comm comm, PetscInt napp, const PetscInt from_array[], const PetscInt to_array[], AO ao, PetscInt *aomap_loc)
263: {
264:   AO_MemoryScalable *aomems  = (AO_MemoryScalable *)ao->data;
265:   PetscLayout        map     = aomems->map;
266:   PetscInt           n_local = map->n, i, j;
267:   PetscMPIInt        rank, size, tag;
268:   PetscInt          *owner, *start, *sizes, nsends, nreceives;
269:   PetscInt           nmax, count, *sindices, *rindices, idx, lastidx;
270:   PetscInt          *owners = aomems->map->range;
271:   MPI_Request       *send_waits, *recv_waits;
272:   MPI_Status         recv_status;
273:   PetscMPIInt        widx;
274:   PetscInt          *rbuf;
275:   PetscInt           n = napp, ip, ia;
276:   MPI_Status        *send_status;
277:   PetscCount         nindices;

279:   PetscFunctionBegin;
280:   PetscCall(PetscArrayzero(aomap_loc, n_local));

282:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
283:   PetscCallMPI(MPI_Comm_size(comm, &size));

285:   /*  first count number of contributors (of from_array[]) to each processor */
286:   PetscCall(PetscCalloc1(2 * size, &sizes));
287:   PetscCall(PetscMalloc1(n, &owner));

289:   j       = 0;
290:   lastidx = -1;
291:   for (i = 0; i < n; i++) {
292:     /* if indices are NOT locally sorted, need to start search at the beginning */
293:     if (lastidx > (idx = from_array[i])) j = 0;
294:     lastidx = idx;
295:     for (; j < size; j++) {
296:       if (idx >= owners[j] && idx < owners[j + 1]) {
297:         sizes[2 * j] += 2;    /* num of indices to be sent - in pairs (ip,ia) */
298:         sizes[2 * j + 1] = 1; /* send to proc[j] */
299:         owner[i]         = j;
300:         break;
301:       }
302:     }
303:   }
304:   sizes[2 * rank] = sizes[2 * rank + 1] = 0; /* do not receive from self! */
305:   nsends                                = 0;
306:   for (i = 0; i < size; i++) nsends += sizes[2 * i + 1];

308:   /* inform other processors of number of messages and max length*/
309:   PetscCall(PetscMaxSum(comm, sizes, &nmax, &nreceives));

311:   /* allocate arrays */
312:   PetscCall(PetscObjectGetNewTag((PetscObject)ao, &tag));
313:   PetscCall(PetscMalloc2(nreceives * nmax, &rindices, nreceives, &recv_waits));
314:   PetscCall(PetscMalloc3(2 * n, &sindices, nsends, &send_waits, nsends, &send_status));
315:   PetscCall(PetscMalloc1(size, &start));

317:   /* post receives: */
318:   for (i = 0; i < nreceives; i++) PetscCallMPI(MPIU_Irecv(rindices + nmax * i, nmax, MPIU_INT, MPI_ANY_SOURCE, tag, comm, recv_waits + i));

320:   /* do sends:
321:       1) starts[i] gives the starting index in svalues for stuff going to
322:          the ith processor
323:   */
324:   start[0] = 0;
325:   for (i = 1; i < size; i++) start[i] = start[i - 1] + sizes[2 * i - 2];
326:   for (i = 0; i < n; i++) {
327:     j = owner[i];
328:     if (j != rank) {
329:       ip                   = from_array[i];
330:       ia                   = to_array[i];
331:       sindices[start[j]++] = ip;
332:       sindices[start[j]++] = ia;
333:     } else { /* compute my own map */
334:       ip            = from_array[i] - owners[rank];
335:       ia            = to_array[i];
336:       aomap_loc[ip] = ia;
337:     }
338:   }

340:   start[0] = 0;
341:   for (PetscMPIInt i = 1; i < size; i++) start[i] = start[i - 1] + sizes[2 * i - 2];
342:   count = 0;
343:   for (PetscMPIInt i = 0; i < size; i++) {
344:     if (sizes[2 * i + 1]) {
345:       PetscCallMPI(MPIU_Isend(sindices + start[i], sizes[2 * i], MPIU_INT, i, tag, comm, send_waits + count));
346:       count++;
347:     }
348:   }
349:   PetscCheck(nsends == count, comm, PETSC_ERR_SUP, "nsends %" PetscInt_FMT " != count %" PetscInt_FMT, nsends, count);

351:   /* wait on sends */
352:   if (nsends) PetscCallMPI(MPI_Waitall((PetscMPIInt)nsends, send_waits, send_status));

354:   /* recvs */
355:   count = 0;
356:   for (j = nreceives; j > 0; j--) {
357:     PetscCallMPI(MPI_Waitany((PetscMPIInt)nreceives, recv_waits, &widx, &recv_status));
358:     PetscCallMPI(MPIU_Get_count(&recv_status, MPIU_INT, &nindices));
359:     rbuf = rindices + nmax * widx; /* global index */

361:     /* compute local mapping */
362:     for (i = 0; i < nindices; i += 2) {       /* pack aomap_loc */
363:       ip            = rbuf[i] - owners[rank]; /* local index */
364:       ia            = rbuf[i + 1];
365:       aomap_loc[ip] = ia;
366:     }
367:     count++;
368:   }

370:   PetscCall(PetscFree(start));
371:   PetscCall(PetscFree3(sindices, send_waits, send_status));
372:   PetscCall(PetscFree2(rindices, recv_waits));
373:   PetscCall(PetscFree(owner));
374:   PetscCall(PetscFree(sizes));
375:   PetscFunctionReturn(PETSC_SUCCESS);
376: }

378: PETSC_INTERN PetscErrorCode AOCreate_MemoryScalable(AO ao)
379: {
380:   IS                 isapp = ao->isapp, ispetsc = ao->ispetsc;
381:   const PetscInt    *mypetsc, *myapp;
382:   PetscInt           napp, n_local, N, i, start, *petsc, *lens, *disp;
383:   MPI_Comm           comm;
384:   AO_MemoryScalable *aomems;
385:   PetscLayout        map;
386:   PetscMPIInt        size, rank;

388:   PetscFunctionBegin;
389:   PetscCheck(isapp, PetscObjectComm((PetscObject)ao), PETSC_ERR_ARG_WRONGSTATE, "AOSetIS() must be called before AOSetType()");
390:   /* create special struct aomems */
391:   PetscCall(PetscNew(&aomems));
392:   ao->data   = (void *)aomems;
393:   ao->ops[0] = AOOps_MemoryScalable;
394:   PetscCall(PetscObjectChangeTypeName((PetscObject)ao, AOMEMORYSCALABLE));

396:   /* transmit all local lengths of isapp to all processors */
397:   PetscCall(PetscObjectGetComm((PetscObject)isapp, &comm));
398:   PetscCallMPI(MPI_Comm_size(comm, &size));
399:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
400:   PetscCall(PetscMalloc2(size, &lens, size, &disp));
401:   PetscCall(ISGetLocalSize(isapp, &napp));
402:   PetscCallMPI(MPI_Allgather(&napp, 1, MPIU_INT, lens, 1, MPIU_INT, comm));

404:   N = 0;
405:   for (i = 0; i < size; i++) {
406:     disp[i] = N;
407:     N += lens[i];
408:   }

410:   /* If ispetsc is 0 then use "natural" numbering */
411:   if (napp) {
412:     if (!ispetsc) {
413:       start = disp[rank];
414:       PetscCall(PetscMalloc1(napp + 1, &petsc));
415:       for (i = 0; i < napp; i++) petsc[i] = start + i;
416:     } else {
417:       PetscCall(ISGetIndices(ispetsc, &mypetsc));
418:       petsc = (PetscInt *)mypetsc;
419:     }
420:   } else {
421:     petsc = NULL;
422:   }

424:   /* create a map with global size N - used to determine the local sizes of ao - shall we use local napp instead of N? */
425:   PetscCall(PetscLayoutCreate(comm, &map));
426:   map->bs = 1;
427:   map->N  = N;
428:   PetscCall(PetscLayoutSetUp(map));

430:   ao->N       = N;
431:   ao->n       = map->n;
432:   aomems->map = map;

434:   /* create distributed indices app_loc: petsc->app and petsc_loc: app->petsc */
435:   n_local = map->n;
436:   PetscCall(PetscCalloc2(n_local, &aomems->app_loc, n_local, &aomems->petsc_loc));
437:   PetscCall(ISGetIndices(isapp, &myapp));

439:   PetscCall(AOCreateMemoryScalable_private(comm, napp, petsc, myapp, ao, aomems->app_loc));
440:   PetscCall(AOCreateMemoryScalable_private(comm, napp, myapp, petsc, ao, aomems->petsc_loc));

442:   PetscCall(ISRestoreIndices(isapp, &myapp));
443:   if (napp) {
444:     if (ispetsc) {
445:       PetscCall(ISRestoreIndices(ispetsc, &mypetsc));
446:     } else {
447:       PetscCall(PetscFree(petsc));
448:     }
449:   }
450:   PetscCall(PetscFree2(lens, disp));
451:   PetscFunctionReturn(PETSC_SUCCESS);
452: }

454: /*@
455:   AOCreateMemoryScalable - Creates a memory scalable application ordering using two integer arrays.

457:   Collective

459:   Input Parameters:
460: + comm    - MPI communicator that is to share the `AO`
461: . napp    - size of `myapp` and `mypetsc`
462: . myapp   - integer array that defines an ordering
463: - mypetsc - integer array that defines another ordering (may be `NULL` to indicate the natural ordering, that is 0,1,2,3,...)

465:   Output Parameter:
466: . aoout - the new application ordering

468:   Level: beginner

470:   Note:
471:   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"
472:   in the indices. Use `AOCreateMapping()` or `AOCreateMappingIS()` if you wish to have "holes" in the indices.
473:   Comparing with `AOCreateBasic()`, this routine trades memory with message communication.

475: .seealso: [](sec_ao), [](sec_scatter), `AO`, `AOCreateMemoryScalableIS()`, `AODestroy()`, `AOPetscToApplication()`, `AOApplicationToPetsc()`
476: @*/
477: PetscErrorCode AOCreateMemoryScalable(MPI_Comm comm, PetscInt napp, const PetscInt myapp[], const PetscInt mypetsc[], AO *aoout)
478: {
479:   IS              isapp, ispetsc;
480:   const PetscInt *app = myapp, *petsc = mypetsc;

482:   PetscFunctionBegin;
483:   PetscCall(ISCreateGeneral(comm, napp, app, PETSC_USE_POINTER, &isapp));
484:   if (mypetsc) {
485:     PetscCall(ISCreateGeneral(comm, napp, petsc, PETSC_USE_POINTER, &ispetsc));
486:   } else {
487:     ispetsc = NULL;
488:   }
489:   PetscCall(AOCreateMemoryScalableIS(isapp, ispetsc, aoout));
490:   PetscCall(ISDestroy(&isapp));
491:   if (mypetsc) PetscCall(ISDestroy(&ispetsc));
492:   PetscFunctionReturn(PETSC_SUCCESS);
493: }

495: /*@
496:   AOCreateMemoryScalableIS - Creates a memory scalable application ordering using two index sets.

498:   Collective

500:   Input Parameters:
501: + isapp   - index set that defines an ordering
502: - ispetsc - index set that defines another ordering (may be `NULL` to use the natural ordering)

504:   Output Parameter:
505: . aoout - the new application ordering

507:   Level: beginner

509:   Notes:
510:   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;
511:   that is there cannot be any "holes".

513:   Comparing with `AOCreateBasicIS()`, this routine trades memory with message communication.

515: .seealso: [](sec_ao), [](sec_scatter), `AO`, `AOCreateBasicIS()`, `AOCreateMemoryScalable()`, `AODestroy()`
516: @*/
517: PetscErrorCode AOCreateMemoryScalableIS(IS isapp, IS ispetsc, AO *aoout)
518: {
519:   MPI_Comm comm;
520:   AO       ao;

522:   PetscFunctionBegin;
523:   PetscCall(PetscObjectGetComm((PetscObject)isapp, &comm));
524:   PetscCall(AOCreate(comm, &ao));
525:   PetscCall(AOSetIS(ao, isapp, ispetsc));
526:   PetscCall(AOSetType(ao, AOMEMORYSCALABLE));
527:   PetscCall(AOViewFromOptions(ao, NULL, "-ao_view"));
528:   *aoout = ao;
529:   PetscFunctionReturn(PETSC_SUCCESS);
530: }