Actual source code: stride.c
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 first,step;
11: } IS_Stride;
13: static PetscErrorCode ISCopy_Stride(IS is,IS isy)
14: {
15: IS_Stride *is_stride = (IS_Stride*)is->data,*isy_stride = (IS_Stride*)isy->data;
19: PetscMemcpy(isy_stride,is_stride,sizeof(IS_Stride));
20: return(0);
21: }
23: PetscErrorCode ISDuplicate_Stride(IS is,IS *newIS)
24: {
26: IS_Stride *sub = (IS_Stride*)is->data;
29: ISCreateStride(PetscObjectComm((PetscObject)is),is->map->n,sub->first,sub->step,newIS);
30: return(0);
31: }
33: PetscErrorCode ISInvertPermutation_Stride(IS is,PetscInt nlocal,IS *perm)
34: {
35: PetscBool isident;
39: ISGetInfo(is,IS_IDENTITY,IS_GLOBAL,PETSC_TRUE,&isident);
40: if (isident) {
41: PetscInt rStart, rEnd;
43: PetscLayoutGetRange(is->map, &rStart, &rEnd);
44: ISCreateStride(PETSC_COMM_SELF,PetscMax(rEnd - rStart, 0),rStart,1,perm);
45: } else {
46: IS tmp;
47: const PetscInt *indices,n = is->map->n;
49: ISGetIndices(is,&indices);
50: ISCreateGeneral(PetscObjectComm((PetscObject)is),n,indices,PETSC_COPY_VALUES,&tmp);
51: ISSetPermutation(tmp);
52: ISRestoreIndices(is,&indices);
53: ISInvertPermutation(tmp,nlocal,perm);
54: ISDestroy(&tmp);
55: }
56: return(0);
57: }
59: /*@
60: ISStrideGetInfo - Returns the first index in a stride index set and the stride width.
62: Not Collective
64: Input Parameter:
65: . is - the index set
67: Output Parameters:
68: + first - the first index
69: - step - the stride width
71: Level: intermediate
73: Notes:
74: Returns info on stride index set. This is a pseudo-public function that
75: should not be needed by most users.
77: .seealso: ISCreateStride(), ISGetSize(), ISSTRIDE
78: @*/
79: PetscErrorCode ISStrideGetInfo(IS is,PetscInt *first,PetscInt *step)
80: {
81: IS_Stride *sub;
82: PetscBool flg;
89: PetscObjectTypeCompare((PetscObject)is,ISSTRIDE,&flg);
90: if (!flg) SETERRQ(PetscObjectComm((PetscObject)is),PETSC_ERR_ARG_WRONG,"IS must be of type ISSTRIDE");
92: sub = (IS_Stride*)is->data;
93: if (first) *first = sub->first;
94: if (step) *step = sub->step;
95: return(0);
96: }
98: PetscErrorCode ISDestroy_Stride(IS is)
99: {
103: PetscObjectComposeFunction((PetscObject)is,"ISStrideSetStride_C",NULL);
104: PetscFree(is->data);
105: return(0);
106: }
108: PetscErrorCode ISToGeneral_Stride(IS inis)
109: {
111: const PetscInt *idx;
112: PetscInt n;
115: ISGetLocalSize(inis,&n);
116: ISGetIndices(inis,&idx);
117: ISSetType(inis,ISGENERAL);
118: ISGeneralSetIndices(inis,n,idx,PETSC_OWN_POINTER);
119: return(0);
120: }
122: PetscErrorCode ISLocate_Stride(IS is,PetscInt key,PetscInt *location)
123: {
124: IS_Stride *sub = (IS_Stride*)is->data;
125: PetscInt rem, step;
128: *location = -1;
129: step = sub->step;
130: key -= sub->first;
131: rem = key / step;
132: if ((rem < is->map->n) && !(key % step)) {
133: *location = rem;
134: }
135: return(0);
136: }
138: /*
139: Returns a legitimate index memory even if
140: the stride index set is empty.
141: */
142: PetscErrorCode ISGetIndices_Stride(IS is,const PetscInt *idx[])
143: {
144: IS_Stride *sub = (IS_Stride*)is->data;
146: PetscInt i,**dx = (PetscInt**)idx;
149: PetscMalloc1(is->map->n,(PetscInt**)idx);
150: if (is->map->n) {
151: (*dx)[0] = sub->first;
152: for (i=1; i<is->map->n; i++) (*dx)[i] = (*dx)[i-1] + sub->step;
153: }
154: return(0);
155: }
157: PetscErrorCode ISRestoreIndices_Stride(IS in,const PetscInt *idx[])
158: {
162: PetscFree(*(void**)idx);
163: return(0);
164: }
166: PetscErrorCode ISView_Stride(IS is,PetscViewer viewer)
167: {
168: IS_Stride *sub = (IS_Stride*)is->data;
169: PetscInt i,n = is->map->n;
170: PetscMPIInt rank,size;
171: PetscBool iascii,ibinary;
172: PetscViewerFormat fmt;
173: PetscErrorCode ierr;
176: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
177: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&ibinary);
178: if (iascii) {
179: PetscBool matl, isperm;
181: MPI_Comm_rank(PetscObjectComm((PetscObject)is),&rank);
182: MPI_Comm_size(PetscObjectComm((PetscObject)is),&size);
183: PetscViewerGetFormat(viewer,&fmt);
184: matl = (PetscBool)(fmt == PETSC_VIEWER_ASCII_MATLAB);
185: ISGetInfo(is,IS_PERMUTATION,IS_GLOBAL,PETSC_FALSE,&isperm);
186: if (isperm && !matl) {PetscViewerASCIIPrintf(viewer,"Index set is permutation\n");}
187: if (size == 1) {
188: if (matl) {
189: const char* name;
191: PetscObjectGetName((PetscObject)is,&name);
192: PetscViewerASCIIPrintf(viewer,"%s = [%D : %D : %D];\n",name,sub->first+1,sub->step,sub->first + sub->step*(n-1)+1);
193: } else {
194: PetscViewerASCIIPrintf(viewer,"Number of indices in (stride) set %D\n",n);
195: for (i=0; i<n; i++) {
196: PetscViewerASCIIPrintf(viewer,"%D %D\n",i,sub->first + i*sub->step);
197: }
198: }
199: PetscViewerFlush(viewer);
200: } else {
201: PetscViewerASCIIPushSynchronized(viewer);
202: if (matl) {
203: const char* name;
205: PetscObjectGetName((PetscObject)is,&name);
206: PetscViewerASCIISynchronizedPrintf(viewer,"%s_%d = [%D : %D : %D];\n",name,rank,sub->first+1,sub->step,sub->first + sub->step*(n-1)+1);
207: } else {
208: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of indices in (stride) set %D\n",rank,n);
209: for (i=0; i<n; i++) {
210: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D %D\n",rank,i,sub->first + i*sub->step);
211: }
212: }
213: PetscViewerFlush(viewer);
214: PetscViewerASCIIPopSynchronized(viewer);
215: }
216: } else if (ibinary) {
217: ISView_Binary(is,viewer);
218: }
219: return(0);
220: }
222: PetscErrorCode ISSort_Stride(IS is)
223: {
224: IS_Stride *sub = (IS_Stride*)is->data;
227: if (sub->step >= 0) return(0);
228: sub->first += (is->map->n - 1)*sub->step;
229: sub->step *= -1;
230: return(0);
231: }
233: PetscErrorCode ISSorted_Stride(IS is,PetscBool * flg)
234: {
235: IS_Stride *sub = (IS_Stride*)is->data;
238: if (sub->step >= 0) *flg = PETSC_TRUE;
239: else *flg = PETSC_FALSE;
240: return(0);
241: }
243: static PetscErrorCode ISUniqueLocal_Stride(IS is, PetscBool *flg)
244: {
245: IS_Stride *sub = (IS_Stride*)is->data;
248: if (!(is->map->n) || sub->step != 0) *flg = PETSC_TRUE;
249: else *flg = PETSC_FALSE;
250: return(0);
251: }
253: static PetscErrorCode ISPermutationLocal_Stride(IS is, PetscBool *flg)
254: {
255: IS_Stride *sub = (IS_Stride*)is->data;
258: if (!(is->map->n) || (PetscAbsInt(sub->step) == 1 && is->min == 0)) *flg = PETSC_TRUE;
259: else *flg = PETSC_FALSE;
260: return(0);
261: }
263: static PetscErrorCode ISIntervalLocal_Stride(IS is, PetscBool *flg)
264: {
265: IS_Stride *sub = (IS_Stride*)is->data;
268: if (!(is->map->n) || sub->step == 1) *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,is->map->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+is->map->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: }
309: static struct _ISOps myops = { ISGetIndices_Stride,
310: ISRestoreIndices_Stride,
311: ISInvertPermutation_Stride,
312: ISSort_Stride,
313: ISSort_Stride,
314: ISSorted_Stride,
315: ISDuplicate_Stride,
316: ISDestroy_Stride,
317: ISView_Stride,
318: ISLoad_Default,
319: ISCopy_Stride,
320: ISToGeneral_Stride,
321: ISOnComm_Stride,
322: ISSetBlockSize_Stride,
323: ISContiguousLocal_Stride,
324: ISLocate_Stride,
325: ISSorted_Stride,
326: NULL,
327: ISUniqueLocal_Stride,
328: NULL,
329: ISPermutationLocal_Stride,
330: NULL,
331: ISIntervalLocal_Stride,
332: NULL};
334: /*@
335: ISStrideSetStride - Sets the stride information for a stride index set.
337: Collective on IS
339: Input Parameters:
340: + is - the index set
341: . n - the length of the locally owned portion of the index set
342: . first - the first element of the locally owned portion of the index set
343: - step - the change to the next index
345: Level: beginner
347: .seealso: ISCreateGeneral(), ISCreateBlock(), ISAllGather(), ISSTRIDE, ISCreateStride(), ISStrideGetInfo()
348: @*/
349: PetscErrorCode ISStrideSetStride(IS is,PetscInt n,PetscInt first,PetscInt step)
350: {
354: if (n < 0) SETERRQ1(PetscObjectComm((PetscObject)is), PETSC_ERR_ARG_OUTOFRANGE, "Negative length %D not valid", n);
355: ISClearInfoCache(is,PETSC_FALSE);
356: PetscUseMethod(is,"ISStrideSetStride_C",(IS,PetscInt,PetscInt,PetscInt),(is,n,first,step));
357: return(0);
358: }
360: PetscErrorCode ISStrideSetStride_Stride(IS is,PetscInt n,PetscInt first,PetscInt step)
361: {
363: PetscInt min,max;
364: IS_Stride *sub = (IS_Stride*)is->data;
365: PetscLayout map;
368: PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)is),n,is->map->N,is->map->bs,&map);
369: PetscLayoutDestroy(&is->map);
370: is->map = map;
372: sub->first = first;
373: sub->step = step;
374: if (step > 0) {min = first; max = first + step*(n-1);}
375: else {max = first; min = first + step*(n-1);}
377: is->min = n > 0 ? min : PETSC_MAX_INT;
378: is->max = n > 0 ? max : PETSC_MIN_INT;
379: is->data = (void*)sub;
380: return(0);
381: }
383: /*@
384: ISCreateStride - Creates a data structure for an index set
385: containing a list of evenly spaced integers.
387: Collective
389: Input Parameters:
390: + comm - the MPI communicator
391: . n - the length of the locally owned portion of the index set
392: . first - the first element of the locally owned portion of the index set
393: - step - the change to the next index
395: Output Parameter:
396: . is - the new index set
398: Notes:
399: When the communicator is not MPI_COMM_SELF, the operations on IS are NOT
400: conceptually the same as MPI_Group operations. The IS are the
401: distributed sets of indices and thus certain operations on them are collective.
403: Level: beginner
405: .seealso: ISCreateGeneral(), ISCreateBlock(), ISAllGather(), ISSTRIDE
406: @*/
407: PetscErrorCode ISCreateStride(MPI_Comm comm,PetscInt n,PetscInt first,PetscInt step,IS *is)
408: {
412: ISCreate(comm,is);
413: ISSetType(*is,ISSTRIDE);
414: ISStrideSetStride(*is,n,first,step);
415: return(0);
416: }
418: PETSC_EXTERN PetscErrorCode ISCreate_Stride(IS is)
419: {
421: IS_Stride *sub;
424: PetscNewLog(is,&sub);
425: is->data = (void *) sub;
426: PetscMemcpy(is->ops,&myops,sizeof(myops));
427: PetscObjectComposeFunction((PetscObject)is,"ISStrideSetStride_C",ISStrideSetStride_Stride);
428: return(0);
429: }