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;
17: PetscMemcpy(isy_stride, is_stride, sizeof(IS_Stride));
18: return 0;
19: }
21: PetscErrorCode ISShift_Stride(IS is, PetscInt shift, IS isy)
22: {
23: IS_Stride *is_stride = (IS_Stride *)is->data, *isy_stride = (IS_Stride *)isy->data;
25: isy_stride->first = is_stride->first + shift;
26: isy_stride->step = is_stride->step;
27: return 0;
28: }
30: PetscErrorCode ISDuplicate_Stride(IS is, IS *newIS)
31: {
32: IS_Stride *sub = (IS_Stride *)is->data;
34: ISCreateStride(PetscObjectComm((PetscObject)is), is->map->n, sub->first, sub->step, newIS);
35: return 0;
36: }
38: PetscErrorCode ISInvertPermutation_Stride(IS is, PetscInt nlocal, IS *perm)
39: {
40: PetscBool isident;
42: ISGetInfo(is, IS_IDENTITY, IS_GLOBAL, PETSC_TRUE, &isident);
43: if (isident) {
44: PetscInt rStart, rEnd;
46: PetscLayoutGetRange(is->map, &rStart, &rEnd);
47: ISCreateStride(PETSC_COMM_SELF, PetscMax(rEnd - rStart, 0), rStart, 1, perm);
48: } else {
49: IS tmp;
50: const PetscInt *indices, n = is->map->n;
52: ISGetIndices(is, &indices);
53: ISCreateGeneral(PetscObjectComm((PetscObject)is), n, indices, PETSC_COPY_VALUES, &tmp);
54: ISSetPermutation(tmp);
55: ISRestoreIndices(is, &indices);
56: ISInvertPermutation(tmp, nlocal, perm);
57: ISDestroy(&tmp);
58: }
59: return 0;
60: }
62: /*@
63: ISStrideGetInfo - Returns the first index in a stride index set and the stride width.
65: Not Collective
67: Input Parameter:
68: . is - the index set
70: Output Parameters:
71: + first - the first index
72: - step - the stride width
74: Level: intermediate
76: Notes:
77: Returns info on stride index set. This is a pseudo-public function that
78: should not be needed by most users.
80: .seealso: `ISCreateStride()`, `ISGetSize()`, `ISSTRIDE`
81: @*/
82: PetscErrorCode ISStrideGetInfo(IS is, PetscInt *first, PetscInt *step)
83: {
84: IS_Stride *sub;
85: PetscBool flg;
90: PetscObjectTypeCompare((PetscObject)is, ISSTRIDE, &flg);
93: sub = (IS_Stride *)is->data;
94: if (first) *first = sub->first;
95: if (step) *step = sub->step;
96: return 0;
97: }
99: PetscErrorCode ISDestroy_Stride(IS is)
100: {
101: PetscObjectComposeFunction((PetscObject)is, "ISStrideSetStride_C", NULL);
102: PetscObjectComposeFunction((PetscObject)is, "ISShift_C", NULL);
103: PetscFree(is->data);
104: return 0;
105: }
107: PetscErrorCode ISToGeneral_Stride(IS inis)
108: {
109: const PetscInt *idx;
110: PetscInt n;
112: ISGetLocalSize(inis, &n);
113: ISGetIndices(inis, &idx);
114: ISSetType(inis, ISGENERAL);
115: ISGeneralSetIndices(inis, n, idx, PETSC_OWN_POINTER);
116: return 0;
117: }
119: PetscErrorCode ISLocate_Stride(IS is, PetscInt key, PetscInt *location)
120: {
121: IS_Stride *sub = (IS_Stride *)is->data;
122: PetscInt rem, step;
124: *location = -1;
125: step = sub->step;
126: key -= sub->first;
127: rem = key / step;
128: if ((rem < is->map->n) && !(key % step)) *location = rem;
129: return 0;
130: }
132: /*
133: Returns a legitimate index memory even if
134: the stride index set is empty.
135: */
136: PetscErrorCode ISGetIndices_Stride(IS is, const PetscInt *idx[])
137: {
138: IS_Stride *sub = (IS_Stride *)is->data;
139: PetscInt i, **dx = (PetscInt **)idx;
141: PetscMalloc1(is->map->n, (PetscInt **)idx);
142: if (is->map->n) {
143: (*dx)[0] = sub->first;
144: for (i = 1; i < is->map->n; i++) (*dx)[i] = (*dx)[i - 1] + sub->step;
145: }
146: return 0;
147: }
149: PetscErrorCode ISRestoreIndices_Stride(IS in, const PetscInt *idx[])
150: {
151: PetscFree(*(void **)idx);
152: return 0;
153: }
155: PetscErrorCode ISView_Stride(IS is, PetscViewer viewer)
156: {
157: IS_Stride *sub = (IS_Stride *)is->data;
158: PetscInt i, n = is->map->n;
159: PetscMPIInt rank, size;
160: PetscBool iascii, ibinary;
161: PetscViewerFormat fmt;
163: PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);
164: PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &ibinary);
165: if (iascii) {
166: PetscBool matl, isperm;
168: MPI_Comm_rank(PetscObjectComm((PetscObject)is), &rank);
169: MPI_Comm_size(PetscObjectComm((PetscObject)is), &size);
170: PetscViewerGetFormat(viewer, &fmt);
171: matl = (PetscBool)(fmt == PETSC_VIEWER_ASCII_MATLAB);
172: ISGetInfo(is, IS_PERMUTATION, IS_GLOBAL, PETSC_FALSE, &isperm);
173: if (isperm && !matl) PetscViewerASCIIPrintf(viewer, "Index set is permutation\n");
174: if (size == 1) {
175: if (matl) {
176: const char *name;
178: PetscObjectGetName((PetscObject)is, &name);
179: PetscViewerASCIIPrintf(viewer, "%s = [%" PetscInt_FMT " : %" PetscInt_FMT " : %" PetscInt_FMT "];\n", name, sub->first + 1, sub->step, sub->first + sub->step * (n - 1) + 1);
180: } else {
181: PetscViewerASCIIPrintf(viewer, "Number of indices in (stride) set %" PetscInt_FMT "\n", n);
182: for (i = 0; i < n; i++) PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT "\n", i, sub->first + i * sub->step);
183: }
184: PetscViewerFlush(viewer);
185: } else {
186: PetscViewerASCIIPushSynchronized(viewer);
187: if (matl) {
188: const char *name;
190: PetscObjectGetName((PetscObject)is, &name);
191: PetscViewerASCIISynchronizedPrintf(viewer, "%s_%d = [%" PetscInt_FMT " : %" PetscInt_FMT " : %" PetscInt_FMT "];\n", name, rank, sub->first + 1, sub->step, sub->first + sub->step * (n - 1) + 1);
192: } else {
193: PetscViewerASCIISynchronizedPrintf(viewer, "[%d] Number of indices in (stride) set %" PetscInt_FMT "\n", rank, n);
194: for (i = 0; i < n; i++) PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %" PetscInt_FMT " %" PetscInt_FMT "\n", rank, i, sub->first + i * sub->step);
195: }
196: PetscViewerFlush(viewer);
197: PetscViewerASCIIPopSynchronized(viewer);
198: }
199: } else if (ibinary) ISView_Binary(is, viewer);
200: return 0;
201: }
203: PetscErrorCode ISSort_Stride(IS is)
204: {
205: IS_Stride *sub = (IS_Stride *)is->data;
207: if (sub->step >= 0) return 0;
208: sub->first += (is->map->n - 1) * sub->step;
209: sub->step *= -1;
210: return 0;
211: }
213: PetscErrorCode ISSorted_Stride(IS is, PetscBool *flg)
214: {
215: IS_Stride *sub = (IS_Stride *)is->data;
217: if (sub->step >= 0) *flg = PETSC_TRUE;
218: else *flg = PETSC_FALSE;
219: return 0;
220: }
222: static PetscErrorCode ISUniqueLocal_Stride(IS is, PetscBool *flg)
223: {
224: IS_Stride *sub = (IS_Stride *)is->data;
226: if (!(is->map->n) || sub->step != 0) *flg = PETSC_TRUE;
227: else *flg = PETSC_FALSE;
228: return 0;
229: }
231: static PetscErrorCode ISPermutationLocal_Stride(IS is, PetscBool *flg)
232: {
233: IS_Stride *sub = (IS_Stride *)is->data;
235: if (!(is->map->n) || (PetscAbsInt(sub->step) == 1 && is->min == 0)) *flg = PETSC_TRUE;
236: else *flg = PETSC_FALSE;
237: return 0;
238: }
240: static PetscErrorCode ISIntervalLocal_Stride(IS is, PetscBool *flg)
241: {
242: IS_Stride *sub = (IS_Stride *)is->data;
244: if (!(is->map->n) || sub->step == 1) *flg = PETSC_TRUE;
245: else *flg = PETSC_FALSE;
246: return 0;
247: }
249: static PetscErrorCode ISOnComm_Stride(IS is, MPI_Comm comm, PetscCopyMode mode, IS *newis)
250: {
251: IS_Stride *sub = (IS_Stride *)is->data;
253: ISCreateStride(comm, is->map->n, sub->first, sub->step, newis);
254: return 0;
255: }
257: static PetscErrorCode ISSetBlockSize_Stride(IS is, PetscInt bs)
258: {
259: IS_Stride *sub = (IS_Stride *)is->data;
262: PetscLayoutSetBlockSize(is->map, bs);
263: return 0;
264: }
266: static PetscErrorCode ISContiguousLocal_Stride(IS is, PetscInt gstart, PetscInt gend, PetscInt *start, PetscBool *contig)
267: {
268: IS_Stride *sub = (IS_Stride *)is->data;
270: if (sub->step == 1 && sub->first >= gstart && sub->first + is->map->n <= gend) {
271: *start = sub->first - gstart;
272: *contig = PETSC_TRUE;
273: } else {
274: *start = -1;
275: *contig = PETSC_FALSE;
276: }
277: return 0;
278: }
280: static struct _ISOps myops = {PetscDesignatedInitializer(getindices, ISGetIndices_Stride), PetscDesignatedInitializer(restoreindices, ISRestoreIndices_Stride), PetscDesignatedInitializer(invertpermutation, ISInvertPermutation_Stride), PetscDesignatedInitializer(sort, ISSort_Stride), PetscDesignatedInitializer(sortremovedups, ISSort_Stride), PetscDesignatedInitializer(sorted, ISSorted_Stride), PetscDesignatedInitializer(duplicate, ISDuplicate_Stride), PetscDesignatedInitializer(destroy, ISDestroy_Stride), PetscDesignatedInitializer(view, ISView_Stride), PetscDesignatedInitializer(load, ISLoad_Default), PetscDesignatedInitializer(copy, ISCopy_Stride), PetscDesignatedInitializer(togeneral, ISToGeneral_Stride), PetscDesignatedInitializer(oncomm, ISOnComm_Stride), PetscDesignatedInitializer(setblocksize, ISSetBlockSize_Stride), PetscDesignatedInitializer(contiguous, ISContiguousLocal_Stride), PetscDesignatedInitializer(locate, ISLocate_Stride), PetscDesignatedInitializer(sortedlocal, ISSorted_Stride), NULL, PetscDesignatedInitializer(uniquelocal, ISUniqueLocal_Stride), NULL, PetscDesignatedInitializer(permlocal, ISPermutationLocal_Stride), NULL, PetscDesignatedInitializer(intervallocal, ISIntervalLocal_Stride), NULL};
282: /*@
283: ISStrideSetStride - Sets the stride information for a stride index set.
285: Collective on IS
287: Input Parameters:
288: + is - the index set
289: . n - the length of the locally owned portion of the index set
290: . first - the first element of the locally owned portion of the index set
291: - step - the change to the next index
293: Level: beginner
295: .seealso: `ISCreateGeneral()`, `ISCreateBlock()`, `ISAllGather()`, `ISSTRIDE`, `ISCreateStride()`, `ISStrideGetInfo()`
296: @*/
297: PetscErrorCode ISStrideSetStride(IS is, PetscInt n, PetscInt first, PetscInt step)
298: {
300: ISClearInfoCache(is, PETSC_FALSE);
301: PetscUseMethod(is, "ISStrideSetStride_C", (IS, PetscInt, PetscInt, PetscInt), (is, n, first, step));
302: return 0;
303: }
305: PetscErrorCode ISStrideSetStride_Stride(IS is, PetscInt n, PetscInt first, PetscInt step)
306: {
307: PetscInt min, max;
308: IS_Stride *sub = (IS_Stride *)is->data;
309: PetscLayout map;
311: PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)is), n, is->map->N, is->map->bs, &map);
312: PetscLayoutDestroy(&is->map);
313: is->map = map;
315: sub->first = first;
316: sub->step = step;
317: if (step > 0) {
318: min = first;
319: max = first + step * (n - 1);
320: } else {
321: max = first;
322: min = first + step * (n - 1);
323: }
325: is->min = n > 0 ? min : PETSC_MAX_INT;
326: is->max = n > 0 ? max : PETSC_MIN_INT;
327: is->data = (void *)sub;
328: return 0;
329: }
331: /*@
332: ISCreateStride - Creates a data structure for an index set
333: containing a list of evenly spaced integers.
335: Collective
337: Input Parameters:
338: + comm - the MPI communicator
339: . n - the length of the locally owned portion of the index set
340: . first - the first element of the locally owned portion of the index set
341: - step - the change to the next index
343: Output Parameter:
344: . is - the new index set
346: Notes:
347: When the communicator is not MPI_COMM_SELF, the operations on IS are NOT
348: conceptually the same as MPI_Group operations. The IS are the
349: distributed sets of indices and thus certain operations on them are collective.
351: Level: beginner
353: .seealso: `ISCreateGeneral()`, `ISCreateBlock()`, `ISAllGather()`, `ISSTRIDE`
354: @*/
355: PetscErrorCode ISCreateStride(MPI_Comm comm, PetscInt n, PetscInt first, PetscInt step, IS *is)
356: {
357: ISCreate(comm, is);
358: ISSetType(*is, ISSTRIDE);
359: ISStrideSetStride(*is, n, first, step);
360: return 0;
361: }
363: PETSC_EXTERN PetscErrorCode ISCreate_Stride(IS is)
364: {
365: IS_Stride *sub;
367: PetscNew(&sub);
368: is->data = (void *)sub;
369: PetscMemcpy(is->ops, &myops, sizeof(myops));
370: PetscObjectComposeFunction((PetscObject)is, "ISStrideSetStride_C", ISStrideSetStride_Stride);
371: PetscObjectComposeFunction((PetscObject)is, "ISShift_C", ISShift_Stride);
372: return 0;
373: }