Actual source code: stride.c
1: #define PETSCVEC_DLL
2: /*
3: Index sets of evenly space integers, defined by a
4: start, stride and length.
5: */
6: #include private/isimpl.h
7: #include petscvec.h
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: static PetscErrorCode ISCopy_Stride(IS is,IS isy)
32: {
33: IS_Stride *is_stride = (IS_Stride*)is->data,*isy_stride = (IS_Stride*)isy->data;
37: PetscMemcpy(isy_stride,is_stride,sizeof(IS_Stride));
38: return(0);
39: }
43: PetscErrorCode ISDuplicate_Stride(IS is,IS *newIS)
44: {
46: IS_Stride *sub = (IS_Stride*)is->data;
49: ISCreateStride(((PetscObject)is)->comm,sub->n,sub->first,sub->step,newIS);
50: return(0);
51: }
55: PetscErrorCode ISInvertPermutation_Stride(IS is,PetscInt nlocal,IS *perm)
56: {
57: IS_Stride *isstride = (IS_Stride*)is->data;
61: if (is->isidentity) {
62: ISCreateStride(PETSC_COMM_SELF,isstride->n,0,1,perm);
63: } else {
64: IS tmp;
65: const PetscInt *indices,n = isstride->n;
66: ISGetIndices(is,&indices);
67: ISCreateGeneral(((PetscObject)is)->comm,n,indices,&tmp);
68: ISSetPermutation(tmp);
69: ISRestoreIndices(is,&indices);
70: ISInvertPermutation(tmp,nlocal,perm);
71: ISDestroy(tmp);
72: }
73: return(0);
74: }
75:
78: /*@
79: ISStrideGetInfo - Returns the first index in a stride index set and
80: the stride width.
82: Not Collective
84: Input Parameter:
85: . is - the index set
87: Output Parameters:
88: . first - the first index
89: . step - the stride width
91: Level: intermediate
93: Notes:
94: Returns info on stride index set. This is a pseudo-public function that
95: should not be needed by most users.
97: Concepts: index sets^getting information
98: Concepts: IS^getting information
100: .seealso: ISCreateStride(), ISGetSize()
101: @*/
102: PetscErrorCode ISStrideGetInfo(IS is,PetscInt *first,PetscInt *step)
103: {
104: IS_Stride *sub;
111: sub = (IS_Stride*)is->data;
112: if (((PetscObject)is)->type != IS_STRIDE) return(0);
113: if (first) *first = sub->first;
114: if (step) *step = sub->step;
115: return(0);
116: }
120: /*@
121: ISStride - Determines if an IS is based on a stride.
123: Not Collective
125: Input Parameter:
126: . is - the index set
128: Output Parameters:
129: . flag - either PETSC_TRUE or PETSC_FALSE
131: Level: intermediate
133: Concepts: index sets^is it stride
134: Concepts: IS^is it stride
136: .seealso: ISCreateStride(), ISGetSize()
137: @*/
138: PetscErrorCode ISStride(IS is,PetscTruth *flag)
139: {
144: if (((PetscObject)is)->type != IS_STRIDE) *flag = PETSC_FALSE;
145: else *flag = PETSC_TRUE;
147: return(0);
148: }
152: PetscErrorCode ISDestroy_Stride(IS is)
153: {
157: PetscFree(is->data);
158: return(0);
159: }
161: /*
162: Returns a legitimate index memory even if
163: the stride index set is empty.
164: */
167: PetscErrorCode ISGetIndices_Stride(IS in,const PetscInt *idx[])
168: {
169: IS_Stride *sub = (IS_Stride*)in->data;
171: PetscInt i,**dx = (PetscInt**)idx;
174: PetscMalloc(sub->n*sizeof(PetscInt),idx);
175: if (sub->n) {
176: (*dx)[0] = sub->first;
177: for (i=1; i<sub->n; i++) (*dx)[i] = (*dx)[i-1] + sub->step;
178: }
179: return(0);
180: }
184: PetscErrorCode ISRestoreIndices_Stride(IS in,const PetscInt *idx[])
185: {
189: PetscFree(*(void**)idx);
190: return(0);
191: }
195: PetscErrorCode ISGetSize_Stride(IS is,PetscInt *size)
196: {
197: IS_Stride *sub = (IS_Stride *)is->data;
200: *size = sub->N;
201: return(0);
202: }
206: PetscErrorCode ISGetLocalSize_Stride(IS is,PetscInt *size)
207: {
208: IS_Stride *sub = (IS_Stride *)is->data;
211: *size = sub->n;
212: return(0);
213: }
217: PetscErrorCode ISView_Stride(IS is,PetscViewer viewer)
218: {
219: IS_Stride *sub = (IS_Stride *)is->data;
220: PetscInt i,n = sub->n;
221: PetscMPIInt rank,size;
222: PetscTruth iascii;
226: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
227: if (iascii) {
228: MPI_Comm_rank(((PetscObject)is)->comm,&rank);
229: MPI_Comm_size(((PetscObject)is)->comm,&size);
230: if (size == 1) {
231: if (is->isperm) {
232: PetscViewerASCIISynchronizedPrintf(viewer,"Index set is permutation\n");
233: }
234: PetscViewerASCIISynchronizedPrintf(viewer,"Number of indices in (stride) set %D\n",n);
235: for (i=0; i<n; i++) {
236: PetscViewerASCIISynchronizedPrintf(viewer,"%D %D\n",i,sub->first + i*sub->step);
237: }
238: } else {
239: if (is->isperm) {
240: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Index set is permutation\n",rank);
241: }
242: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of indices in (stride) set %D\n",rank,n);
243: for (i=0; i<n; i++) {
244: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D %D\n",rank,i,sub->first + i*sub->step);
245: }
246: }
247: PetscViewerFlush(viewer);
248: } else {
249: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for this object",((PetscObject)viewer)->type_name);
250: }
251: return(0);
252: }
253:
256: PetscErrorCode ISSort_Stride(IS is)
257: {
258: IS_Stride *sub = (IS_Stride*)is->data;
261: if (sub->step >= 0) return(0);
262: sub->first += (sub->n - 1)*sub->step;
263: sub->step *= -1;
264: return(0);
265: }
269: PetscErrorCode ISSorted_Stride(IS is,PetscTruth* flg)
270: {
271: IS_Stride *sub = (IS_Stride*)is->data;
274: if (sub->step >= 0) *flg = PETSC_TRUE;
275: else *flg = PETSC_FALSE;
276: return(0);
277: }
279: static struct _ISOps myops = { ISGetSize_Stride,
280: ISGetLocalSize_Stride,
281: ISGetIndices_Stride,
282: ISRestoreIndices_Stride,
283: ISInvertPermutation_Stride,
284: ISSort_Stride,
285: ISSorted_Stride,
286: ISDuplicate_Stride,
287: ISDestroy_Stride,
288: ISView_Stride,
289: ISIdentity_Stride,
290: ISCopy_Stride };
294: /*@
295: ISCreateStride - Creates a data structure for an index set
296: containing a list of evenly spaced integers.
298: Collective on MPI_Comm
300: Input Parameters:
301: + comm - the MPI communicator
302: . n - the length of the locally owned portion of the index set
303: . first - the first element of the locally owned portion of the index set
304: - step - the change to the next index
306: Output Parameter:
307: . is - the new index set
309: Notes:
310: When the communicator is not MPI_COMM_SELF, the operations on IS are NOT
311: conceptually the same as MPI_Group operations. The IS are the
312: distributed sets of indices and thus certain operations on them are collective.
314: Level: beginner
316: Concepts: IS^stride
317: Concepts: index sets^stride
318: Concepts: stride^index set
320: .seealso: ISCreateGeneral(), ISCreateBlock(), ISAllGather()
321: @*/
322: PetscErrorCode ISCreateStride(MPI_Comm comm,PetscInt n,PetscInt first,PetscInt step,IS *is)
323: {
325: PetscInt min,max;
326: IS Nindex;
327: IS_Stride *sub;
328: PetscTruth flg = PETSC_FALSE;
332: *is = PETSC_NULL;
333: if (n < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Number of indices < 0");
334: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
335: ISInitializePackage(PETSC_NULL);
336: #endif
338: PetscHeaderCreate(Nindex,_p_IS,struct _ISOps,IS_COOKIE,IS_STRIDE,"IS",comm,ISDestroy,ISView);
339: PetscNewLog(Nindex,IS_Stride,&sub);
340: sub->n = n;
341: MPI_Allreduce(&n,&sub->N,1,MPIU_INT,MPI_SUM,comm);
342: sub->first = first;
343: sub->step = step;
344: if (step > 0) {min = first; max = first + step*(n-1);}
345: else {max = first; min = first + step*(n-1);}
347: Nindex->min = min;
348: Nindex->max = max;
349: Nindex->data = (void*)sub;
350: PetscMemcpy(Nindex->ops,&myops,sizeof(myops));
352: if ((!first && step == 1) || (first == max && step == -1 && !min)) {
353: Nindex->isperm = PETSC_TRUE;
354: } else {
355: Nindex->isperm = PETSC_FALSE;
356: }
357: PetscOptionsGetTruth(PETSC_NULL,"-is_view",&flg,PETSC_NULL);
358: if (flg) {
359: PetscViewer viewer;
360: PetscViewerASCIIGetStdout(((PetscObject)Nindex)->comm,&viewer);
361: ISView(Nindex,viewer);
362: }
363: *is = Nindex;
364: return(0);
365: }