Actual source code: stride.c

petsc-3.4.5 2014-06-29
  2: /*
  3:        Index sets of evenly space integers, defined by a
  4:     start, stride and length.
  5: */
  6: #include <petsc-private/isimpl.h>             /*I   "petscis.h"   I*/
  7: #include <petscvec.h>
  8: #include <petscviewer.h>

 10: typedef struct {
 11:   PetscInt N,n,first,step;
 12: } IS_Stride;

 16: PetscErrorCode ISIdentity_Stride(IS is,PetscBool  *ident)
 17: {
 18:   IS_Stride *is_stride = (IS_Stride*)is->data;

 21:   is->isidentity = PETSC_FALSE;
 22:   *ident         = PETSC_FALSE;
 23:   if (is_stride->first != 0) return(0);
 24:   if (is_stride->step  != 1) return(0);
 25:   *ident         = PETSC_TRUE;
 26:   is->isidentity = PETSC_TRUE;
 27:   return(0);
 28: }

 32: static PetscErrorCode ISCopy_Stride(IS is,IS isy)
 33: {
 34:   IS_Stride      *is_stride = (IS_Stride*)is->data,*isy_stride = (IS_Stride*)isy->data;

 38:   PetscMemcpy(isy_stride,is_stride,sizeof(IS_Stride));
 39:   return(0);
 40: }

 44: PetscErrorCode ISDuplicate_Stride(IS is,IS *newIS)
 45: {
 47:   IS_Stride      *sub = (IS_Stride*)is->data;

 50:   ISCreateStride(PetscObjectComm((PetscObject)is),sub->n,sub->first,sub->step,newIS);
 51:   return(0);
 52: }

 56: PetscErrorCode ISInvertPermutation_Stride(IS is,PetscInt nlocal,IS *perm)
 57: {
 58:   IS_Stride      *isstride = (IS_Stride*)is->data;

 62:   if (is->isidentity) {
 63:     ISCreateStride(PETSC_COMM_SELF,isstride->n,0,1,perm);
 64:   } else {
 65:     IS             tmp;
 66:     const PetscInt *indices,n = isstride->n;
 67:     ISGetIndices(is,&indices);
 68:     ISCreateGeneral(PetscObjectComm((PetscObject)is),n,indices,PETSC_COPY_VALUES,&tmp);
 69:     ISSetPermutation(tmp);
 70:     ISRestoreIndices(is,&indices);
 71:     ISInvertPermutation(tmp,nlocal,perm);
 72:     ISDestroy(&tmp);
 73:   }
 74:   return(0);
 75: }

 79: /*@
 80:    ISStrideGetInfo - Returns the first index in a stride index set and
 81:    the stride width.

 83:    Not Collective

 85:    Input Parameter:
 86: .  is - the index set

 88:    Output Parameters:
 89: .  first - the first index
 90: .  step - the stride width

 92:    Level: intermediate

 94:    Notes:
 95:    Returns info on stride index set. This is a pseudo-public function that
 96:    should not be needed by most users.

 98:    Concepts: index sets^getting information
 99:    Concepts: IS^getting information

101: .seealso: ISCreateStride(), ISGetSize()
102: @*/
103: PetscErrorCode  ISStrideGetInfo(IS is,PetscInt *first,PetscInt *step)
104: {
105:   IS_Stride *sub;


112:   sub = (IS_Stride*)is->data;
113:   if (first) *first = sub->first;
114:   if (step)  *step  = sub->step;
115:   return(0);
116: }

120: PetscErrorCode ISDestroy_Stride(IS is)
121: {

125:   PetscObjectComposeFunction((PetscObject)is,"ISStrideSetStride_C",0);
126:   PetscFree(is->data);
127:   return(0);
128: }

132: PetscErrorCode  ISToGeneral_Stride(IS inis)
133: {
135:   const PetscInt *idx;
136:   PetscInt       n;

139:   ISGetLocalSize(inis,&n);
140:   ISGetIndices(inis,&idx);
141:   ISSetType(inis,ISGENERAL);
142:   ISGeneralSetIndices(inis,n,idx,PETSC_OWN_POINTER);
143:   return(0);
144: }


147: /*
148:      Returns a legitimate index memory even if
149:    the stride index set is empty.
150: */
153: PetscErrorCode ISGetIndices_Stride(IS in,const PetscInt *idx[])
154: {
155:   IS_Stride      *sub = (IS_Stride*)in->data;
157:   PetscInt       i,**dx = (PetscInt**)idx;

160:   PetscMalloc(sub->n*sizeof(PetscInt),idx);
161:   if (sub->n) {
162:     (*dx)[0] = sub->first;
163:     for (i=1; i<sub->n; i++) (*dx)[i] = (*dx)[i-1] + sub->step;
164:   }
165:   return(0);
166: }

170: PetscErrorCode ISRestoreIndices_Stride(IS in,const PetscInt *idx[])
171: {

175:   PetscFree(*(void**)idx);
176:   return(0);
177: }

181: PetscErrorCode ISGetSize_Stride(IS is,PetscInt *size)
182: {
183:   IS_Stride *sub = (IS_Stride*)is->data;

186:   *size = sub->N;
187:   return(0);
188: }

192: PetscErrorCode ISGetLocalSize_Stride(IS is,PetscInt *size)
193: {
194:   IS_Stride *sub = (IS_Stride*)is->data;

197:   *size = sub->n;
198:   return(0);
199: }

203: PetscErrorCode ISView_Stride(IS is,PetscViewer viewer)
204: {
205:   IS_Stride      *sub = (IS_Stride*)is->data;
206:   PetscInt       i,n = sub->n;
207:   PetscMPIInt    rank,size;
208:   PetscBool      iascii;

212:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
213:   if (iascii) {
214:     MPI_Comm_rank(PetscObjectComm((PetscObject)is),&rank);
215:     MPI_Comm_size(PetscObjectComm((PetscObject)is),&size);
216:     if (size == 1) {
217:       if (is->isperm) {
218:         PetscViewerASCIIPrintf(viewer,"Index set is permutation\n");
219:       }
220:       PetscViewerASCIIPrintf(viewer,"Number of indices in (stride) set %D\n",n);
221:       for (i=0; i<n; i++) {
222:         PetscViewerASCIIPrintf(viewer,"%D %D\n",i,sub->first + i*sub->step);
223:       }
224:       PetscViewerFlush(viewer);
225:     } else {
226:       PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);
227:       if (is->isperm) {
228:         PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Index set is permutation\n",rank);
229:       }
230:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of indices in (stride) set %D\n",rank,n);
231:       for (i=0; i<n; i++) {
232:         PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D %D\n",rank,i,sub->first + i*sub->step);
233:       }
234:       PetscViewerFlush(viewer);
235:       PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);
236:     }
237:   }
238:   return(0);
239: }

243: PetscErrorCode ISSort_Stride(IS is)
244: {
245:   IS_Stride *sub = (IS_Stride*)is->data;

248:   if (sub->step >= 0) return(0);
249:   sub->first += (sub->n - 1)*sub->step;
250:   sub->step  *= -1;
251:   return(0);
252: }

256: PetscErrorCode ISSorted_Stride(IS is,PetscBool * flg)
257: {
258:   IS_Stride *sub = (IS_Stride*)is->data;

261:   if (sub->step >= 0) *flg = PETSC_TRUE;
262:   else *flg = PETSC_FALSE;
263:   return(0);
264: }

268: static PetscErrorCode ISOnComm_Stride(IS is,MPI_Comm comm,PetscCopyMode mode,IS *newis)
269: {
271:   IS_Stride      *sub = (IS_Stride*)is->data;

274:   ISCreateStride(comm,sub->n,sub->first,sub->step,newis);
275:   return(0);
276: }

280: static PetscErrorCode ISSetBlockSize_Stride(IS is,PetscInt bs)
281: {
282:   IS_Stride *sub = (IS_Stride*)is->data;

285:   if (sub->step != 1 && bs != 1) SETERRQ2(PetscObjectComm((PetscObject)is),PETSC_ERR_ARG_SIZ,"ISSTRIDE has stride %D, cannot be blocked of size %D",sub->step,bs);
286:   is->bs = bs;
287:   return(0);
288: }

292: static PetscErrorCode ISContiguousLocal_Stride(IS is,PetscInt gstart,PetscInt gend,PetscInt *start,PetscBool *contig)
293: {
294:   IS_Stride *sub = (IS_Stride*)is->data;

297:   if (sub->step == 1 && sub->first >= gstart && sub->first+sub->n <= gend) {
298:     *start  = sub->first - gstart;
299:     *contig = PETSC_TRUE;
300:   } else {
301:     *start  = -1;
302:     *contig = PETSC_FALSE;
303:   }
304:   return(0);
305: }


308: static struct _ISOps myops = { ISGetSize_Stride,
309:                                ISGetLocalSize_Stride,
310:                                ISGetIndices_Stride,
311:                                ISRestoreIndices_Stride,
312:                                ISInvertPermutation_Stride,
313:                                ISSort_Stride,
314:                                ISSorted_Stride,
315:                                ISDuplicate_Stride,
316:                                ISDestroy_Stride,
317:                                ISView_Stride,
318:                                ISIdentity_Stride,
319:                                ISCopy_Stride,
320:                                ISToGeneral_Stride,
321:                                ISOnComm_Stride,
322:                                ISSetBlockSize_Stride,
323:                                ISContiguousLocal_Stride};


328: /*@
329:    ISStrideSetStride - Sets the stride information for a stride index set.

331:    Collective on IS

333:    Input Parameters:
334: +  is - the index set
335: .  n - the length of the locally owned portion of the index set
336: .  first - the first element of the locally owned portion of the index set
337: -  step - the change to the next index

339:    Level: beginner

341:   Concepts: IS^stride
342:   Concepts: index sets^stride
343:   Concepts: stride^index set

345: .seealso: ISCreateGeneral(), ISCreateBlock(), ISAllGather()
346: @*/
347: PetscErrorCode  ISStrideSetStride(IS is,PetscInt n,PetscInt first,PetscInt step)
348: {

352:   if (n < 0) SETERRQ1(PetscObjectComm((PetscObject)is), PETSC_ERR_ARG_OUTOFRANGE, "Negative length %d not valid", n);
353:   PetscUseMethod(is,"ISStrideSetStride_C",(IS,PetscInt,PetscInt,PetscInt),(is,n,first,step));
354:   return(0);
355: }

359: PetscErrorCode  ISStrideSetStride_Stride(IS is,PetscInt n,PetscInt first,PetscInt step)
360: {
362:   PetscInt       min,max;
363:   IS_Stride      *sub = (IS_Stride*)is->data;

366:   sub->n     = n;
367:   MPI_Allreduce(&n,&sub->N,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)is));
368:   sub->first = first;
369:   sub->step  = step;
370:   if (step > 0) {min = first; max = first + step*(n-1);}
371:   else          {max = first; min = first + step*(n-1);}

373:   is->min  = n > 0 ? min : -1;
374:   is->max  = n > 0 ? max : -1;
375:   is->data = (void*)sub;

377:   if ((!first && step == 1) || (first == max && step == -1 && !min)) is->isperm = PETSC_TRUE;
378:   else is->isperm = PETSC_FALSE;
379:   return(0);
380: }

384: /*@
385:    ISCreateStride - Creates a data structure for an index set
386:    containing a list of evenly spaced integers.

388:    Collective on MPI_Comm

390:    Input Parameters:
391: +  comm - the MPI communicator
392: .  n - the length of the locally owned portion of the index set
393: .  first - the first element of the locally owned portion of the index set
394: -  step - the change to the next index

396:    Output Parameter:
397: .  is - the new index set

399:    Notes:
400:    When the communicator is not MPI_COMM_SELF, the operations on IS are NOT
401:    conceptually the same as MPI_Group operations. The IS are the
402:    distributed sets of indices and thus certain operations on them are collective.

404:    Level: beginner

406:   Concepts: IS^stride
407:   Concepts: index sets^stride
408:   Concepts: stride^index set

410: .seealso: ISCreateGeneral(), ISCreateBlock(), ISAllGather()
411: @*/
412: PetscErrorCode  ISCreateStride(MPI_Comm comm,PetscInt n,PetscInt first,PetscInt step,IS *is)
413: {

417:   ISCreate(comm,is);
418:   ISSetType(*is,ISSTRIDE);
419:   ISStrideSetStride(*is,n,first,step);
420:   return(0);
421: }

425: PETSC_EXTERN PetscErrorCode ISCreate_Stride(IS is)
426: {
428:   IS_Stride      *sub;

431:   PetscMemcpy(is->ops,&myops,sizeof(myops));
432:   PetscNewLog(is,IS_Stride,&sub);
433:   is->bs   = 1;
434:   is->data = sub;
435:   PetscObjectComposeFunction((PetscObject)is,"ISStrideSetStride_C",ISStrideSetStride_Stride);
436:   return(0);
437: }