Actual source code: pmap.c

petsc-3.8.4 2018-03-24
Report Typos and Errors

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

  6:  #include <petscis.h>
  7:  #include <petscsf.h>

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

 12:   Collective on MPI_Comm

 14:   Input Parameters:
 15: + comm - the MPI communicator
 16: - map - pointer to the map

 18:   Level: advanced

 20:   Notes:
 21:   Typical calling sequence
 22: .vb
 23:        PetscLayoutCreate(MPI_Comm,PetscLayout *);
 24:        PetscLayoutSetBlockSize(PetscLayout,1);
 25:        PetscLayoutSetSize(PetscLayout,N) // or PetscLayoutSetLocalSize(PetscLayout,n);
 26:        PetscLayoutSetUp(PetscLayout);
 27: .ve
 28:   Optionally use any of the following:

 30: + PetscLayoutGetSize(PetscLayout,PetscInt *);
 31: . PetscLayoutGetLocalSize(PetscLayout,PetscInt *);
 32: . PetscLayoutGetRange(PetscLayout,PetscInt *rstart,PetscInt *rend);
 33: . PetscLayoutGetRanges(PetscLayout,const PetscInt *range[]);
 34: - PetscLayoutDestroy(PetscLayout*);

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

 39: .seealso: PetscLayoutSetLocalSize(), PetscLayoutSetSize(), PetscLayoutGetSize(), PetscLayoutGetLocalSize(), PetscLayout, PetscLayoutDestroy(),
 40:           PetscLayoutGetRange(), PetscLayoutGetRanges(), PetscLayoutSetBlockSize(), PetscLayoutGetBlockSize(), PetscLayoutSetUp()

 42: @*/
 43: PetscErrorCode PetscLayoutCreate(MPI_Comm comm,PetscLayout *map)
 44: {

 48:   PetscNew(map);

 50:   (*map)->comm   = comm;
 51:   (*map)->bs     = -1;
 52:   (*map)->n      = -1;
 53:   (*map)->N      = -1;
 54:   (*map)->range  = NULL;
 55:   (*map)->rstart = 0;
 56:   (*map)->rend   = 0;
 57:   return(0);
 58: }

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

 63:   Collective on MPI_Comm

 65:   Input Parameters:
 66: . map - the PetscLayout

 68:   Level: developer

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

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

 77: @*/
 78: PetscErrorCode PetscLayoutDestroy(PetscLayout *map)
 79: {

 83:   if (!*map) return(0);
 84:   if (!(*map)->refcnt--) {
 85:     PetscFree((*map)->range);
 86:     ISLocalToGlobalMappingDestroy(&(*map)->mapping);
 87:     PetscFree((*map));
 88:   }
 89:   *map = NULL;
 90:   return(0);
 91: }

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

 97:   Collective on MPI_Comm

 99:   Input Parameters:
100: . map - pointer to the map

102:   Level: developer

104:   Notes: Typical calling sequence
105: $ PetscLayoutCreate(MPI_Comm,PetscLayout *);
106: $ PetscLayoutSetBlockSize(PetscLayout,1);
107: $ PetscLayoutSetSize(PetscLayout,n) or PetscLayoutSetLocalSize(PetscLayout,N); or both
108: $ PetscLayoutSetUp(PetscLayout);
109: $ PetscLayoutGetSize(PetscLayout,PetscInt *);


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

114: .seealso: PetscLayoutSetLocalSize(), PetscLayoutSetSize(), PetscLayoutGetSize(), PetscLayoutGetLocalSize(), PetscLayout, PetscLayoutDestroy(),
115:           PetscLayoutGetRange(), PetscLayoutGetRanges(), PetscLayoutSetBlockSize(), PetscLayoutGetBlockSize(), PetscLayoutCreate()
116: @*/
117: PetscErrorCode PetscLayoutSetUp(PetscLayout map)
118: {
119:   PetscMPIInt    rank,size;
120:   PetscInt       p;

124:   if ((map->n >= 0) && (map->N >= 0) && (map->range)) return(0);

126:   if (map->n > 0 && map->bs > 1) {
127:     if (map->n % map->bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Local matrix size %D must be divisible by blocksize %D",map->n,map->bs);
128:   }
129:   if (map->N > 0 && map->bs > 1) {
130:     if (map->N % map->bs) SETERRQ2(map->comm,PETSC_ERR_PLIB,"Global matrix size %D must be divisible by blocksize %D",map->N,map->bs);
131:   }

133:   MPI_Comm_size(map->comm, &size);
134:   MPI_Comm_rank(map->comm, &rank);
135:   if (map->n > 0) map->n = map->n/PetscAbs(map->bs);
136:   if (map->N > 0) map->N = map->N/PetscAbs(map->bs);
137:   PetscSplitOwnership(map->comm,&map->n,&map->N);
138:   map->n = map->n*PetscAbs(map->bs);
139:   map->N = map->N*PetscAbs(map->bs);
140:   if (!map->range) {
141:     PetscMalloc1(size+1, &map->range);
142:   }
143:   MPI_Allgather(&map->n, 1, MPIU_INT, map->range+1, 1, MPIU_INT, map->comm);

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

148:   map->rstart = map->range[rank];
149:   map->rend   = map->range[rank+1];
150:   return(0);
151: }

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

156:   Collective on PetscLayout

158:   Input Parameter:
159: . in - input PetscLayout to be duplicated

161:   Output Parameter:
162: . out - the copy

164:   Level: developer

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

168: .seealso: PetscLayoutCreate(), PetscLayoutDestroy(), PetscLayoutSetUp(), PetscLayoutReference()
169: @*/
170: PetscErrorCode PetscLayoutDuplicate(PetscLayout in,PetscLayout *out)
171: {
172:   PetscMPIInt    size;
174:   MPI_Comm       comm = in->comm;

177:   PetscLayoutDestroy(out);
178:   PetscLayoutCreate(comm,out);
179:   MPI_Comm_size(comm,&size);
180:   PetscMemcpy(*out,in,sizeof(struct _n_PetscLayout));
181:   PetscMalloc1(size+1,&(*out)->range);
182:   PetscMemcpy((*out)->range,in->range,(size+1)*sizeof(PetscInt));

184:   (*out)->refcnt = 0;
185:   return(0);
186: }

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

191:   Collective on PetscLayout

193:   Input Parameter:
194: . in - input PetscLayout to be copied

196:   Output Parameter:
197: . out - the reference location

199:   Level: developer

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

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

205: .seealso: PetscLayoutCreate(), PetscLayoutDestroy(), PetscLayoutSetUp(), PetscLayoutDuplicate()
206: @*/
207: PetscErrorCode PetscLayoutReference(PetscLayout in,PetscLayout *out)
208: {

212:   in->refcnt++;
213:   PetscLayoutDestroy(out);
214:   *out = in;
215:   return(0);
216: }

218: /*@
219:   PetscLayoutSetISLocalToGlobalMapping - sets a ISLocalGlobalMapping into a PetscLayout

221:   Collective on PetscLayout

223:   Input Parameter:
224: + in - input PetscLayout
225: - ltog - the local to global mapping


228:   Level: developer

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

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

234: .seealso: PetscLayoutCreate(), PetscLayoutDestroy(), PetscLayoutSetUp(), PetscLayoutDuplicate()
235: @*/
236: PetscErrorCode PetscLayoutSetISLocalToGlobalMapping(PetscLayout in,ISLocalToGlobalMapping ltog)
237: {
239:   PetscInt       bs;

242:   ISLocalToGlobalMappingGetBlockSize(ltog,&bs);
243:   if (in->bs > 0 && in->bs != bs) SETERRQ2(in->comm,PETSC_ERR_PLIB,"Blocksize of layout %D must match that of mapping %D",in->bs,bs);
244:   PetscObjectReference((PetscObject)ltog);
245:   ISLocalToGlobalMappingDestroy(&in->mapping);
246:   in->mapping = ltog;
247:   return(0);
248: }

250: /*@
251:   PetscLayoutSetLocalSize - Sets the local size for a PetscLayout object.

253:   Collective on PetscLayout

255:   Input Parameters:
256: + map - pointer to the map
257: - n - the local size

259:   Level: developer

261:   Notes:
262:   Call this after the call to PetscLayoutCreate()

264: .seealso: PetscLayoutCreate(), PetscLayoutSetSize(), PetscLayoutGetSize(), PetscLayoutGetLocalSize(), PetscLayoutSetUp()
265:           PetscLayoutGetRange(), PetscLayoutGetRanges(), PetscLayoutSetBlockSize(), PetscLayoutGetBlockSize()
266: @*/
267: PetscErrorCode PetscLayoutSetLocalSize(PetscLayout map,PetscInt n)
268: {
270:   if (map->bs > 1 && n % map->bs) SETERRQ2(map->comm,PETSC_ERR_ARG_INCOMP,"Local size %D not compatible with block size %D",n,map->bs);
271:   map->n = n;
272:   return(0);
273: }

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

278:     Not Collective

280:    Input Parameters:
281: .    map - pointer to the map

283:    Output Parameters:
284: .    n - the local size

286:    Level: developer

288:     Notes:
289:        Call this after the call to PetscLayoutSetUp()

291:     Fortran Notes:
292:       Not available from Fortran

294: .seealso: PetscLayoutCreate(), PetscLayoutSetSize(), PetscLayoutGetSize(), PetscLayoutGetLocalSize(), PetscLayoutSetUp()
295:           PetscLayoutGetRange(), PetscLayoutGetRanges(), PetscLayoutSetBlockSize(), PetscLayoutGetBlockSize()

297: @*/
298: PetscErrorCode  PetscLayoutGetLocalSize(PetscLayout map,PetscInt *n)
299: {
301:   *n = map->n;
302:   return(0);
303: }

305: /*@
306:   PetscLayoutSetSize - Sets the global size for a PetscLayout object.

308:   Logically Collective on PetscLayout

310:   Input Parameters:
311: + map - pointer to the map
312: - n - the global size

314:   Level: developer

316:   Notes:
317:   Call this after the call to PetscLayoutCreate()

319: .seealso: PetscLayoutCreate(), PetscLayoutSetLocalSize(), PetscLayoutGetLocalSize(), PetscLayoutGetSize(), PetscLayoutSetUp()
320:           PetscLayoutGetRange(), PetscLayoutGetRanges(), PetscLayoutSetBlockSize(), PetscLayoutGetBlockSize()
321: @*/
322: PetscErrorCode PetscLayoutSetSize(PetscLayout map,PetscInt n)
323: {
325:   map->N = n;
326:   return(0);
327: }

329: /*@
330:   PetscLayoutGetSize - Gets the global size for a PetscLayout object.

332:   Not Collective

334:   Input Parameters:
335: . map - pointer to the map

337:   Output Parameters:
338: . n - the global size

340:   Level: developer

342:   Notes:
343:   Call this after the call to PetscLayoutSetUp()

345: .seealso: PetscLayoutCreate(), PetscLayoutSetLocalSize(), PetscLayoutGetLocalSize(), PetscLayoutSetSize(), PetscLayoutSetUp()
346:           PetscLayoutGetRange(), PetscLayoutGetRanges(), PetscLayoutSetBlockSize(), PetscLayoutGetBlockSize()
347: @*/
348: PetscErrorCode PetscLayoutGetSize(PetscLayout map,PetscInt *n)
349: {
351:   *n = map->N;
352:   return(0);
353: }

355: /*@
356:   PetscLayoutSetBlockSize - Sets the block size for a PetscLayout object.

358:   Logically Collective on PetscLayout

360:   Input Parameters:
361: + map - pointer to the map
362: - bs - the size

364:   Level: developer

366:   Notes:
367:   Call this after the call to PetscLayoutCreate()

369: .seealso: PetscLayoutCreate(), PetscLayoutSetLocalSize(), PetscLayoutGetLocalSize(), PetscLayoutGetBlockSize(),
370:           PetscLayoutGetRange(), PetscLayoutGetRanges(), PetscLayoutSetSize(), PetscLayoutGetSize(), PetscLayoutSetUp()
371: @*/
372: PetscErrorCode PetscLayoutSetBlockSize(PetscLayout map,PetscInt bs)
373: {
375:   if (bs < 0) return(0);
376:   if (map->n > 0 && map->n % bs) SETERRQ2(map->comm,PETSC_ERR_ARG_INCOMP,"Local size %D not compatible with block size %D",map->n,bs);
377:   if (map->mapping) {
378:     PetscInt       obs;

381:     ISLocalToGlobalMappingGetBlockSize(map->mapping,&obs);
382:     if (obs > 1) {
383:       ISLocalToGlobalMappingSetBlockSize(map->mapping,bs);
384:     }
385:   }
386:   map->bs = bs;
387:   return(0);
388: }

390: /*@
391:   PetscLayoutGetBlockSize - Gets the block size for a PetscLayout object.

393:   Not Collective

395:   Input Parameters:
396: . map - pointer to the map

398:   Output Parameters:
399: . bs - the size

401:   Level: developer

403:   Notes:
404:   Call this after the call to PetscLayoutSetUp()

406: .seealso: PetscLayoutCreate(), PetscLayoutSetLocalSize(), PetscLayoutGetLocalSize(), PetscLayoutSetSize(), PetscLayoutSetUp()
407:           PetscLayoutGetRange(), PetscLayoutGetRanges(), PetscLayoutSetBlockSize(), PetscLayoutGetSize()
408: @*/
409: PetscErrorCode PetscLayoutGetBlockSize(PetscLayout map,PetscInt *bs)
410: {
412:   *bs = PetscAbs(map->bs);
413:   return(0);
414: }

416: /*@
417:   PetscLayoutGetRange - gets the range of values owned by this process

419:   Not Collective

421:   Input Parameters:
422: . map - pointer to the map

424:   Output Parameters:
425: + rstart - first index owned by this process
426: - rend   - one more than the last index owned by this process

428:   Level: developer

430:   Notes:
431:   Call this after the call to PetscLayoutSetUp()

433: .seealso: PetscLayoutCreate(), PetscLayoutSetLocalSize(), PetscLayoutGetLocalSize(), PetscLayoutSetSize(),
434:           PetscLayoutGetSize(), PetscLayoutGetRanges(), PetscLayoutSetBlockSize(), PetscLayoutGetSize(), PetscLayoutSetUp()
435: @*/
436: PetscErrorCode PetscLayoutGetRange(PetscLayout map,PetscInt *rstart,PetscInt *rend)
437: {
439:   if (rstart) *rstart = map->rstart;
440:   if (rend)   *rend   = map->rend;
441:   return(0);
442: }

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

447:     Not Collective

449:    Input Parameters:
450: .    map - pointer to the map

452:    Output Parameters:
453: .    range - start of each processors range of indices (the final entry is one more then the
454:              last index on the last process)

456:    Level: developer

458:     Notes:
459:        Call this after the call to PetscLayoutSetUp()

461:     Fortran Notes:
462:       Not available from Fortran

464: .seealso: PetscLayoutCreate(), PetscLayoutSetLocalSize(), PetscLayoutGetLocalSize(), PetscLayoutSetSize(),
465:           PetscLayoutGetSize(), PetscLayoutGetRange(), PetscLayoutSetBlockSize(), PetscLayoutGetSize(), PetscLayoutSetUp()

467: @*/
468: PetscErrorCode  PetscLayoutGetRanges(PetscLayout map,const PetscInt *range[])
469: {
471:   *range = map->range;
472:   return(0);
473: }

475: /*@C
476:    PetscSFSetGraphLayout - Set a parallel star forest via global indices and a PetscLayout

478:    Collective

480:    Input Arguments:
481: +  sf - star forest
482: .  layout - PetscLayout defining the global space
483: .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
484: .  ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage
485: -  iremote - remote locations of root vertices for each leaf on the current process

487:    Level: intermediate

489: .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph()
490: @*/
491: PetscErrorCode PetscSFSetGraphLayout(PetscSF sf,PetscLayout layout,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscInt *iremote)
492: {
494:   PetscInt       i,nroots;
495:   PetscSFNode    *remote;

498:   PetscLayoutGetLocalSize(layout,&nroots);
499:   PetscMalloc1(nleaves,&remote);
500:   for (i=0; i<nleaves; i++) {
501:     PetscInt owner = -1;
502:     PetscLayoutFindOwner(layout,iremote[i],&owner);
503:     remote[i].rank  = owner;
504:     remote[i].index = iremote[i] - layout->range[owner];
505:   }
506:   PetscSFSetGraph(sf,nroots,nleaves,ilocal,localmode,remote,PETSC_OWN_POINTER);
507:   return(0);
508: }

510: /*@
511:   PetscLayoutCompare - Compares two layouts

513:   Not Collective

515:   Input Parameters:
516: + mapa - pointer to the first map
517: - mapb - pointer to the second map

519:   Output Parameters:
520: . congruent - PETSC_TRUE if the two layouts are congruent, PETSC_FALSE otherwise

522:   Level: beginner

524:   Notes:

526: .seealso: PetscLayoutCreate(), PetscLayoutSetLocalSize(), PetscLayoutGetLocalSize(), PetscLayoutGetBlockSize(),
527:           PetscLayoutGetRange(), PetscLayoutGetRanges(), PetscLayoutSetSize(), PetscLayoutGetSize(), PetscLayoutSetUp()
528: @*/
529: PetscErrorCode PetscLayoutCompare(PetscLayout mapa,PetscLayout mapb,PetscBool *congruent)
530: {
532:   PetscMPIInt    sizea,sizeb;

535:   *congruent = PETSC_FALSE;
536:   MPI_Comm_size(mapa->comm,&sizea);
537:   MPI_Comm_size(mapb->comm,&sizeb);
538:   if (mapa->N == mapb->N && mapa->range && mapb->range && sizea == sizeb) {
539:     PetscMemcmp(mapa->range,mapb->range,(sizea+1)*sizeof(PetscInt),congruent);
540:   }
541:   return(0);
542: }