Actual source code: stride.c

petsc-master 2019-07-21
Report Typos and Errors

  2: /*
  3:        Index sets of evenly space integers, defined by a
  4:     start, stride and length.
  5: */
  6:  #include <petsc/private/isimpl.h>
  7:  #include <petscviewer.h>

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

 13: PetscErrorCode ISIdentity_Stride(IS is,PetscBool  *ident)
 14: {
 15:   IS_Stride *is_stride = (IS_Stride*)is->data;

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

 27: static PetscErrorCode ISCopy_Stride(IS is,IS isy)
 28: {
 29:   IS_Stride      *is_stride = (IS_Stride*)is->data,*isy_stride = (IS_Stride*)isy->data;

 33:   PetscMemcpy(isy_stride,is_stride,sizeof(IS_Stride));
 34:   return(0);
 35: }

 37: PetscErrorCode ISDuplicate_Stride(IS is,IS *newIS)
 38: {
 40:   IS_Stride      *sub = (IS_Stride*)is->data;

 43:   ISCreateStride(PetscObjectComm((PetscObject)is),sub->n,sub->first,sub->step,newIS);
 44:   return(0);
 45: }

 47: PetscErrorCode ISInvertPermutation_Stride(IS is,PetscInt nlocal,IS *perm)
 48: {
 49:   IS_Stride      *isstride = (IS_Stride*)is->data;

 53:   if (is->isidentity) {
 54:     ISCreateStride(PETSC_COMM_SELF,isstride->n,0,1,perm);
 55:   } else {
 56:     IS             tmp;
 57:     const PetscInt *indices,n = isstride->n;
 58:     ISGetIndices(is,&indices);
 59:     ISCreateGeneral(PetscObjectComm((PetscObject)is),n,indices,PETSC_COPY_VALUES,&tmp);
 60:     ISSetPermutation(tmp);
 61:     ISRestoreIndices(is,&indices);
 62:     ISInvertPermutation(tmp,nlocal,perm);
 63:     ISDestroy(&tmp);
 64:   }
 65:   return(0);
 66: }

 68: /*@
 69:    ISStrideGetInfo - Returns the first index in a stride index set and
 70:    the stride width.

 72:    Not Collective

 74:    Input Parameter:
 75: .  is - the index set

 77:    Output Parameters:
 78: +  first - the first index
 79: -  step - the stride width

 81:    Level: intermediate

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


 88: .seealso: ISCreateStride(), ISGetSize()
 89: @*/
 90: PetscErrorCode  ISStrideGetInfo(IS is,PetscInt *first,PetscInt *step)
 91: {
 92:   IS_Stride      *sub;
 93:   PetscBool      flg;

100:   PetscObjectTypeCompare((PetscObject)is,ISSTRIDE,&flg);
101:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)is),PETSC_ERR_ARG_WRONG,"IS must be of type ISSTRIDE");

103:   sub = (IS_Stride*)is->data;
104:   if (first) *first = sub->first;
105:   if (step)  *step  = sub->step;
106:   return(0);
107: }

109: PetscErrorCode ISDestroy_Stride(IS is)
110: {

114:   PetscObjectComposeFunction((PetscObject)is,"ISStrideSetStride_C",NULL);
115:   PetscFree(is->data);
116:   return(0);
117: }

119: PetscErrorCode  ISToGeneral_Stride(IS inis)
120: {
122:   const PetscInt *idx;
123:   PetscInt       n;

126:   ISGetLocalSize(inis,&n);
127:   ISGetIndices(inis,&idx);
128:   ISSetType(inis,ISGENERAL);
129:   ISGeneralSetIndices(inis,n,idx,PETSC_OWN_POINTER);
130:   return(0);
131: }

133: PetscErrorCode ISLocate_Stride(IS is,PetscInt key,PetscInt *location)
134: {
135:   IS_Stride      *sub = (IS_Stride*)is->data;
136:   PetscInt       rem, step;

139:   *location = -1;
140:   step      = sub->step;
141:   key      -= sub->first;
142:   rem       = key / step;
143:   if ((rem < sub->n) && !(key % step)) {
144:     *location = rem;
145:   }
146:   return(0);
147: }

149: /*
150:      Returns a legitimate index memory even if
151:    the stride index set is empty.
152: */
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:   PetscMalloc1(sub->n,(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: }

168: PetscErrorCode ISRestoreIndices_Stride(IS in,const PetscInt *idx[])
169: {

173:   PetscFree(*(void**)idx);
174:   return(0);
175: }

177: PetscErrorCode ISGetSize_Stride(IS is,PetscInt *size)
178: {
179:   IS_Stride *sub = (IS_Stride*)is->data;

182:   *size = sub->N;
183:   return(0);
184: }

186: PetscErrorCode ISGetLocalSize_Stride(IS is,PetscInt *size)
187: {
188:   IS_Stride *sub = (IS_Stride*)is->data;

191:   *size = sub->n;
192:   return(0);
193: }

195: PetscErrorCode ISView_Stride(IS is,PetscViewer viewer)
196: {
197:   IS_Stride         *sub = (IS_Stride*)is->data;
198:   PetscInt          i,n = sub->n;
199:   PetscMPIInt       rank,size;
200:   PetscBool         iascii;
201:   PetscViewerFormat fmt;
202:   PetscErrorCode    ierr;

205:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
206:   if (iascii) {
207:     PetscBool matl;

209:     MPI_Comm_rank(PetscObjectComm((PetscObject)is),&rank);
210:     MPI_Comm_size(PetscObjectComm((PetscObject)is),&size);
211:     PetscViewerGetFormat(viewer,&fmt);
212:     matl = (PetscBool)(fmt == PETSC_VIEWER_ASCII_MATLAB);
213:     if (size == 1) {
214:       if (matl) {
215:         const char* name;

217:         PetscObjectGetName((PetscObject)is,&name);
218:         PetscViewerASCIIPrintf(viewer,"%s = [%D : %D : %D];\n",name,sub->first+1,sub->step,sub->first + sub->step*(n-1)+1);
219:       } else {
220:         if (is->isperm) {
221:           PetscViewerASCIIPrintf(viewer,"Index set is permutation\n");
222:         }
223:         PetscViewerASCIIPrintf(viewer,"Number of indices in (stride) set %D\n",n);
224:         for (i=0; i<n; i++) {
225:           PetscViewerASCIIPrintf(viewer,"%D %D\n",i,sub->first + i*sub->step);
226:         }
227:       }
228:       PetscViewerFlush(viewer);
229:     } else {
230:       PetscViewerASCIIPushSynchronized(viewer);
231:       if (matl) {
232:         const char* name;

234:         PetscObjectGetName((PetscObject)is,&name);
235:         PetscViewerASCIISynchronizedPrintf(viewer,"%s_%d = [%D : %D : %D];\n",name,rank,sub->first+1,sub->step,sub->first + sub->step*(n-1)+1);
236:       } else {
237:         if (is->isperm) {
238:           PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Index set is permutation\n",rank);
239:         }
240:         PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of indices in (stride) set %D\n",rank,n);
241:         for (i=0; i<n; i++) {
242:           PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D %D\n",rank,i,sub->first + i*sub->step);
243:         }
244:       }
245:       PetscViewerFlush(viewer);
246:       PetscViewerASCIIPopSynchronized(viewer);
247:     }
248:   }
249:   return(0);
250: }

252: PetscErrorCode ISSort_Stride(IS is)
253: {
254:   IS_Stride *sub = (IS_Stride*)is->data;

257:   if (sub->step >= 0) return(0);
258:   sub->first += (sub->n - 1)*sub->step;
259:   sub->step  *= -1;
260:   return(0);
261: }

263: PetscErrorCode ISSorted_Stride(IS is,PetscBool * flg)
264: {
265:   IS_Stride *sub = (IS_Stride*)is->data;

268:   if (sub->step >= 0) *flg = PETSC_TRUE;
269:   else *flg = PETSC_FALSE;
270:   return(0);
271: }

273: static PetscErrorCode ISOnComm_Stride(IS is,MPI_Comm comm,PetscCopyMode mode,IS *newis)
274: {
276:   IS_Stride      *sub = (IS_Stride*)is->data;

279:   ISCreateStride(comm,sub->n,sub->first,sub->step,newis);
280:   return(0);
281: }

283: static PetscErrorCode ISSetBlockSize_Stride(IS is,PetscInt bs)
284: {
285:   IS_Stride     *sub = (IS_Stride*)is->data;

289:   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);
290:   PetscLayoutSetBlockSize(is->map, bs);
291:   return(0);
292: }

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

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


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


331: /*@
332:    ISStrideSetStride - Sets the stride information for a stride index set.

334:    Collective on IS

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

342:    Level: beginner

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

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

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

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

370:   is->min  = n > 0 ? min : PETSC_MAX_INT;
371:   is->max  = n > 0 ? max : PETSC_MIN_INT;
372:   is->data = (void*)sub;

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

380: /*@
381:    ISCreateStride - Creates a data structure for an index set
382:    containing a list of evenly spaced integers.

384:    Collective

386:    Input Parameters:
387: +  comm - the MPI communicator
388: .  n - the length of the locally owned portion of the index set
389: .  first - the first element of the locally owned portion of the index set
390: -  step - the change to the next index

392:    Output Parameter:
393: .  is - the new index set

395:    Notes:
396:    When the communicator is not MPI_COMM_SELF, the operations on IS are NOT
397:    conceptually the same as MPI_Group operations. The IS are the
398:    distributed sets of indices and thus certain operations on them are collective.

400:    Level: beginner

402: .seealso: ISCreateGeneral(), ISCreateBlock(), ISAllGather()
403: @*/
404: PetscErrorCode  ISCreateStride(MPI_Comm comm,PetscInt n,PetscInt first,PetscInt step,IS *is)
405: {

409:   ISCreate(comm,is);
410:   ISSetType(*is,ISSTRIDE);
411:   ISStrideSetStride(*is,n,first,step);
412:   return(0);
413: }

415: PETSC_EXTERN PetscErrorCode ISCreate_Stride(IS is)
416: {
418:   IS_Stride      *sub;

421:   PetscNewLog(is,&sub);
422:   is->data = (void *) sub;
423:   PetscMemcpy(is->ops,&myops,sizeof(myops));
424:   PetscObjectComposeFunction((PetscObject)is,"ISStrideSetStride_C",ISStrideSetStride_Stride);
425:   return(0);
426: }