Actual source code: stride.c
1: /*$Id: stride.c,v 1.104 2001/06/21 21:15:53 bsmith Exp $*/
2: /*
3: Index sets of evenly space integers, defined by a
4: start, stride and length.
5: */
6: #include src/vec/is/isimpl.h
8: EXTERN int VecInitializePackage(char *);
10: typedef struct {
11: int N,n,first,step;
12: } IS_Stride;
14: int ISIdentity_Stride(IS is,PetscTruth *ident)
15: {
16: IS_Stride *is_stride = (IS_Stride*)is->data;
19: is->isidentity = PETSC_FALSE;
20: *ident = PETSC_FALSE;
21: if (is_stride->first != 0) return(0);
22: if (is_stride->step != 1) return(0);
23: *ident = PETSC_TRUE;
24: is->isidentity = PETSC_TRUE;
25: return(0);
26: }
28: int ISDuplicate_Stride(IS is,IS *newIS)
29: {
30: int ierr;
31: IS_Stride *sub = (IS_Stride*)is->data;
34: ISCreateStride(is->comm,sub->n,sub->first,sub->step,newIS);
35: return(0);
36: }
38: int ISInvertPermutation_Stride(IS is,int nlocal,IS *perm)
39: {
40: IS_Stride *isstride = (IS_Stride*)is->data;
41: int ierr;
44: if (is->isidentity) {
45: ISCreateStride(PETSC_COMM_SELF,isstride->n,0,1,perm);
46: } else {
47: IS tmp;
48: int *indices,n = isstride->n;
49: ISGetIndices(is,&indices);
50: ISCreateGeneral(is->comm,n,indices,&tmp);
51: ISRestoreIndices(is,&indices);
52: ISInvertPermutation(tmp,nlocal,perm);
53: ISDestroy(tmp);
54: }
55: return(0);
56: }
57:
58: /*@
59: ISStrideGetInfo - Returns the first index in a stride index set and
60: 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: Concepts: index sets^getting information
78: Concepts: IS^getting information
80: .seealso: ISCreateStride(), ISGetSize()
81: @*/
82: int ISStrideGetInfo(IS is,int *first,int *step)
83: {
84: IS_Stride *sub;
91: sub = (IS_Stride*)is->data;
92: if (is->type != IS_STRIDE) return(0);
93: if (first) *first = sub->first;
94: if (step) *step = sub->step;
95: return(0);
96: }
98: /*@C
99: ISStride - Determines if an IS is based on a stride.
101: Not Collective
103: Input Parameter:
104: . is - the index set
106: Output Parameters:
107: . flag - either PETSC_TRUE or PETSC_FALSE
109: Level: intermediate
111: Concepts: index sets^is it stride
112: Concepts: IS^is it stride
114: .seealso: ISCreateStride(), ISGetSize()
115: @*/
116: int ISStride(IS is,PetscTruth *flag)
117: {
122: if (is->type != IS_STRIDE) *flag = PETSC_FALSE;
123: else *flag = PETSC_TRUE;
125: return(0);
126: }
128: int ISDestroy_Stride(IS is)
129: {
133: PetscFree(is->data);
134: PetscLogObjectDestroy(is);
135: PetscHeaderDestroy(is); return(0);
136: }
138: /*
139: Returns a legitimate index memory even if
140: the stride index set is empty.
141: */
142: int ISGetIndices_Stride(IS in,int **idx)
143: {
144: IS_Stride *sub = (IS_Stride*)in->data;
145: int i,ierr;
148: ierr = PetscMalloc((sub->n+1)*sizeof(int),idx);
149: (*idx)[0] = sub->first;
150: for (i=1; i<sub->n; i++) (*idx)[i] = (*idx)[i-1] + sub->step;
151: return(0);
152: }
154: int ISRestoreIndices_Stride(IS in,int **idx)
155: {
159: PetscFree(*idx);
160: return(0);
161: }
163: int ISGetSize_Stride(IS is,int *size)
164: {
165: IS_Stride *sub = (IS_Stride *)is->data;
168: *size = sub->N;
169: return(0);
170: }
172: int ISGetLocalSize_Stride(IS is,int *size)
173: {
174: IS_Stride *sub = (IS_Stride *)is->data;
177: *size = sub->n;
178: return(0);
179: }
181: int ISView_Stride(IS is,PetscViewer viewer)
182: {
183: IS_Stride *sub = (IS_Stride *)is->data;
184: int i,n = sub->n,ierr,rank,size;
185: PetscTruth isascii;
188: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
189: if (isascii) {
190: MPI_Comm_rank(is->comm,&rank);
191: MPI_Comm_size(is->comm,&size);
192: if (size == 1) {
193: if (is->isperm) {
194: PetscViewerASCIISynchronizedPrintf(viewer,"Index set is permutationn");
195: }
196: PetscViewerASCIISynchronizedPrintf(viewer,"Number of indices in (stride) set %dn",n);
197: for (i=0; i<n; i++) {
198: PetscViewerASCIISynchronizedPrintf(viewer,"%d %dn",i,sub->first + i*sub->step);
199: }
200: } else {
201: if (is->isperm) {
202: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Index set is permutationn",rank);
203: }
204: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of indices in (stride) set %dn",rank,n);
205: for (i=0; i<n; i++) {
206: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %d %dn",rank,i,sub->first + i*sub->step);
207: }
208: }
209: PetscViewerFlush(viewer);
210: } else {
211: SETERRQ1(1,"Viewer type %s not supported for this object",((PetscObject)viewer)->type_name);
212: }
213: return(0);
214: }
215:
216: int ISSort_Stride(IS is)
217: {
218: IS_Stride *sub = (IS_Stride*)is->data;
221: if (sub->step >= 0) return(0);
222: sub->first += (sub->n - 1)*sub->step;
223: sub->step *= -1;
224: return(0);
225: }
227: int ISSorted_Stride(IS is,PetscTruth* flg)
228: {
229: IS_Stride *sub = (IS_Stride*)is->data;
232: if (sub->step >= 0) *flg = PETSC_TRUE;
233: else *flg = PETSC_FALSE;
234: return(0);
235: }
237: static struct _ISOps myops = { ISGetSize_Stride,
238: ISGetLocalSize_Stride,
239: ISGetIndices_Stride,
240: ISRestoreIndices_Stride,
241: ISInvertPermutation_Stride,
242: ISSort_Stride,
243: ISSorted_Stride,
244: ISDuplicate_Stride,
245: ISDestroy_Stride,
246: ISView_Stride,
247: ISIdentity_Stride };
249: /*@C
250: ISCreateStride - Creates a data structure for an index set
251: containing a list of evenly spaced integers.
253: Collective on MPI_Comm
255: Input Parameters:
256: + comm - the MPI communicator
257: . n - the length of the index set
258: . first - the first element of the index set
259: - step - the change to the next index
261: Output Parameter:
262: . is - the new index set
264: Notes:
265: When the communicator is not MPI_COMM_SELF, the operations on IS are NOT
266: conceptually the same as MPI_Group operations. The IS are the
267: distributed sets of indices and thus certain operations on them are collective.
269: Level: beginner
271: Concepts: IS^stride
272: Concepts: index sets^stride
273: Concepts: stride^index set
275: .seealso: ISCreateGeneral(), ISCreateBlock(), ISAllGather()
276: @*/
277: int ISCreateStride(MPI_Comm comm,int n,int first,int step,IS *is)
278: {
279: int min,max,ierr;
280: IS Nindex;
281: IS_Stride *sub;
282: PetscTruth flg;
286: *is = PETSC_NULL;
287: if (n < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Number of indices < 0");
288: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
289: VecInitializePackage(PETSC_NULL);
290: #endif
292: PetscHeaderCreate(Nindex,_p_IS,struct _ISOps,IS_COOKIE,IS_STRIDE,"IS",comm,ISDestroy,ISView);
293: PetscLogObjectCreate(Nindex);
294: PetscLogObjectMemory(Nindex,sizeof(IS_Stride) + sizeof(struct _p_IS));
295: ierr = PetscNew(IS_Stride,&sub);
296: sub->n = n;
297: MPI_Allreduce(&n,&sub->N,1,MPI_INT,MPI_SUM,comm);
298: sub->first = first;
299: sub->step = step;
300: if (step > 0) {min = first; max = first + step*(n-1);}
301: else {max = first; min = first + step*(n-1);}
303: Nindex->min = min;
304: Nindex->max = max;
305: Nindex->data = (void*)sub;
306: PetscMemcpy(Nindex->ops,&myops,sizeof(myops));
308: if ((first == 0 && step == 1) || (first == max && step == -1 && min == 0)) {
309: Nindex->isperm = PETSC_TRUE;
310: } else {
311: Nindex->isperm = PETSC_FALSE;
312: }
313: PetscOptionsHasName(PETSC_NULL,"-is_view",&flg);
314: if (flg) {
315: ISView(Nindex,PETSC_VIEWER_STDOUT_(Nindex->comm));
316: }
317: *is = Nindex;
318: return(0);
319: }