Actual source code: block.c

  1: /*$Id: block.c,v 1.54 2001/03/23 23:21:10 balay Exp $*/
  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 src/vec/is/isimpl.h
 7:  #include petscsys.h

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

 16: int ISDestroy_Block(IS is)
 17: {
 18:   IS_Block *is_block = (IS_Block*)is->data;

 22:   PetscFree(is_block->idx);
 23:   PetscFree(is_block);
 24:   PetscLogObjectDestroy(is);
 25:   PetscHeaderDestroy(is); return(0);
 26: }

 28: int ISGetIndices_Block(IS in,int **idx)
 29: {
 30:   IS_Block *sub = (IS_Block*)in->data;
 31:   int      i,j,k,bs = sub->bs,n = sub->n,*ii,*jj,ierr;

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

 50: int ISRestoreIndices_Block(IS in,int **idx)
 51: {
 52:   IS_Block *sub = (IS_Block*)in->data;
 53:   int      ierr;

 56:   if (sub->bs != 1) {
 57:     PetscFree(*idx);
 58:   } else {
 59:     if (*idx !=  sub->idx) {
 60:       SETERRQ(PETSC_ERR_ARG_WRONG,"Must restore with value from ISGetIndices()");
 61:     }
 62:   }
 63:   return(0);
 64: }

 66: int ISGetSize_Block(IS is,int *size)
 67: {
 68:   IS_Block *sub = (IS_Block *)is->data;

 71:   *size = sub->bs*sub->N;
 72:   return(0);
 73: }

 75: int ISGetLocalSize_Block(IS is,int *size)
 76: {
 77:   IS_Block *sub = (IS_Block *)is->data;

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

 84: int ISInvertPermutation_Block(IS is,int nlocal,IS *isout)
 85: {
 86:   IS_Block *sub = (IS_Block *)is->data;
 87:   int      i,ierr,*ii,n = sub->n,*idx = sub->idx,size;

 90:   MPI_Comm_size(is->comm,&size);
 91:   if (size == 1) {
 92:     PetscMalloc((n+1)*sizeof(int),&ii);
 93:     for (i=0; i<n; i++) {
 94:       ii[idx[i]] = i;
 95:     }
 96:     ISCreateBlock(PETSC_COMM_SELF,sub->bs,n,ii,isout);
 97:     ISSetPermutation(*isout);
 98:     PetscFree(ii);
 99:   } else {
100:     SETERRQ(1,"No inversion written yet for block IS");
101:   }
102:   return(0);
103: }

105: int ISView_Block(IS is, PetscViewer viewer)
106: {
107:   IS_Block    *sub = (IS_Block *)is->data;
108:   int         i,n = sub->n,*idx = sub->idx,ierr;
109:   PetscTruth  isascii;

112:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
113:   if (isascii) {
114:     if (is->isperm) {
115:       PetscViewerASCIISynchronizedPrintf(viewer,"Block Index set is permutationn");
116:     }
117:     PetscViewerASCIISynchronizedPrintf(viewer,"Block size %dn",sub->bs);
118:     PetscViewerASCIISynchronizedPrintf(viewer,"Number of block indices in set %dn",n);
119:     PetscViewerASCIISynchronizedPrintf(viewer,"The first indices of each block aren");
120:     for (i=0; i<n; i++) {
121:       PetscViewerASCIISynchronizedPrintf(viewer,"%d %dn",i,idx[i]);
122:     }
123:     PetscViewerFlush(viewer);
124:   } else {
125:     SETERRQ1(1,"Viewer type %s not supported for this object",((PetscObject)viewer)->type_name);
126:   }
127:   return(0);
128: }

130: int ISSort_Block(IS is)
131: {
132:   IS_Block *sub = (IS_Block *)is->data;
133:   int      ierr;

136:   if (sub->sorted) return(0);
137:   PetscSortInt(sub->n,sub->idx);
138:   sub->sorted = PETSC_TRUE;
139:   return(0);
140: }

142: int ISSorted_Block(IS is,PetscTruth *flg)
143: {
144:   IS_Block *sub = (IS_Block *)is->data;

147:   *flg = sub->sorted;
148:   return(0);
149: }

151: int ISDuplicate_Block(IS is,IS *newIS)
152: {
153:   int      ierr;
154:   IS_Block *sub = (IS_Block *)is->data;

157:   ISCreateBlock(is->comm,sub->bs,sub->n,sub->idx,newIS);
158:   return(0);
159: }

161: int ISIdentity_Block(IS is,PetscTruth *ident)
162: {
163:   IS_Block *is_block = (IS_Block*)is->data;
164:   int      i,n = is_block->n,*idx = is_block->idx,bs = is_block->bs;

167:   is->isidentity = PETSC_TRUE;
168:   *ident         = PETSC_TRUE;
169:   for (i=0; i<n; i++) {
170:     if (idx[i] != bs*i) {
171:       is->isidentity = PETSC_FALSE;
172:       *ident         = PETSC_FALSE;
173:       return(0);
174:     }
175:   }
176:   return(0);
177: }

179: static struct _ISOps myops = { ISGetSize_Block,
180:                                ISGetLocalSize_Block,
181:                                ISGetIndices_Block,
182:                                ISRestoreIndices_Block,
183:                                ISInvertPermutation_Block,
184:                                ISSort_Block,
185:                                ISSorted_Block,
186:                                ISDuplicate_Block,
187:                                ISDestroy_Block,
188:                                ISView_Block,
189:                                ISIdentity_Block };
190: /*@C
191:    ISCreateBlock - Creates a data structure for an index set containing
192:    a list of integers. The indices are relative to entries, not blocks. 

194:    Collective on MPI_Comm

196:    Input Parameters:
197: +  n - the length of the index set (the number of blocks)
198: .  bs - number of elements in each block
199: .  idx - the list of integers
200: -  comm - the MPI communicator

202:    Output Parameter:
203: .  is - the new index set

205:    Notes:
206:    When the communicator is not MPI_COMM_SELF, the operations on the 
207:    index sets, IS, are NOT conceptually the same as MPI_Group operations. 
208:    The index sets are then distributed sets of indices and thus certain operations
209:    on them are collective. 

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

215:    Level: beginner

217:   Concepts: IS^block
218:   Concepts: index sets^block
219:   Concepts: block^index set

221: .seealso: ISCreateStride(), ISCreateGeneral(), ISAllGather()
222: @*/
223: int ISCreateBlock(MPI_Comm comm,int bs,int n,const int idx[],IS *is)
224: {
225:   int        i,min,max,ierr;
226:   IS         Nindex;
227:   IS_Block   *sub;
228:   PetscTruth sorted = PETSC_TRUE;

231:   *is = 0;
232:   PetscHeaderCreate(Nindex,_p_IS,struct _ISOps,IS_COOKIE,IS_BLOCK,"IS",comm,ISDestroy,ISView);
233:   PetscLogObjectCreate(Nindex);
234:   PetscNew(IS_Block,&sub);
235:   PetscLogObjectMemory(Nindex,sizeof(IS_Block)+n*sizeof(int)+sizeof(struct _p_IS));
236:   ierr   = PetscMalloc((n+1)*sizeof(int),&sub->idx);
237:   sub->n = n;
238:   MPI_Allreduce(&n,&sub->N,1,MPI_INT,MPI_SUM,comm);
239:   for (i=1; i<n; i++) {
240:     if (idx[i] < idx[i-1]) {sorted = PETSC_FALSE; break;}
241:   }
242:   if (n) {min = max = idx[0];} else {min = max = 0;}
243:   for (i=1; i<n; i++) {
244:     if (idx[i] < min) min = idx[i];
245:     if (idx[i] > max) max = idx[i];
246:   }
247:   PetscMemcpy(sub->idx,idx,n*sizeof(int));
248:   sub->sorted     = sorted;
249:   sub->bs         = bs;
250:   Nindex->min     = min;
251:   Nindex->max     = max;
252:   Nindex->data    = (void*)sub;
253:   PetscMemcpy(Nindex->ops,&myops,sizeof(myops));
254:   Nindex->isperm  = PETSC_FALSE;
255:   *is = Nindex; return(0);
256: }


259: /*@C
260:    ISBlockGetIndices - Gets the indices associated with each block.

262:    Not Collective

264:    Input Parameter:
265: .  is - the index set

267:    Output Parameter:
268: .  idx - the integer indices

270:    Level: intermediate

272:    Concepts: IS^block
273:    Concepts: index sets^getting indices
274:    Concepts: index sets^block

276: .seealso: ISGetIndices(), ISBlockRestoreIndices()
277: @*/
278: int ISBlockGetIndices(IS in,int *idx[])
279: {
280:   IS_Block *sub;

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

287:   sub = (IS_Block*)in->data;
288:   *idx = sub->idx;
289:   return(0);
290: }

292: /*@C
293:    ISBlockRestoreIndices - Restores the indices associated with each block.

295:    Not Collective

297:    Input Parameter:
298: .  is - the index set

300:    Output Parameter:
301: .  idx - the integer indices

303:    Level: intermediate

305:    Concepts: IS^block
306:    Concepts: index sets^getting indices
307:    Concepts: index sets^block

309: .seealso: ISRestoreIndices(), ISBlockGetIndices()
310: @*/
311: int ISBlockRestoreIndices(IS is,int *idx[])
312: {
316:   if (is->type != IS_BLOCK) SETERRQ(PETSC_ERR_ARG_WRONG,"Not a block index set");
317:   return(0);
318: }

320: /*@
321:    ISBlockGetBlockSize - Returns the number of elements in a block.

323:    Not Collective

325:    Input Parameter:
326: .  is - the index set

328:    Output Parameter:
329: .  size - the number of elements in a block

331:    Level: intermediate

333:    Concepts: IS^block size
334:    Concepts: index sets^block size

336: .seealso: ISBlockGetSize(), ISGetSize(), ISBlock(), ISCreateBlock()
337: @*/
338: int ISBlockGetBlockSize(IS is,int *size)
339: {
340:   IS_Block *sub;

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

347:   sub = (IS_Block *)is->data;
348:   *size = sub->bs;
349:   return(0);
350: }

352: /*@C
353:    ISBlock - Checks whether an index set is blocked.

355:    Not Collective

357:    Input Parameter:
358: .  is - the index set

360:    Output Parameter:
361: .  flag - PETSC_TRUE if a block index set, else PETSC_FALSE

363:    Level: intermediate

365:    Concepts: IS^block
366:    Concepts: index sets^block

368: .seealso: ISBlockGetSize(), ISGetSize(), ISBlockGetBlockSize(), ISCreateBlock()
369: @*/
370: int ISBlock(IS is,PetscTruth *flag)
371: {
375:   if (is->type != IS_BLOCK) *flag = PETSC_FALSE;
376:   else                          *flag = PETSC_TRUE;
377:   return(0);
378: }

380: /*@
381:    ISBlockGetSize - Returns the number of blocks in the index set.

383:    Not Collective

385:    Input Parameter:
386: .  is - the index set

388:    Output Parameter:
389: .  size - the number of blocks

391:    Level: intermediate

393:    Concepts: IS^block sizes
394:    Concepts: index sets^block sizes

396: .seealso: ISBlockGetBlockSize(), ISGetSize(), ISBlock(), ISCreateBlock()
397: @*/
398: int ISBlockGetSize(IS is,int *size)
399: {
400:   IS_Block *sub;

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

407:   sub = (IS_Block *)is->data;
408:   *size = sub->n;
409:   return(0);
410: }