Actual source code: stride.c
1: /*
2: Index sets of evenly space integers, defined by a
3: start, stride and length.
4: */
5: #include src/vec/is/isimpl.h
7: EXTERN PetscErrorCode VecInitializePackage(char *);
9: typedef struct {
10: PetscInt N,n,first,step;
11: } IS_Stride;
15: PetscErrorCode ISIdentity_Stride(IS is,PetscTruth *ident)
16: {
17: IS_Stride *is_stride = (IS_Stride*)is->data;
20: is->isidentity = PETSC_FALSE;
21: *ident = PETSC_FALSE;
22: if (is_stride->first != 0) return(0);
23: if (is_stride->step != 1) return(0);
24: *ident = PETSC_TRUE;
25: is->isidentity = PETSC_TRUE;
26: return(0);
27: }
31: PetscErrorCode ISDuplicate_Stride(IS is,IS *newIS)
32: {
34: IS_Stride *sub = (IS_Stride*)is->data;
37: ISCreateStride(is->comm,sub->n,sub->first,sub->step,newIS);
38: return(0);
39: }
43: PetscErrorCode ISInvertPermutation_Stride(IS is,PetscInt nlocal,IS *perm)
44: {
45: IS_Stride *isstride = (IS_Stride*)is->data;
49: if (is->isidentity) {
50: ISCreateStride(PETSC_COMM_SELF,isstride->n,0,1,perm);
51: } else {
52: IS tmp;
53: PetscInt *indices,n = isstride->n;
54: ISGetIndices(is,&indices);
55: ISCreateGeneral(is->comm,n,indices,&tmp);
56: ISRestoreIndices(is,&indices);
57: ISInvertPermutation(tmp,nlocal,perm);
58: ISDestroy(tmp);
59: }
60: return(0);
61: }
62:
65: /*@
66: ISStrideGetInfo - Returns the first index in a stride index set and
67: the stride width.
69: Not Collective
71: Input Parameter:
72: . is - the index set
74: Output Parameters:
75: . first - the first index
76: . step - the stride width
78: Level: intermediate
80: Notes:
81: Returns info on stride index set. This is a pseudo-public function that
82: should not be needed by most users.
84: Concepts: index sets^getting information
85: Concepts: IS^getting information
87: .seealso: ISCreateStride(), ISGetSize()
88: @*/
89: PetscErrorCode ISStrideGetInfo(IS is,PetscInt *first,PetscInt *step)
90: {
91: IS_Stride *sub;
98: sub = (IS_Stride*)is->data;
99: if (is->type != IS_STRIDE) return(0);
100: if (first) *first = sub->first;
101: if (step) *step = sub->step;
102: return(0);
103: }
107: /*@C
108: ISStride - Determines if an IS is based on a stride.
110: Not Collective
112: Input Parameter:
113: . is - the index set
115: Output Parameters:
116: . flag - either PETSC_TRUE or PETSC_FALSE
118: Level: intermediate
120: Concepts: index sets^is it stride
121: Concepts: IS^is it stride
123: .seealso: ISCreateStride(), ISGetSize()
124: @*/
125: PetscErrorCode ISStride(IS is,PetscTruth *flag)
126: {
131: if (is->type != IS_STRIDE) *flag = PETSC_FALSE;
132: else *flag = PETSC_TRUE;
134: return(0);
135: }
139: PetscErrorCode ISDestroy_Stride(IS is)
140: {
144: PetscFree(is->data);
145: PetscLogObjectDestroy(is);
146: PetscHeaderDestroy(is); return(0);
147: }
149: /*
150: Returns a legitimate index memory even if
151: the stride index set is empty.
152: */
155: PetscErrorCode ISGetIndices_Stride(IS in,PetscInt **idx)
156: {
157: IS_Stride *sub = (IS_Stride*)in->data;
159: PetscInt i;
162: PetscMalloc(sub->n*sizeof(PetscInt),idx);
163: if (sub->n) {
164: (*idx)[0] = sub->first;
165: for (i=1; i<sub->n; i++) (*idx)[i] = (*idx)[i-1] + sub->step;
166: }
167: return(0);
168: }
172: PetscErrorCode ISRestoreIndices_Stride(IS in,PetscInt **idx)
173: {
177: PetscFree(*idx);
178: return(0);
179: }
183: PetscErrorCode ISGetSize_Stride(IS is,PetscInt *size)
184: {
185: IS_Stride *sub = (IS_Stride *)is->data;
188: *size = sub->N;
189: return(0);
190: }
194: PetscErrorCode ISGetLocalSize_Stride(IS is,PetscInt *size)
195: {
196: IS_Stride *sub = (IS_Stride *)is->data;
199: *size = sub->n;
200: return(0);
201: }
205: PetscErrorCode ISView_Stride(IS is,PetscViewer viewer)
206: {
207: IS_Stride *sub = (IS_Stride *)is->data;
208: PetscInt i,n = sub->n;
209: PetscMPIInt rank,size;
210: PetscTruth iascii;
214: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
215: if (iascii) {
216: MPI_Comm_rank(is->comm,&rank);
217: MPI_Comm_size(is->comm,&size);
218: if (size == 1) {
219: if (is->isperm) {
220: PetscViewerASCIISynchronizedPrintf(viewer,"Index set is permutation\n");
221: }
222: PetscViewerASCIISynchronizedPrintf(viewer,"Number of indices in (stride) set %D\n",n);
223: for (i=0; i<n; i++) {
224: PetscViewerASCIISynchronizedPrintf(viewer,"%D %D\n",i,sub->first + i*sub->step);
225: }
226: } else {
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: }
235: PetscViewerFlush(viewer);
236: } else {
237: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for this object",((PetscObject)viewer)->type_name);
238: }
239: return(0);
240: }
241:
244: PetscErrorCode ISSort_Stride(IS is)
245: {
246: IS_Stride *sub = (IS_Stride*)is->data;
249: if (sub->step >= 0) return(0);
250: sub->first += (sub->n - 1)*sub->step;
251: sub->step *= -1;
252: return(0);
253: }
257: PetscErrorCode ISSorted_Stride(IS is,PetscTruth* flg)
258: {
259: IS_Stride *sub = (IS_Stride*)is->data;
262: if (sub->step >= 0) *flg = PETSC_TRUE;
263: else *flg = PETSC_FALSE;
264: return(0);
265: }
267: static struct _ISOps myops = { ISGetSize_Stride,
268: ISGetLocalSize_Stride,
269: ISGetIndices_Stride,
270: ISRestoreIndices_Stride,
271: ISInvertPermutation_Stride,
272: ISSort_Stride,
273: ISSorted_Stride,
274: ISDuplicate_Stride,
275: ISDestroy_Stride,
276: ISView_Stride,
277: ISIdentity_Stride };
281: /*@C
282: ISCreateStride - Creates a data structure for an index set
283: containing a list of evenly spaced integers.
285: Collective on MPI_Comm
287: Input Parameters:
288: + comm - the MPI communicator
289: . n - the length of the index set
290: . first - the first element of the index set
291: - step - the change to the next index
293: Output Parameter:
294: . is - the new index set
296: Notes:
297: When the communicator is not MPI_COMM_SELF, the operations on IS are NOT
298: conceptually the same as MPI_Group operations. The IS are the
299: distributed sets of indices and thus certain operations on them are collective.
301: Level: beginner
303: Concepts: IS^stride
304: Concepts: index sets^stride
305: Concepts: stride^index set
307: .seealso: ISCreateGeneral(), ISCreateBlock(), ISAllGather()
308: @*/
309: PetscErrorCode ISCreateStride(MPI_Comm comm,PetscInt n,PetscInt first,PetscInt step,IS *is)
310: {
312: PetscInt min,max;
313: IS Nindex;
314: IS_Stride *sub;
315: PetscTruth flg;
319: *is = PETSC_NULL;
320: if (n < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Number of indices < 0");
321: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
322: VecInitializePackage(PETSC_NULL);
323: #endif
325: PetscHeaderCreate(Nindex,_p_IS,struct _ISOps,IS_COOKIE,IS_STRIDE,"IS",comm,ISDestroy,ISView);
326: PetscLogObjectCreate(Nindex);
327: PetscLogObjectMemory(Nindex,sizeof(IS_Stride) + sizeof(struct _p_IS));
328: PetscNew(IS_Stride,&sub);
329: sub->n = n;
330: MPI_Allreduce(&n,&sub->N,1,MPIU_INT,MPI_SUM,comm);
331: sub->first = first;
332: sub->step = step;
333: if (step > 0) {min = first; max = first + step*(n-1);}
334: else {max = first; min = first + step*(n-1);}
336: Nindex->min = min;
337: Nindex->max = max;
338: Nindex->data = (void*)sub;
339: PetscMemcpy(Nindex->ops,&myops,sizeof(myops));
341: if ((!first && step == 1) || (first == max && step == -1 && !min)) {
342: Nindex->isperm = PETSC_TRUE;
343: } else {
344: Nindex->isperm = PETSC_FALSE;
345: }
346: PetscOptionsHasName(PETSC_NULL,"-is_view",&flg);
347: if (flg) {
348: ISView(Nindex,PETSC_VIEWER_STDOUT_(Nindex->comm));
349: }
350: *is = Nindex;
351: return(0);
352: }