Actual source code: pmap.c


  2: /*
  3:    This file contains routines for basic map object implementation.
  4: */

  6: #include <petsc/private/isimpl.h>

  8: /*@
  9:   PetscLayoutCreate - Allocates PetscLayout space and sets the PetscLayout contents to the default.

 11:   Collective

 13:   Input Parameters:
 14: . comm - the MPI communicator

 16:   Output Parameters:
 17: . map - the new PetscLayout

 19:   Level: advanced

 21:   Notes:
 22:   Typical calling sequence
 23: .vb
 24:        PetscLayoutCreate(MPI_Comm,PetscLayout *);
 25:        PetscLayoutSetBlockSize(PetscLayout,bs);
 26:        PetscLayoutSetSize(PetscLayout,N); // or PetscLayoutSetLocalSize(PetscLayout,n);
 27:        PetscLayoutSetUp(PetscLayout);
 28: .ve
 29:   Alternatively,
 30: .vb
 31:       PetscLayoutCreateFromSizes(comm,n,N,bs,&layout);
 32: .ve

 34:   Optionally use any of the following
 35: .vb
 36:   PetscLayoutGetSize(PetscLayout,PetscInt *);
 37:   PetscLayoutGetLocalSize(PetscLayout,PetscInt *);
 38:   PetscLayoutGetRange(PetscLayout,PetscInt *rstart,PetscInt *rend);
 39:   PetscLayoutGetRanges(PetscLayout,const PetscInt *range[]);
 40:   PetscLayoutDestroy(PetscLayout*);
 41: .ve

 43:   The PetscLayout object and methods are intended to be used in the PETSc Vec and Mat implementations; it is often not needed in
 44:   user codes unless you really gain something in their use.

 46: .seealso: PetscLayoutSetLocalSize(), PetscLayoutSetSize(), PetscLayoutGetSize(), PetscLayoutGetLocalSize(), PetscLayout, PetscLayoutDestroy(),
 47:           PetscLayoutGetRange(), PetscLayoutGetRanges(), PetscLayoutSetBlockSize(), PetscLayoutGetBlockSize(), PetscLayoutSetUp(),
 48:           PetscLayoutCreateFromSizes()

 50: @*/
 51: PetscErrorCode PetscLayoutCreate(MPI_Comm comm,PetscLayout *map)
 52: {
 53:   PetscNew(map);
 54:   MPI_Comm_size(comm, &(*map)->size);
 55:   (*map)->comm        = comm;
 56:   (*map)->bs          = -1;
 57:   (*map)->n           = -1;
 58:   (*map)->N           = -1;
 59:   (*map)->range       = NULL;
 60:   (*map)->range_alloc = PETSC_TRUE;
 61:   (*map)->rstart      = 0;
 62:   (*map)->rend        = 0;
 63:   (*map)->setupcalled = PETSC_FALSE;
 64:   (*map)->oldn        = -1;
 65:   (*map)->oldN        = -1;
 66:   (*map)->oldbs       = -1;
 67:   return 0;
 68: }

 70: /*@
 71:   PetscLayoutCreateFromSizes - Allocates PetscLayout space, sets the layout sizes, and sets the layout up.

 73:   Collective

 75:   Input Parameters:
 76: + comm  - the MPI communicator
 77: . n     - the local size (or PETSC_DECIDE)
 78: . N     - the global size (or PETSC_DECIDE)
 79: - bs    - the block size (or PETSC_DECIDE)

 81:   Output Parameters:
 82: . map - the new PetscLayout

 84:   Level: advanced

 86:   Notes:
 87: $ PetscLayoutCreateFromSizes(comm,n,N,bs,&layout);
 88:   is a shorthand for
 89: .vb
 90:   PetscLayoutCreate(comm,&layout);
 91:   PetscLayoutSetLocalSize(layout,n);
 92:   PetscLayoutSetSize(layout,N);
 93:   PetscLayoutSetBlockSize(layout,bs);
 94:   PetscLayoutSetUp(layout);
 95: .ve

 97: .seealso: PetscLayoutCreate(), PetscLayoutSetLocalSize(), PetscLayoutSetSize(), PetscLayoutGetSize(), PetscLayoutGetLocalSize(), PetscLayout, PetscLayoutDestroy(),
 98:           PetscLayoutGetRange(), PetscLayoutGetRanges(), PetscLayoutSetBlockSize(), PetscLayoutGetBlockSize(), PetscLayoutSetUp(), PetscLayoutCreateFromRanges()

100: @*/
101: PetscErrorCode PetscLayoutCreateFromSizes(MPI_Comm comm,PetscInt n,PetscInt N,PetscInt bs,PetscLayout *map)
102: {
103:   PetscLayoutCreate(comm, map);
104:   PetscLayoutSetLocalSize(*map, n);
105:   PetscLayoutSetSize(*map, N);
106:   PetscLayoutSetBlockSize(*map, bs);
107:   PetscLayoutSetUp(*map);
108:   return 0;
109: }

111: /*@
112:   PetscLayoutDestroy - Frees a map object and frees its range if that exists.

114:   Collective

116:   Input Parameters:
117: . map - the PetscLayout

119:   Level: developer

121:   Note:
122:   The PetscLayout object and methods are intended to be used in the PETSc Vec and Mat implementations; it is
123:   recommended they not be used in user codes unless you really gain something in their use.

125: .seealso: PetscLayoutSetLocalSize(), PetscLayoutSetSize(), PetscLayoutGetSize(), PetscLayoutGetLocalSize(), PetscLayout, PetscLayoutCreate(),
126:           PetscLayoutGetRange(), PetscLayoutGetRanges(), PetscLayoutSetBlockSize(), PetscLayoutGetBlockSize(), PetscLayoutSetUp()

128: @*/
129: PetscErrorCode PetscLayoutDestroy(PetscLayout *map)
130: {
131:   if (!*map) return 0;
132:   if (!(*map)->refcnt--) {
133:     if ((*map)->range_alloc) PetscFree((*map)->range);
134:     ISLocalToGlobalMappingDestroy(&(*map)->mapping);
135:     PetscFree((*map));
136:   }
137:   *map = NULL;
138:   return 0;
139: }

141: /*@
142:   PetscLayoutCreateFromRanges - Creates a new PetscLayout with the given ownership ranges and sets it up.

144:   Collective

146:   Input Parameters:
147: + comm  - the MPI communicator
148: . range - the array of ownership ranges for each rank with length commsize+1
149: . mode  - the copy mode for range
150: - bs    - the block size (or PETSC_DECIDE)

152:   Output Parameters:
153: . newmap - the new PetscLayout

155:   Level: developer

157: .seealso: PetscLayoutCreate(), PetscLayoutSetLocalSize(), PetscLayoutSetSize(), PetscLayoutGetSize(), PetscLayoutGetLocalSize(), PetscLayout, PetscLayoutDestroy(),
158:           PetscLayoutGetRange(), PetscLayoutGetRanges(), PetscLayoutSetBlockSize(), PetscLayoutGetBlockSize(), PetscLayoutSetUp(), PetscLayoutCreateFromSizes()

160: @*/
161: PetscErrorCode PetscLayoutCreateFromRanges(MPI_Comm comm,const PetscInt range[],PetscCopyMode mode,PetscInt bs,PetscLayout *newmap)
162: {
163:   PetscLayout    map;
164:   PetscMPIInt    rank;

166:   MPI_Comm_rank(comm, &rank);
167:   PetscLayoutCreate(comm, &map);
168:   PetscLayoutSetBlockSize(map, bs);
169:   switch (mode) {
170:     case PETSC_COPY_VALUES:
171:       PetscMalloc1(map->size+1, &map->range);
172:       PetscArraycpy(map->range, range, map->size+1);
173:       break;
174:     case PETSC_USE_POINTER:
175:       map->range_alloc = PETSC_FALSE;
176:     default:
177:       map->range = (PetscInt*) range;
178:       break;
179:   }
180:   map->rstart = map->range[rank];
181:   map->rend   = map->range[rank+1];
182:   map->n      = map->rend - map->rstart;
183:   map->N      = map->range[map->size];
184:   if (PetscDefined(USE_DEBUG)) {  /* just check that n, N and bs are consistent */
185:     PetscInt tmp;
186:     MPIU_Allreduce(&map->n,&tmp,1,MPIU_INT,MPI_SUM,map->comm);
188:     if (map->bs > 1) {
190:     }
191:     if (map->bs > 1) {
193:     }
194:   }
195:   /* lock the layout */
196:   map->setupcalled = PETSC_TRUE;
197:   map->oldn = map->n;
198:   map->oldN = map->N;
199:   map->oldbs = map->bs;
200:   *newmap = map;
201:   return 0;
202: }

204: /*@
205:   PetscLayoutSetUp - given a map where you have set either the global or local
206:                      size sets up the map so that it may be used.

208:   Collective

210:   Input Parameters:
211: . map - pointer to the map

213:   Level: developer

215:   Notes:
216:     Typical calling sequence
217: $ PetscLayoutCreate(MPI_Comm,PetscLayout *);
218: $ PetscLayoutSetBlockSize(PetscLayout,1);
219: $ PetscLayoutSetSize(PetscLayout,n) or PetscLayoutSetLocalSize(PetscLayout,N); or both
220: $ PetscLayoutSetUp(PetscLayout);
221: $ PetscLayoutGetSize(PetscLayout,PetscInt *);

223:   If range exists, and local size is not set, everything gets computed from the range.

225:   If the local size, global size are already set and range exists then this does nothing.

227: .seealso: PetscLayoutSetLocalSize(), PetscLayoutSetSize(), PetscLayoutGetSize(), PetscLayoutGetLocalSize(), PetscLayout, PetscLayoutDestroy(),
228:           PetscLayoutGetRange(), PetscLayoutGetRanges(), PetscLayoutSetBlockSize(), PetscLayoutGetBlockSize(), PetscLayoutCreate()
229: @*/
230: PetscErrorCode PetscLayoutSetUp(PetscLayout map)
231: {
232:   PetscMPIInt    rank;
233:   PetscInt       p;

236:   if (map->setupcalled) return 0;

238:   if (map->n > 0 && map->bs > 1) {
240:   }
241:   if (map->N > 0 && map->bs > 1) {
243:   }

245:   MPI_Comm_rank(map->comm, &rank);
246:   if (map->n > 0) map->n = map->n/PetscAbs(map->bs);
247:   if (map->N > 0) map->N = map->N/PetscAbs(map->bs);
248:   PetscSplitOwnership(map->comm,&map->n,&map->N);
249:   map->n = map->n*PetscAbs(map->bs);
250:   map->N = map->N*PetscAbs(map->bs);
251:   if (!map->range) {
252:     PetscMalloc1(map->size+1, &map->range);
253:   }
254:   MPI_Allgather(&map->n, 1, MPIU_INT, map->range+1, 1, MPIU_INT, map->comm);

256:   map->range[0] = 0;
257:   for (p = 2; p <= map->size; p++) map->range[p] += map->range[p-1];

259:   map->rstart = map->range[rank];
260:   map->rend   = map->range[rank+1];

262:   /* lock the layout */
263:   map->setupcalled = PETSC_TRUE;
264:   map->oldn = map->n;
265:   map->oldN = map->N;
266:   map->oldbs = map->bs;
267:   return 0;
268: }

270: /*@
271:   PetscLayoutDuplicate - creates a new PetscLayout with the same information as a given one. If the PetscLayout already exists it is destroyed first.

273:   Collective on PetscLayout

275:   Input Parameter:
276: . in - input PetscLayout to be duplicated

278:   Output Parameter:
279: . out - the copy

281:   Level: developer

283:   Notes:
284:     PetscLayoutSetUp() does not need to be called on the resulting PetscLayout

286: .seealso: PetscLayoutCreate(), PetscLayoutDestroy(), PetscLayoutSetUp(), PetscLayoutReference()
287: @*/
288: PetscErrorCode PetscLayoutDuplicate(PetscLayout in,PetscLayout *out)
289: {
290:   MPI_Comm       comm = in->comm;

292:   PetscLayoutDestroy(out);
293:   PetscLayoutCreate(comm,out);
294:   PetscMemcpy(*out,in,sizeof(struct _n_PetscLayout));
295:   if (in->range) {
296:     PetscMalloc1((*out)->size+1,&(*out)->range);
297:     PetscArraycpy((*out)->range,in->range,(*out)->size+1);
298:   }
299:   (*out)->refcnt = 0;
300:   return 0;
301: }

303: /*@
304:   PetscLayoutReference - Causes a PETSc Vec or Mat to share a PetscLayout with one that already exists. Used by Vec/MatDuplicate_XXX()

306:   Collective on PetscLayout

308:   Input Parameter:
309: . in - input PetscLayout to be copied

311:   Output Parameter:
312: . out - the reference location

314:   Level: developer

316:   Notes:
317:     PetscLayoutSetUp() does not need to be called on the resulting PetscLayout

319:   If the out location already contains a PetscLayout it is destroyed

321: .seealso: PetscLayoutCreate(), PetscLayoutDestroy(), PetscLayoutSetUp(), PetscLayoutDuplicate()
322: @*/
323: PetscErrorCode PetscLayoutReference(PetscLayout in,PetscLayout *out)
324: {
325:   in->refcnt++;
326:   PetscLayoutDestroy(out);
327:   *out = in;
328:   return 0;
329: }

331: /*@
332:   PetscLayoutSetISLocalToGlobalMapping - sets a ISLocalGlobalMapping into a PetscLayout

334:   Collective on PetscLayout

336:   Input Parameters:
337: + in - input PetscLayout
338: - ltog - the local to global mapping

340:   Level: developer

342:   Notes:
343:     PetscLayoutSetUp() does not need to be called on the resulting PetscLayout

345:   If the ltog location already contains a PetscLayout it is destroyed

347: .seealso: PetscLayoutCreate(), PetscLayoutDestroy(), PetscLayoutSetUp(), PetscLayoutDuplicate()
348: @*/
349: PetscErrorCode PetscLayoutSetISLocalToGlobalMapping(PetscLayout in,ISLocalToGlobalMapping ltog)
350: {
351:   if (ltog) {
352:     PetscInt bs;

354:     ISLocalToGlobalMappingGetBlockSize(ltog,&bs);
356:     PetscObjectReference((PetscObject)ltog);
357:   }
358:   ISLocalToGlobalMappingDestroy(&in->mapping);
359:   in->mapping = ltog;
360:   return 0;
361: }

363: /*@
364:   PetscLayoutSetLocalSize - Sets the local size for a PetscLayout object.

366:   Collective on PetscLayout

368:   Input Parameters:
369: + map - pointer to the map
370: - n - the local size

372:   Level: developer

374:   Notes:
375:   Call this after the call to PetscLayoutCreate()

377: .seealso: PetscLayoutCreate(), PetscLayoutSetSize(), PetscLayoutGetSize(), PetscLayoutGetLocalSize(), PetscLayoutSetUp()
378:           PetscLayoutGetRange(), PetscLayoutGetRanges(), PetscLayoutSetBlockSize(), PetscLayoutGetBlockSize()
379: @*/
380: PetscErrorCode PetscLayoutSetLocalSize(PetscLayout map,PetscInt n)
381: {
383:   map->n = n;
384:   return 0;
385: }

387: /*@C
388:      PetscLayoutGetLocalSize - Gets the local size for a PetscLayout object.

390:     Not Collective

392:    Input Parameters:
393: .    map - pointer to the map

395:    Output Parameters:
396: .    n - the local size

398:    Level: developer

400:     Notes:
401:        Call this after the call to PetscLayoutSetUp()

403:     Fortran Notes:
404:       Not available from Fortran

406: .seealso: PetscLayoutCreate(), PetscLayoutSetSize(), PetscLayoutGetSize(), PetscLayoutGetLocalSize(), PetscLayoutSetUp()
407:           PetscLayoutGetRange(), PetscLayoutGetRanges(), PetscLayoutSetBlockSize(), PetscLayoutGetBlockSize()

409: @*/
410: PetscErrorCode  PetscLayoutGetLocalSize(PetscLayout map,PetscInt *n)
411: {
412:   *n = map->n;
413:   return 0;
414: }

416: /*@
417:   PetscLayoutSetSize - Sets the global size for a PetscLayout object.

419:   Logically Collective on PetscLayout

421:   Input Parameters:
422: + map - pointer to the map
423: - n - the global size

425:   Level: developer

427:   Notes:
428:   Call this after the call to PetscLayoutCreate()

430: .seealso: PetscLayoutCreate(), PetscLayoutSetLocalSize(), PetscLayoutGetLocalSize(), PetscLayoutGetSize(), PetscLayoutSetUp()
431:           PetscLayoutGetRange(), PetscLayoutGetRanges(), PetscLayoutSetBlockSize(), PetscLayoutGetBlockSize()
432: @*/
433: PetscErrorCode PetscLayoutSetSize(PetscLayout map,PetscInt n)
434: {
435:   map->N = n;
436:   return 0;
437: }

439: /*@
440:   PetscLayoutGetSize - Gets the global size for a PetscLayout object.

442:   Not Collective

444:   Input Parameters:
445: . map - pointer to the map

447:   Output Parameters:
448: . n - the global size

450:   Level: developer

452:   Notes:
453:   Call this after the call to PetscLayoutSetUp()

455: .seealso: PetscLayoutCreate(), PetscLayoutSetLocalSize(), PetscLayoutGetLocalSize(), PetscLayoutSetSize(), PetscLayoutSetUp()
456:           PetscLayoutGetRange(), PetscLayoutGetRanges(), PetscLayoutSetBlockSize(), PetscLayoutGetBlockSize()
457: @*/
458: PetscErrorCode PetscLayoutGetSize(PetscLayout map,PetscInt *n)
459: {
460:   *n = map->N;
461:   return 0;
462: }

464: /*@
465:   PetscLayoutSetBlockSize - Sets the block size for a PetscLayout object.

467:   Logically Collective on PetscLayout

469:   Input Parameters:
470: + map - pointer to the map
471: - bs - the size

473:   Level: developer

475:   Notes:
476:   Call this after the call to PetscLayoutCreate()

478: .seealso: PetscLayoutCreate(), PetscLayoutSetLocalSize(), PetscLayoutGetLocalSize(), PetscLayoutGetBlockSize(),
479:           PetscLayoutGetRange(), PetscLayoutGetRanges(), PetscLayoutSetSize(), PetscLayoutGetSize(), PetscLayoutSetUp()
480: @*/
481: PetscErrorCode PetscLayoutSetBlockSize(PetscLayout map,PetscInt bs)
482: {
483:   if (bs < 0) return 0;
485:   if (map->mapping) {
486:     PetscInt       obs;

488:     ISLocalToGlobalMappingGetBlockSize(map->mapping,&obs);
489:     if (obs > 1) {
490:       ISLocalToGlobalMappingSetBlockSize(map->mapping,bs);
491:     }
492:   }
493:   map->bs = bs;
494:   return 0;
495: }

497: /*@
498:   PetscLayoutGetBlockSize - Gets the block size for a PetscLayout object.

500:   Not Collective

502:   Input Parameters:
503: . map - pointer to the map

505:   Output Parameters:
506: . bs - the size

508:   Level: developer

510:   Notes:
511:   Call this after the call to PetscLayoutSetUp()

513: .seealso: PetscLayoutCreate(), PetscLayoutSetLocalSize(), PetscLayoutGetLocalSize(), PetscLayoutSetSize(), PetscLayoutSetUp()
514:           PetscLayoutGetRange(), PetscLayoutGetRanges(), PetscLayoutSetBlockSize(), PetscLayoutGetSize()
515: @*/
516: PetscErrorCode PetscLayoutGetBlockSize(PetscLayout map,PetscInt *bs)
517: {
518:   *bs = PetscAbs(map->bs);
519:   return 0;
520: }

522: /*@
523:   PetscLayoutGetRange - gets the range of values owned by this process

525:   Not Collective

527:   Input Parameter:
528: . map - pointer to the map

530:   Output Parameters:
531: + rstart - first index owned by this process
532: - rend   - one more than the last index owned by this process

534:   Level: developer

536:   Notes:
537:   Call this after the call to PetscLayoutSetUp()

539: .seealso: PetscLayoutCreate(), PetscLayoutSetLocalSize(), PetscLayoutGetLocalSize(), PetscLayoutSetSize(),
540:           PetscLayoutGetSize(), PetscLayoutGetRanges(), PetscLayoutSetBlockSize(), PetscLayoutGetSize(), PetscLayoutSetUp()
541: @*/
542: PetscErrorCode PetscLayoutGetRange(PetscLayout map,PetscInt *rstart,PetscInt *rend)
543: {
544:   if (rstart) *rstart = map->rstart;
545:   if (rend)   *rend   = map->rend;
546:   return 0;
547: }

549: /*@C
550:      PetscLayoutGetRanges - gets the range of values owned by all processes

552:     Not Collective

554:    Input Parameters:
555: .    map - pointer to the map

557:    Output Parameters:
558: .    range - start of each processors range of indices (the final entry is one more than the
559:              last index on the last process)

561:    Level: developer

563:     Notes:
564:        Call this after the call to PetscLayoutSetUp()

566:     Fortran Notes:
567:       Not available from Fortran

569: .seealso: PetscLayoutCreate(), PetscLayoutSetLocalSize(), PetscLayoutGetLocalSize(), PetscLayoutSetSize(),
570:           PetscLayoutGetSize(), PetscLayoutGetRange(), PetscLayoutSetBlockSize(), PetscLayoutGetSize(), PetscLayoutSetUp()

572: @*/
573: PetscErrorCode  PetscLayoutGetRanges(PetscLayout map,const PetscInt *range[])
574: {
575:   *range = map->range;
576:   return 0;
577: }

579: /*@
580:   PetscLayoutCompare - Compares two layouts

582:   Not Collective

584:   Input Parameters:
585: + mapa - pointer to the first map
586: - mapb - pointer to the second map

588:   Output Parameters:
589: . congruent - PETSC_TRUE if the two layouts are congruent, PETSC_FALSE otherwise

591:   Level: beginner

593:   Notes:

595: .seealso: PetscLayoutCreate(), PetscLayoutSetLocalSize(), PetscLayoutGetLocalSize(), PetscLayoutGetBlockSize(),
596:           PetscLayoutGetRange(), PetscLayoutGetRanges(), PetscLayoutSetSize(), PetscLayoutGetSize(), PetscLayoutSetUp()
597: @*/
598: PetscErrorCode PetscLayoutCompare(PetscLayout mapa,PetscLayout mapb,PetscBool *congruent)
599: {
600:   *congruent = PETSC_FALSE;
601:   if (mapa->N == mapb->N && mapa->range && mapb->range && mapa->size == mapb->size) {
602:     PetscArraycmp(mapa->range,mapb->range,mapa->size+1,congruent);
603:   }
604:   return 0;
605: }