Actual source code: pmap.c

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

 6:  #include private/vecimpl.h
  7: /*@C
  8:      PetscMapInitialize - given a map where you have set either the global or local
  9:            size sets up the map so that it may be used.

 11:     Collective on MPI_Comm

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

 17:    Level: intermediate

 19:     Notes:
 20:        You must call PetscMapSetBlockSize() and either PetscMapSetSize() or PetscMapSetLocalSize()
 21:        before calling this routine.

 23:        Unlike regular PETSc objects you work with a pointer to the object instead of 
 24:      the object directly.

 26:     Fortran Notes: 
 27:       Not available from Fortran

 29: .seealso: PetscMapSetLocalSize(), PetscMapSetSize(), PetscMapGetSize(), PetscMapGetLocalSize(), PetscMap,
 30:           PetscMapGetLocalRange(), PetscMapGetGlobalRange(), PetscMapSetBlockSize(), PetscMapGetBlockSize()

 32: @*/
 35: PetscErrorCode  PetscMapInitialize(MPI_Comm comm,PetscMap *map)
 36: {
 37:   PetscMPIInt    rank,size;
 38:   PetscInt       p;

 42:   MPI_Comm_size(comm, &size);
 43:   MPI_Comm_rank(comm, &rank);
 44:   if (map->bs <=0) {SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"BlockSize not yet set");}
 45:   if (map->n > 0) map->n = map->n/map->bs;
 46:   if (map->N > 0) map->N = map->N/map->bs;
 47:   PetscSplitOwnership(comm,&map->n,&map->N);
 48:   map->n = map->n*map->bs;
 49:   map->N = map->N*map->bs;
 50:   if (!map->range) {
 51:     PetscMalloc((size+1)*sizeof(PetscInt), &map->range);
 52:   }
 53:   MPI_Allgather(&map->n, 1, MPIU_INT, map->range+1, 1, MPIU_INT, comm);

 55:   map->range[0] = 0;
 56:   for(p = 2; p <= size; p++) {
 57:     map->range[p] += map->range[p-1];
 58:   }

 60:   map->rstart = map->range[rank];
 61:   map->rend   = map->range[rank+1];
 62:   return(0);
 63: }

 67: PetscErrorCode  PetscMapCopy(MPI_Comm comm,PetscMap *in,PetscMap *out)
 68: {
 69:   PetscMPIInt    size;
 71:   PetscInt       *range = out->range;

 74:   MPI_Comm_size(comm,&size);
 75:   PetscMemcpy(out,in,sizeof(PetscMap));
 76:   if (!range) {
 77:     PetscMalloc((size+1)*sizeof(PetscInt),&out->range);
 78:   } else {
 79:     out->range = range;
 80:   }
 81:   PetscMemcpy(out->range,in->range,(size+1)*sizeof(PetscInt));
 82:   return(0);
 83: }

 85: /*@C
 86:      PetscMapSetLocalSize - Sets the local size for a PetscMap object.

 88:     Collective on PetscMap

 90:    Input Parameters:
 91: +    map - pointer to the map
 92: -    n - the local size

 94:    Level: intermediate

 96:     Notes:
 97:        Call this after the call to PetscMapInitialize()

 99:        Unlike regular PETSc objects you work with a pointer to the object instead of 
100:      the object directly.

102:     Fortran Notes: 
103:       Not available from Fortran

105: .seealso: PetscMapInitialize(), PetscMapSetSize(), PetscMapGetSize(), PetscMapGetLocalSize(),
106:           PetscMapGetLocalRange(), PetscMapGetGlobalRange(), PetscMapSetBlockSize(), PetscMapGetBlockSize()

108: @*/
111: PetscErrorCode  PetscMapSetLocalSize(PetscMap *map,PetscInt n)
112: {
114:   map->n = n;
115:   return(0);
116: }

118: /*@C
119:      PetscMapGetLocalSize - Gets the local size for a PetscMap object.

121:     Not Collective

123:    Input Parameters:
124: .    map - pointer to the map

126:    Output Parameters:
127: .    n - the local size

129:    Level: intermediate

131:     Notes:
132:        Call this after the call to PetscMapInitialize()

134:        Unlike regular PETSc objects you work with a pointer to the object instead of 
135:      the object directly.

137:     Fortran Notes: 
138:       Not available from Fortran

140: .seealso: PetscMapInitialize(), PetscMapSetSize(), PetscMapGetSize(), PetscMapGetLocalSize(),
141:           PetscMapGetLocalRange(), PetscMapGetGlobalRange(), PetscMapSetBlockSize(), PetscMapGetBlockSize()

143: @*/
146: PetscErrorCode  PetscMapGetLocalSize(PetscMap *map,PetscInt *n)
147: {
149:   *n = map->n;
150:   return(0);
151: }

153: /*@C
154:      PetscMapSetSize - Sets the global size for a PetscMap object.

156:     Collective on PetscMap

158:    Input Parameters:
159: +    map - pointer to the map
160: -    n - the global size

162:    Level: intermediate

164:     Notes:
165:        Call this after the call to PetscMapInitialize()

167:        Unlike regular PETSc objects you work with a pointer to the object instead of 
168:      the object directly.

170:     Fortran Notes: 
171:       Not available from Fortran

173: .seealso: PetscMapInitialize(), PetscMapSetLocalSize(), PetscMapGetLocalSize(), PetscMapGetSize(),
174:           PetscMapGetLocalRange(), PetscMapGetGlobalRange(), PetscMapSetBlockSize(), PetscMapGetBlockSize()

176: @*/
179: PetscErrorCode  PetscMapSetSize(PetscMap *map,PetscInt n)
180: {
182:   map->N = n;
183:   return(0);
184: }

186: /*@C
187:      PetscMapGetSize - Gets the global size for a PetscMap object.

189:     Not Collective

191:    Input Parameters:
192: .    map - pointer to the map

194:    Output Parameters:
195: .    n - the global size

197:    Level: intermediate

199:     Notes:
200:        Call this after the call to PetscMapInitialize()

202:        Unlike regular PETSc objects you work with a pointer to the object instead of 
203:      the object directly.

205:     Fortran Notes: 
206:       Not available from Fortran

208: .seealso: PetscMapInitialize(), PetscMapSetLocalSize(), PetscMapGetLocalSize(), PetscMapSetSize(),
209:           PetscMapGetLocalRange(), PetscMapGetGlobalRange(), PetscMapSetBlockSize(), PetscMapGetBlockSize()

211: @*/
214: PetscErrorCode  PetscMapGetSize(PetscMap *map,PetscInt *n)
215: {
217:   *n = map->n;
218:   return(0);
219: }

221: /*@C
222:      PetscMapSetBlockSize - Sets the block size for a PetscMap object.

224:     Collective on PetscMap

226:    Input Parameters:
227: +    map - pointer to the map
228: -    bs - the size

230:    Level: intermediate

232:     Notes:
233:        Call this after the call to PetscMapInitialize()

235:        Unlike regular PETSc objects you work with a pointer to the object instead of 
236:      the object directly.

238:     Fortran Notes: 
239:       Not available from Fortran

241: .seealso: PetscMapInitialize(), PetscMapSetLocalSize(), PetscMapGetLocalSize(), PetscMapGetBlockSize(),
242:           PetscMapGetLocalRange(), PetscMapGetGlobalRange(), PetscMapSetSize(), PetscMapGetSize()

244: @*/
247: PetscErrorCode  PetscMapSetBlockSize(PetscMap *map,PetscInt bs)
248: {
250:   map->bs = bs;
251:   return(0);
252: }

254: /*@C
255:      PetscMapGetBlockSize - Gets the block size for a PetscMap object.

257:     Not Collective

259:    Input Parameters:
260: .    map - pointer to the map

262:    Output Parameters:
263: .    bs - the size

265:    Level: intermediate

267:     Notes:
268:        Call this after the call to PetscMapInitialize()

270:        Unlike regular PETSc objects you work with a pointer to the object instead of 
271:      the object directly.

273:     Fortran Notes: 
274:       Not available from Fortran

276: .seealso: PetscMapInitialize(), PetscMapSetLocalSize(), PetscMapGetLocalSize(), PetscMapSetSize(),
277:           PetscMapGetLocalRange(), PetscMapGetGlobalRange(), PetscMapSetBlockSize(), PetscMapGetSize()

279: @*/
282: PetscErrorCode  PetscMapGetBlockSize(PetscMap *map,PetscInt *bs)
283: {
285:   *bs = map->bs;
286:   return(0);
287: }


290: /*@C
291:      PetscMapGetLocalRange - gets the range of values owned by this process

293:     Not Collective

295:    Input Parameters:
296: .    map - pointer to the map

298:    Output Parameters:
299: +    rstart - first index owned by this process
300: -    rend - one more than the last index owned by this process

302:    Level: intermediate

304:     Notes:
305:        Call this after the call to PetscMapInitialize()

307:        Unlike regular PETSc objects you work with a pointer to the object instead of 
308:      the object directly.

310:     Fortran Notes: 
311:       Not available from Fortran

313: .seealso: PetscMapInitialize(), PetscMapSetLocalSize(), PetscMapGetLocalSize(), PetscMapSetSize(),
314:           PetscMapGetSize(), PetscMapGetGlobalRange(), PetscMapSetBlockSize(), PetscMapGetSize()

316: @*/
319: PetscErrorCode  PetscMapGetLocalRange(PetscMap *map,PetscInt *rstart,PetscInt *rend)
320: {
322:   if (rstart) *rstart = map->rstart;
323:   if (rend)   *rend   = map->rend;
324:   return(0);
325: }

327: /*@C
328:      PetscMapGetGlobalRange - gets the range of values owned by all processes

330:     Not Collective

332:    Input Parameters:
333: .    map - pointer to the map

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

339:    Level: intermediate

341:     Notes:
342:        Call this after the call to PetscMapInitialize()

344:        Unlike regular PETSc objects you work with a pointer to the object instead of 
345:      the object directly.

347:     Fortran Notes: 
348:       Not available from Fortran

350: .seealso: PetscMapInitialize(), PetscMapSetLocalSize(), PetscMapGetLocalSize(), PetscMapSetSize(),
351:           PetscMapGetSize(), PetscMapGetLocalRange(), PetscMapSetBlockSize(), PetscMapGetSize()

353: @*/
356: PetscErrorCode  PetscMapGetGlobalRange(PetscMap *map,const PetscInt *range[])
357: {
359:   *range = map->range;
360:   return(0);
361: }