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: EXTERN int VecInitializePackage(char *);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

149:   *flg = sub->sorted;
150:   return(0);
151: }

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

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

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

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

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

196:    Collective on MPI_Comm

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

204:    Output Parameter:
205: .  is - the new index set

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

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

217:    Level: beginner

219:   Concepts: IS^block
220:   Concepts: index sets^block
221:   Concepts: block^index set

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

234:   if (n < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"length < 0");
236:   *is = PETSC_NULL;
237: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
238:   VecInitializePackage(PETSC_NULL);
239: #endif

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


268: /*@C
269:    ISBlockGetIndices - Gets the indices associated with each block.

271:    Not Collective

273:    Input Parameter:
274: .  is - the index set

276:    Output Parameter:
277: .  idx - the integer indices

279:    Level: intermediate

281:    Concepts: IS^block
282:    Concepts: index sets^getting indices
283:    Concepts: index sets^block

285: .seealso: ISGetIndices(), ISBlockRestoreIndices()
286: @*/
287: int ISBlockGetIndices(IS in,int *idx[])
288: {
289:   IS_Block *sub;

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

296:   sub = (IS_Block*)in->data;
297:   *idx = sub->idx;
298:   return(0);
299: }

301: /*@C
302:    ISBlockRestoreIndices - Restores the indices associated with each block.

304:    Not Collective

306:    Input Parameter:
307: .  is - the index set

309:    Output Parameter:
310: .  idx - the integer indices

312:    Level: intermediate

314:    Concepts: IS^block
315:    Concepts: index sets^getting indices
316:    Concepts: index sets^block

318: .seealso: ISRestoreIndices(), ISBlockGetIndices()
319: @*/
320: int ISBlockRestoreIndices(IS is,int *idx[])
321: {
325:   if (is->type != IS_BLOCK) SETERRQ(PETSC_ERR_ARG_WRONG,"Not a block index set");
326:   return(0);
327: }

329: /*@
330:    ISBlockGetBlockSize - Returns the number of elements in a block.

332:    Not Collective

334:    Input Parameter:
335: .  is - the index set

337:    Output Parameter:
338: .  size - the number of elements in a block

340:    Level: intermediate

342:    Concepts: IS^block size
343:    Concepts: index sets^block size

345: .seealso: ISBlockGetSize(), ISGetSize(), ISBlock(), ISCreateBlock()
346: @*/
347: int ISBlockGetBlockSize(IS is,int *size)
348: {
349:   IS_Block *sub;

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

356:   sub = (IS_Block *)is->data;
357:   *size = sub->bs;
358:   return(0);
359: }

361: /*@C
362:    ISBlock - Checks whether an index set is blocked.

364:    Not Collective

366:    Input Parameter:
367: .  is - the index set

369:    Output Parameter:
370: .  flag - PETSC_TRUE if a block index set, else PETSC_FALSE

372:    Level: intermediate

374:    Concepts: IS^block
375:    Concepts: index sets^block

377: .seealso: ISBlockGetSize(), ISGetSize(), ISBlockGetBlockSize(), ISCreateBlock()
378: @*/
379: int ISBlock(IS is,PetscTruth *flag)
380: {
384:   if (is->type != IS_BLOCK) *flag = PETSC_FALSE;
385:   else                          *flag = PETSC_TRUE;
386:   return(0);
387: }

389: /*@
390:    ISBlockGetSize - Returns the number of blocks in the index set.

392:    Not Collective

394:    Input Parameter:
395: .  is - the index set

397:    Output Parameter:
398: .  size - the number of blocks

400:    Level: intermediate

402:    Concepts: IS^block sizes
403:    Concepts: index sets^block sizes

405: .seealso: ISBlockGetBlockSize(), ISGetSize(), ISBlock(), ISCreateBlock()
406: @*/
407: int ISBlockGetSize(IS is,int *size)
408: {
409:   IS_Block *sub;

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

416:   sub = (IS_Block *)is->data;
417:   *size = sub->n;
418:   return(0);
419: }