Actual source code: block.c

  1: #define PETSCVEC_DLL
  2: /*
  3:      Provides the functions for index sets (IS) defined by a list of integers.
  4:    These are for blocks of data, each block is indicated with a single integer.
  5: */
 6:  #include private/isimpl.h
 7:  #include petscvec.h

  9: typedef struct {
 10:   PetscInt        N,n;            /* number of blocks */
 11:   PetscTruth      sorted;       /* are the blocks sorted? */
 12:   PetscInt        *idx;
 13:   PetscInt        bs;           /* blocksize */
 14: } IS_Block;

 18: PetscErrorCode ISDestroy_Block(IS is)
 19: {
 20:   IS_Block       *is_block = (IS_Block*)is->data;

 24:   PetscFree(is_block->idx);
 25:   PetscFree(is_block);
 26:   return(0);
 27: }

 31: PetscErrorCode ISGetIndices_Block(IS in,const PetscInt *idx[])
 32: {
 33:   IS_Block       *sub = (IS_Block*)in->data;
 35:   PetscInt       i,j,k,bs = sub->bs,n = sub->n,*ii,*jj;

 38:   if (sub->bs == 1) {
 39:     *idx = sub->idx;
 40:   } else {
 41:     PetscMalloc(sub->bs*sub->n*sizeof(PetscInt),&jj);
 42:     *idx = jj;
 43:     k    = 0;
 44:     ii   = sub->idx;
 45:     for (i=0; i<n; i++) {
 46:       for (j=0; j<bs; j++) {
 47:         jj[k++] = ii[i] + j;
 48:       }
 49:     }
 50:   }
 51:   return(0);
 52: }

 56: PetscErrorCode ISRestoreIndices_Block(IS in,const PetscInt *idx[])
 57: {
 58:   IS_Block       *sub = (IS_Block*)in->data;

 62:   if (sub->bs != 1) {
 63:     PetscFree(*(void**)idx);
 64:   } else {
 65:     if (*idx !=  sub->idx) {
 66:       SETERRQ(PETSC_ERR_ARG_WRONG,"Must restore with value from ISGetIndices()");
 67:     }
 68:   }
 69:   return(0);
 70: }

 74: PetscErrorCode ISGetSize_Block(IS is,PetscInt *size)
 75: {
 76:   IS_Block *sub = (IS_Block *)is->data;

 79:   *size = sub->bs*sub->N;
 80:   return(0);
 81: }

 85: PetscErrorCode ISGetLocalSize_Block(IS is,PetscInt *size)
 86: {
 87:   IS_Block *sub = (IS_Block *)is->data;

 90:   *size = sub->bs*sub->n;
 91:   return(0);
 92: }

 96: PetscErrorCode ISInvertPermutation_Block(IS is,PetscInt nlocal,IS *isout)
 97: {
 98:   IS_Block       *sub = (IS_Block *)is->data;
 99:   PetscInt       i,*ii,n = sub->n,*idx = sub->idx;
100:   PetscMPIInt    size;

104:   MPI_Comm_size(((PetscObject)is)->comm,&size);
105:   if (size == 1) {
106:     PetscMalloc(n*sizeof(PetscInt),&ii);
107:     for (i=0; i<n; i++) {
108:       ii[idx[i]] = i;
109:     }
110:     ISCreateBlock(PETSC_COMM_SELF,sub->bs,n,ii,isout);
111:     ISSetPermutation(*isout);
112:     PetscFree(ii);
113:   } else {
114:     SETERRQ(PETSC_ERR_SUP,"No inversion written yet for block IS");
115:   }
116:   return(0);
117: }

121: PetscErrorCode ISView_Block(IS is, PetscViewer viewer)
122: {
123:   IS_Block       *sub = (IS_Block *)is->data;
125:   PetscInt       i,n = sub->n,*idx = sub->idx;
126:   PetscTruth     iascii;

129:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
130:   if (iascii) {
131:     if (is->isperm) {
132:       PetscViewerASCIISynchronizedPrintf(viewer,"Block Index set is permutation\n");
133:     }
134:     PetscViewerASCIISynchronizedPrintf(viewer,"Block size %D\n",sub->bs);
135:     PetscViewerASCIISynchronizedPrintf(viewer,"Number of block indices in set %D\n",n);
136:     PetscViewerASCIISynchronizedPrintf(viewer,"The first indices of each block are\n");
137:     for (i=0; i<n; i++) {
138:       PetscViewerASCIISynchronizedPrintf(viewer,"Block %D Index %D\n",i,idx[i]);
139:     }
140:     PetscViewerFlush(viewer);
141:   } else {
142:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for this object",((PetscObject)viewer)->type_name);
143:   }
144:   return(0);
145: }

149: PetscErrorCode ISSort_Block(IS is)
150: {
151:   IS_Block       *sub = (IS_Block *)is->data;

155:   if (sub->sorted) return(0);
156:   PetscSortInt(sub->n,sub->idx);
157:   sub->sorted = PETSC_TRUE;
158:   return(0);
159: }

163: PetscErrorCode ISSorted_Block(IS is,PetscTruth *flg)
164: {
165:   IS_Block *sub = (IS_Block *)is->data;

168:   *flg = sub->sorted;
169:   return(0);
170: }

174: PetscErrorCode ISDuplicate_Block(IS is,IS *newIS)
175: {
177:   IS_Block       *sub = (IS_Block *)is->data;

180:   ISCreateBlock(((PetscObject)is)->comm,sub->bs,sub->n,sub->idx,newIS);
181:   return(0);
182: }

186: PetscErrorCode ISIdentity_Block(IS is,PetscTruth *ident)
187: {
188:   IS_Block *is_block = (IS_Block*)is->data;
189:   PetscInt i,n = is_block->n,*idx = is_block->idx,bs = is_block->bs;

192:   is->isidentity = PETSC_TRUE;
193:   *ident         = PETSC_TRUE;
194:   for (i=0; i<n; i++) {
195:     if (idx[i] != bs*i) {
196:       is->isidentity = PETSC_FALSE;
197:       *ident         = PETSC_FALSE;
198:       return(0);
199:     }
200:   }
201:   return(0);
202: }

206: static PetscErrorCode ISCopy_Block(IS is,IS isy)
207: {
208:   IS_Block *is_block = (IS_Block*)is->data,*isy_block = (IS_Block*)isy->data;

212:   if (is_block->n != isy_block->n || is_block->N != isy_block->N || is_block->bs != isy_block->bs) SETERRQ(PETSC_ERR_ARG_INCOMP,"Index sets incompatible");
213:   isy_block->sorted = is_block->sorted;
214:   PetscMemcpy(isy_block->idx,is_block->idx,is_block->n*sizeof(PetscInt));
215:   return(0);
216: }

218: static struct _ISOps myops = { ISGetSize_Block,
219:                                ISGetLocalSize_Block,
220:                                ISGetIndices_Block,
221:                                ISRestoreIndices_Block,
222:                                ISInvertPermutation_Block,
223:                                ISSort_Block,
224:                                ISSorted_Block,
225:                                ISDuplicate_Block,
226:                                ISDestroy_Block,
227:                                ISView_Block,
228:                                ISIdentity_Block,
229:                                ISCopy_Block };
232: /*@
233:    ISCreateBlock - Creates a data structure for an index set containing
234:    a list of integers. The indices are relative to entries, not blocks. 

236:    Collective on MPI_Comm

238:    Input Parameters:
239: +  n - the length of the index set (the number of blocks)
240: .  bs - number of elements in each block
241: .  idx - the list of integers
242: -  comm - the MPI communicator

244:    Output Parameter:
245: .  is - the new index set

247:    Notes:
248:    When the communicator is not MPI_COMM_SELF, the operations on the 
249:    index sets, IS, are NOT conceptually the same as MPI_Group operations. 
250:    The index sets are then distributed sets of indices and thus certain operations
251:    on them are collective. 

253:    Example:
254:    If you wish to index the values {0,1,4,5}, then use
255:    a block size of 2 and idx of {0,4}.

257:    Level: beginner

259:   Concepts: IS^block
260:   Concepts: index sets^block
261:   Concepts: block^index set

263: .seealso: ISCreateStride(), ISCreateGeneral(), ISAllGather()
264: @*/
265: PetscErrorCode  ISCreateBlock(MPI_Comm comm,PetscInt bs,PetscInt n,const PetscInt idx[],IS *is)
266: {
268:   PetscInt       i,min,max;
269:   IS             Nindex;
270:   IS_Block       *sub;
271:   PetscTruth     sorted = PETSC_TRUE;

275:   if (n < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"length < 0");
277:   *is = PETSC_NULL;
278: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
279:   ISInitializePackage(PETSC_NULL);
280: #endif

282:   PetscHeaderCreate(Nindex,_p_IS,struct _ISOps,IS_COOKIE,IS_BLOCK,"IS",comm,ISDestroy,ISView);
283:   PetscNewLog(Nindex,IS_Block,&sub);
284:   PetscMalloc(n*sizeof(PetscInt),&sub->idx);
285:   PetscLogObjectMemory(Nindex,n*sizeof(PetscInt));
286:   sub->n = n;
287:   MPI_Allreduce(&n,&sub->N,1,MPIU_INT,MPI_SUM,comm);
288:   for (i=1; i<n; i++) {
289:     if (idx[i] < idx[i-1]) {sorted = PETSC_FALSE; break;}
290:   }
291:   if (n) {min = max = idx[0];} else {min = max = 0;}
292:   for (i=1; i<n; i++) {
293:     if (idx[i] < min) min = idx[i];
294:     if (idx[i] > max) max = idx[i];
295:   }
296:   PetscMemcpy(sub->idx,idx,n*sizeof(PetscInt));
297:   sub->sorted     = sorted;
298:   sub->bs         = bs;
299:   Nindex->min     = min;
300:   Nindex->max     = max;
301:   Nindex->data    = (void*)sub;
302:   PetscMemcpy(Nindex->ops,&myops,sizeof(myops));
303:   Nindex->isperm  = PETSC_FALSE;
304:   *is = Nindex; return(0);
305: }


310: /*@C
311:    ISBlockGetIndices - Gets the indices associated with each block.

313:    Not Collective

315:    Input Parameter:
316: .  is - the index set

318:    Output Parameter:
319: .  idx - the integer indices

321:    Level: intermediate

323:    Concepts: IS^block
324:    Concepts: index sets^getting indices
325:    Concepts: index sets^block

327: .seealso: ISGetIndices(), ISBlockRestoreIndices()
328: @*/
329: PetscErrorCode  ISBlockGetIndices(IS in,const PetscInt *idx[])
330: {
331:   IS_Block *sub;

336:   if (((PetscObject)in)->type != IS_BLOCK) SETERRQ(PETSC_ERR_ARG_WRONG,"Not a block index set");

338:   sub = (IS_Block*)in->data;
339:   *idx = sub->idx;
340:   return(0);
341: }

345: /*@C
346:    ISBlockRestoreIndices - Restores the indices associated with each block.

348:    Not Collective

350:    Input Parameter:
351: .  is - the index set

353:    Output Parameter:
354: .  idx - the integer indices

356:    Level: intermediate

358:    Concepts: IS^block
359:    Concepts: index sets^getting indices
360:    Concepts: index sets^block

362: .seealso: ISRestoreIndices(), ISBlockGetIndices()
363: @*/
364: PetscErrorCode  ISBlockRestoreIndices(IS is,const PetscInt *idx[])
365: {
369:   if (((PetscObject)is)->type != IS_BLOCK) SETERRQ(PETSC_ERR_ARG_WRONG,"Not a block index set");
370:   return(0);
371: }

375: /*@
376:    ISBlockGetBlockSize - Returns the number of elements in a block.

378:    Not Collective

380:    Input Parameter:
381: .  is - the index set

383:    Output Parameter:
384: .  size - the number of elements in a block

386:    Level: intermediate

388:    Concepts: IS^block size
389:    Concepts: index sets^block size

391: .seealso: ISBlockGetSize(), ISGetSize(), ISBlock(), ISCreateBlock()
392: @*/
393: PetscErrorCode  ISBlockGetBlockSize(IS is,PetscInt *size)
394: {
395:   IS_Block *sub;

400:   if (((PetscObject)is)->type != IS_BLOCK) SETERRQ(PETSC_ERR_ARG_WRONG,"Not a block index set");

402:   sub = (IS_Block *)is->data;
403:   *size = sub->bs;
404:   return(0);
405: }

409: /*@
410:    ISBlock - Checks whether an index set is blocked.

412:    Not Collective

414:    Input Parameter:
415: .  is - the index set

417:    Output Parameter:
418: .  flag - PETSC_TRUE if a block index set, else PETSC_FALSE

420:    Level: intermediate

422:    Concepts: IS^block
423:    Concepts: index sets^block

425: .seealso: ISBlockGetSize(), ISGetSize(), ISBlockGetBlockSize(), ISCreateBlock()
426: @*/
427: PetscErrorCode  ISBlock(IS is,PetscTruth *flag)
428: {
432:   if (((PetscObject)is)->type != IS_BLOCK) *flag = PETSC_FALSE;
433:   else                          *flag = PETSC_TRUE;
434:   return(0);
435: }

439: /*@
440:    ISBlockGetLocalSize - Returns the local number of blocks in the index set.

442:    Not Collective

444:    Input Parameter:
445: .  is - the index set

447:    Output Parameter:
448: .  size - the local number of blocks

450:    Level: intermediate

452:    Concepts: IS^block sizes
453:    Concepts: index sets^block sizes

455: .seealso: ISBlockGetBlockSize(), ISBlockGetSize(), ISGetSize(), ISBlock(), ISCreateBlock()
456: @*/
457: PetscErrorCode  ISBlockGetLocalSize(IS is,PetscInt *size)
458: {
459:   IS_Block *sub;

464:   if (((PetscObject)is)->type != IS_BLOCK) SETERRQ(PETSC_ERR_ARG_WRONG,"Not a block index set");

466:   sub = (IS_Block *)is->data;
467:   *size = sub->n;
468:   return(0);
469: }

473: /*@
474:    ISBlockGetSize - Returns the global number of blocks in the index set.

476:    Not Collective

478:    Input Parameter:
479: .  is - the index set

481:    Output Parameter:
482: .  size - the global number of blocks

484:    Level: intermediate

486:    Concepts: IS^block sizes
487:    Concepts: index sets^block sizes

489: .seealso: ISBlockGetBlockSize(), ISBlockGetLocalSize(), ISGetSize(), ISBlock(), ISCreateBlock()
490: @*/
491: PetscErrorCode  ISBlockGetSize(IS is,PetscInt *size)
492: {
493:   IS_Block *sub;

498:   if (((PetscObject)is)->type != IS_BLOCK) SETERRQ(PETSC_ERR_ARG_WRONG,"Not a block index set");

500:   sub = (IS_Block *)is->data;
501:   *size = sub->N;
502:   return(0);
503: }