Actual source code: stride.c

  1: /*$Id: stride.c,v 1.103 2001/03/23 23:21:09 balay 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"             /*I   "petscis.h"   I*/

  8: typedef struct {
  9:   int N,n,first,step;
 10: } IS_Stride;

 12: int ISIdentity_Stride(IS is,PetscTruth *ident)
 13: {
 14:   IS_Stride *is_stride = (IS_Stride*)is->data;

 17:   is->isidentity = PETSC_FALSE;
 18:   *ident         = PETSC_FALSE;
 19:   if (is_stride->first != 0) return(0);
 20:   if (is_stride->step  != 1) return(0);
 21:   *ident          = PETSC_TRUE;
 22:   is->isidentity  = PETSC_TRUE;
 23:   return(0);
 24: }

 26: int ISDuplicate_Stride(IS is,IS *newIS)
 27: {
 28:   int       ierr;
 29:   IS_Stride *sub = (IS_Stride*)is->data;

 32:   ISCreateStride(is->comm,sub->n,sub->first,sub->step,newIS);
 33:   return(0);
 34: }

 36: int ISInvertPermutation_Stride(IS is,int nlocal,IS *perm)
 37: {
 38:   IS_Stride *isstride = (IS_Stride*)is->data;
 39:   int       ierr;

 42:   if (is->isidentity) {
 43:     ISCreateStride(PETSC_COMM_SELF,isstride->n,0,1,perm);
 44:     ISSetPermutation(*perm);
 45:   } else {
 46:     IS  tmp;
 47:     int *indices,n = isstride->n;
 48:     ISGetIndices(is,&indices);
 49:     ISCreateGeneral(is->comm,n,indices,&tmp);
 50:     ISRestoreIndices(is,&indices);
 51:     ISInvertPermutation(tmp,nlocal,perm);
 52:     ISDestroy(tmp);
 53:   }
 54:   return(0);
 55: }
 56: 
 57: /*@
 58:    ISStrideGetInfo - Returns the first index in a stride index set and 
 59:    the stride width.

 61:    Not Collective

 63:    Input Parameter:
 64: .  is - the index set

 66:    Output Parameters:
 67: .  first - the first index
 68: .  step - the stride width

 70:    Level: intermediate

 72:    Notes:
 73:    Returns info on stride index set. This is a pseudo-public function that
 74:    should not be needed by most users.

 76:    Concepts: index sets^getting information
 77:    Concepts: IS^getting information

 79: .seealso: ISCreateStride(), ISGetSize()
 80: @*/
 81: int ISStrideGetInfo(IS is,int *first,int *step)
 82: {
 83:   IS_Stride *sub;


 90:   sub = (IS_Stride*)is->data;
 91:   if (is->type != IS_STRIDE) return(0);
 92:   if (first) *first = sub->first;
 93:   if (step)  *step  = sub->step;
 94:   return(0);
 95: }

 97: /*@C
 98:    ISStride - Determines if an IS is based on a stride.

100:    Not Collective

102:    Input Parameter:
103: .  is - the index set

105:    Output Parameters:
106: .  flag - either PETSC_TRUE or PETSC_FALSE

108:    Level: intermediate

110:    Concepts: index sets^is it stride
111:    Concepts: IS^is it stride

113: .seealso: ISCreateStride(), ISGetSize()
114: @*/
115: int ISStride(IS is,PetscTruth *flag)
116: {

121:   if (is->type != IS_STRIDE) *flag = PETSC_FALSE;
122:   else                       *flag = PETSC_TRUE;

124:   return(0);
125: }

127: int ISDestroy_Stride(IS is)
128: {

132:   PetscFree(is->data);
133:   PetscLogObjectDestroy(is);
134:   PetscHeaderDestroy(is); return(0);
135: }

137: /*
138:      Returns a legitimate index memory even if 
139:    the stride index set is empty.
140: */
141: int ISGetIndices_Stride(IS in,int **idx)
142: {
143:   IS_Stride *sub = (IS_Stride*)in->data;
144:   int       i,ierr;

147:   ierr      = PetscMalloc((sub->n+1)*sizeof(int),idx);
148:   (*idx)[0] = sub->first;
149:   for (i=1; i<sub->n; i++) (*idx)[i] = (*idx)[i-1] + sub->step;
150:   return(0);
151: }

153: int ISRestoreIndices_Stride(IS in,int **idx)
154: {

158:   PetscFree(*idx);
159:   return(0);
160: }

162: int ISGetSize_Stride(IS is,int *size)
163: {
164:   IS_Stride *sub = (IS_Stride *)is->data;

167:   *size = sub->N;
168:   return(0);
169: }

171: int ISGetLocalSize_Stride(IS is,int *size)
172: {
173:   IS_Stride *sub = (IS_Stride *)is->data;

176:   *size = sub->n;
177:   return(0);
178: }

180: int ISView_Stride(IS is,PetscViewer viewer)
181: {
182:   IS_Stride   *sub = (IS_Stride *)is->data;
183:   int         i,n = sub->n,ierr,rank,size;
184:   PetscTruth  isascii;

187:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
188:   if (isascii) {
189:     MPI_Comm_rank(is->comm,&rank);
190:     MPI_Comm_size(is->comm,&size);
191:     if (size == 1) {
192:       if (is->isperm) {
193:         PetscViewerASCIISynchronizedPrintf(viewer,"Index set is permutationn");
194:       }
195:       PetscViewerASCIISynchronizedPrintf(viewer,"Number of indices in (stride) set %dn",n);
196:       for (i=0; i<n; i++) {
197:         PetscViewerASCIISynchronizedPrintf(viewer,"%d %dn",i,sub->first + i*sub->step);
198:       }
199:     } else {
200:       if (is->isperm) {
201:         PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Index set is permutationn",rank);
202:       }
203:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of indices in (stride) set %dn",rank,n);
204:       for (i=0; i<n; i++) {
205:         PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %d %dn",rank,i,sub->first + i*sub->step);
206:       }
207:     }
208:     PetscViewerFlush(viewer);
209:   } else {
210:     SETERRQ1(1,"Viewer type %s not supported for this object",((PetscObject)viewer)->type_name);
211:   }
212:   return(0);
213: }
214: 
215: int ISSort_Stride(IS is)
216: {
217:   IS_Stride *sub = (IS_Stride*)is->data;

220:   if (sub->step >= 0) return(0);
221:   sub->first += (sub->n - 1)*sub->step;
222:   sub->step *= -1;
223:   return(0);
224: }

226: int ISSorted_Stride(IS is,PetscTruth* flg)
227: {
228:   IS_Stride *sub = (IS_Stride*)is->data;

231:   if (sub->step >= 0) *flg = PETSC_TRUE;
232:   else *flg = PETSC_FALSE;
233:   return(0);
234: }

236: static struct _ISOps myops = { ISGetSize_Stride,
237:                                ISGetLocalSize_Stride,
238:                                ISGetIndices_Stride,
239:                                ISRestoreIndices_Stride,
240:                                ISInvertPermutation_Stride,
241:                                ISSort_Stride,
242:                                ISSorted_Stride,
243:                                ISDuplicate_Stride,
244:                                ISDestroy_Stride,
245:                                ISView_Stride,
246:                                ISIdentity_Stride };

248: /*@C
249:    ISCreateStride - Creates a data structure for an index set 
250:    containing a list of evenly spaced integers.

252:    Collective on MPI_Comm

254:    Input Parameters:
255: +  comm - the MPI communicator
256: .  n - the length of the index set
257: .  first - the first element of the index set
258: -  step - the change to the next index

260:    Output Parameter:
261: .  is - the new index set

263:    Notes: 
264:    When the communicator is not MPI_COMM_SELF, the operations on IS are NOT
265:    conceptually the same as MPI_Group operations. The IS are the 
266:    distributed sets of indices and thus certain operations on them are collective. 

268:    Level: beginner

270:   Concepts: IS^stride
271:   Concepts: index sets^stride
272:   Concepts: stride^index set

274: .seealso: ISCreateGeneral(), ISCreateBlock(), ISAllGather()
275: @*/
276: int ISCreateStride(MPI_Comm comm,int n,int first,int step,IS *is)
277: {
278:   int        min,max,ierr;
279:   IS         Nindex;
280:   IS_Stride  *sub;
281:   PetscTruth flg;

284:   *is = 0;
285:   if (n < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Number of indices < 0");

287:   PetscHeaderCreate(Nindex,_p_IS,struct _ISOps,IS_COOKIE,IS_STRIDE,"IS",comm,ISDestroy,ISView);
288:   PetscLogObjectCreate(Nindex);
289:   PetscLogObjectMemory(Nindex,sizeof(IS_Stride) + sizeof(struct _p_IS));
290:   ierr           = PetscNew(IS_Stride,&sub);
291:   sub->n         = n;
292:   MPI_Allreduce(&n,&sub->N,1,MPI_INT,MPI_SUM,comm);
293:   sub->first     = first;
294:   sub->step      = step;
295:   if (step > 0) {min = first; max = first + step*(n-1);}
296:   else          {max = first; min = first + step*(n-1);}

298:   Nindex->min     = min;
299:   Nindex->max     = max;
300:   Nindex->data    = (void*)sub;
301:   PetscMemcpy(Nindex->ops,&myops,sizeof(myops));
302:   Nindex->isperm  = PETSC_FALSE;
303:   PetscOptionsHasName(PETSC_NULL,"-is_view",&flg);
304:   if (flg) {
305:     ISView(Nindex,PETSC_VIEWER_STDOUT_(Nindex->comm));
306:   }
307:   *is = Nindex;
308:   return(0);
309: }