Actual source code: general.c

  1: #define PETSCVEC_DLL
  2: /*
  3:      Provides the functions for index sets (IS) defined by a list of integers.
  4: */
 5:  #include ../src/vec/is/impls/general/general.h
 6:  #include petscvec.h

 10: PetscErrorCode ISDuplicate_General(IS is,IS *newIS)
 11: {
 13:   IS_General     *sub = (IS_General *)is->data;

 16:   ISCreateGeneral(((PetscObject)is)->comm,sub->n,sub->idx,newIS);
 17:   return(0);
 18: }

 22: PetscErrorCode ISDestroy_General(IS is)
 23: {
 24:   IS_General     *is_general = (IS_General*)is->data;

 28:   if (is_general->allocated) {
 29:     PetscFree(is_general->idx);
 30:   }
 31:   PetscFree(is_general);
 32:   PetscHeaderDestroy(is);
 33:   return(0);
 34: }

 38: PetscErrorCode ISIdentity_General(IS is,PetscTruth *ident)
 39: {
 40:   IS_General *is_general = (IS_General*)is->data;
 41:   PetscInt   i,n = is_general->n,*idx = is_general->idx;

 44:   is->isidentity = PETSC_TRUE;
 45:   *ident         = PETSC_TRUE;
 46:   for (i=0; i<n; i++) {
 47:     if (idx[i] != i) {
 48:       is->isidentity = PETSC_FALSE;
 49:       *ident         = PETSC_FALSE;
 50:       break;
 51:     }
 52:   }
 53:   return(0);
 54: }

 58: PetscErrorCode ISGetIndices_General(IS in,const PetscInt *idx[])
 59: {
 60:   IS_General *sub = (IS_General*)in->data;

 63:   *idx = sub->idx;
 64:   return(0);
 65: }

 69: PetscErrorCode ISRestoreIndices_General(IS in,const PetscInt *idx[])
 70: {
 71:   IS_General *sub = (IS_General*)in->data;

 74:   if (*idx != sub->idx) {
 75:     SETERRQ(PETSC_ERR_ARG_WRONG,"Must restore with value from ISGetIndices()");
 76:   }
 77:   return(0);
 78: }

 82: PetscErrorCode ISGetSize_General(IS is,PetscInt *size)
 83: {
 84:   IS_General *sub = (IS_General *)is->data;

 87:   *size = sub->N;
 88:   return(0);
 89: }

 93: PetscErrorCode ISGetLocalSize_General(IS is,PetscInt *size)
 94: {
 95:   IS_General *sub = (IS_General *)is->data;

 98:   *size = sub->n;
 99:   return(0);
100: }

104: PetscErrorCode ISInvertPermutation_General(IS is,PetscInt nlocal,IS *isout)
105: {
106:   IS_General     *sub = (IS_General *)is->data;
107:   PetscInt       i,*ii,n = sub->n,nstart;
108:   const PetscInt *idx = sub->idx;
109:   PetscMPIInt    size;
110:   IS             istmp,nistmp;

114:   MPI_Comm_size(((PetscObject)is)->comm,&size);
115:   if (size == 1) {
116:     PetscMalloc(n*sizeof(PetscInt),&ii);
117:     for (i=0; i<n; i++) {
118:       ii[idx[i]] = i;
119:     }
120:     ISCreateGeneral(PETSC_COMM_SELF,n,ii,isout);
121:     ISSetPermutation(*isout);
122:     PetscFree(ii);
123:   } else {
124:     /* crude, nonscalable get entire IS on each processor */
125:     if (nlocal == PETSC_DECIDE) SETERRQ(PETSC_ERR_SUP,"Do not yet support nlocal of PETSC_DECIDE");
126:     ISAllGather(is,&istmp);
127:     ISSetPermutation(istmp);
128:     ISInvertPermutation(istmp,PETSC_DECIDE,&nistmp);
129:     ISDestroy(istmp);
130:     /* get the part we need */
131:     MPI_Scan(&nlocal,&nstart,1,MPIU_INT,MPI_SUM,((PetscObject)is)->comm);
132: #if defined(PETSC_USE_DEBUG)
133:     {
134:       PetscMPIInt rank;
135:       MPI_Comm_rank(((PetscObject)is)->comm,&rank);
136:       if (rank == size-1) {
137:         if (nstart != sub->N) SETERRQ2(PETSC_ERR_ARG_INCOMP,"Sum of nlocal lengths %d != total IS length %d",nstart,sub->N);
138:       }
139:     }
140: #endif
141:     nstart -= nlocal;
142:     ISGetIndices(nistmp,&idx);
143:     ISCreateGeneral(((PetscObject)is)->comm,nlocal,idx+nstart,isout);
144:     ISRestoreIndices(nistmp,&idx);
145:     ISDestroy(nistmp);
146:   }
147:   return(0);
148: }

152: PetscErrorCode ISView_General(IS is,PetscViewer viewer)
153: {
154:   IS_General     *sub = (IS_General *)is->data;
156:   PetscInt       i,n = sub->n,*idx = sub->idx;
157:   PetscTruth     iascii;

160:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
161:   if (iascii) {
162:     MPI_Comm    comm;
163:     PetscMPIInt rank,size;

165:     PetscObjectGetComm((PetscObject)viewer,&comm);
166:     MPI_Comm_rank(comm,&rank);
167:     MPI_Comm_size(comm,&size);

169:     if (size > 1) {
170:       if (is->isperm) {
171:         PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Index set is permutation\n",rank);
172:       }
173:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of indices in set %D\n",rank,n);
174:       for (i=0; i<n; i++) {
175:         PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D %D\n",rank,i,idx[i]);
176:       }
177:     } else {
178:       if (is->isperm) {
179:         PetscViewerASCIISynchronizedPrintf(viewer,"Index set is permutation\n");
180:       }
181:       PetscViewerASCIISynchronizedPrintf(viewer,"Number of indices in set %D\n",n);
182:       for (i=0; i<n; i++) {
183:         PetscViewerASCIISynchronizedPrintf(viewer,"%D %D\n",i,idx[i]);
184:       }
185:     }
186:     PetscViewerFlush(viewer);
187:   } else {
188:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for this object",((PetscObject)viewer)->type_name);
189:   }
190:   return(0);
191: }

195: PetscErrorCode ISSort_General(IS is)
196: {
197:   IS_General     *sub = (IS_General *)is->data;

201:   if (sub->sorted) return(0);
202:   PetscSortInt(sub->n,sub->idx);
203:   sub->sorted = PETSC_TRUE;
204:   return(0);
205: }

209: PetscErrorCode ISSorted_General(IS is,PetscTruth *flg)
210: {
211:   IS_General *sub = (IS_General *)is->data;

214:   *flg = sub->sorted;
215:   return(0);
216: }

218: static struct _ISOps myops = { ISGetSize_General,
219:                                ISGetLocalSize_General,
220:                                ISGetIndices_General,
221:                                ISRestoreIndices_General,
222:                                ISInvertPermutation_General,
223:                                ISSort_General,
224:                                ISSorted_General,
225:                                ISDuplicate_General,
226:                                ISDestroy_General,
227:                                ISView_General,
228:                                ISIdentity_General };

232: PetscErrorCode ISCreateGeneral_Private(MPI_Comm comm,IS *is)
233: {
235:   IS             Nindex = *is;
236:   IS_General     *sub = (IS_General*)Nindex->data;
237:   PetscInt       n = sub->n,i,min,max;
238:   const PetscInt *idx = sub->idx;
239:   PetscTruth     sorted = PETSC_TRUE;
240:   PetscTruth     flg;

244:   if (n < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"length < 0");
246:   *is = PETSC_NULL;
247: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
248:   ISInitializePackage(PETSC_NULL);
249: #endif

251:   MPI_Allreduce(&n,&sub->N,1,MPIU_INT,MPI_SUM,comm);
252:   for (i=1; i<n; i++) {
253:     if (idx[i] < idx[i-1]) {sorted = PETSC_FALSE; break;}
254:   }
255:   if (n) {min = max = idx[0];} else {min = max = 0;}
256:   for (i=1; i<n; i++) {
257:     if (idx[i] < min) min = idx[i];
258:     if (idx[i] > max) max = idx[i];
259:   }
260:   sub->sorted     = sorted;
261:   Nindex->min     = min;
262:   Nindex->max     = max;
263:   PetscMemcpy(Nindex->ops,&myops,sizeof(myops));
264:   Nindex->isperm     = PETSC_FALSE;
265:   Nindex->isidentity = PETSC_FALSE;
266:   PetscOptionsHasName(PETSC_NULL,"-is_view",&flg);
267:   if (flg) {
268:     PetscViewer viewer;
269:     PetscViewerASCIIGetStdout(((PetscObject)Nindex)->comm,&viewer);
270:     ISView(Nindex,viewer);
271:   }
272:   *is = Nindex;
273:   return(0);
274: }


279: /*@
280:    ISCreateGeneral - Creates a data structure for an index set 
281:    containing a list of integers.

283:    Collective on MPI_Comm

285:    Input Parameters:
286: +  comm - the MPI communicator
287: .  n - the length of the index set
288: -  idx - the list of integers

290:    Output Parameter:
291: .  is - the new index set

293:    Notes:
294:    The index array is copied to internally allocated storage. After the call,
295:    the user can free the index array. Use ISCreateGeneralNC() to use the pointers
296:    passed in and NOT make a copy of the index array.

298:    When the communicator is not MPI_COMM_SELF, the operations on IS are NOT
299:    conceptually the same as MPI_Group operations. The IS are then
300:    distributed sets of indices and thus certain operations on them are
301:    collective.

303:    Level: beginner

305:   Concepts: index sets^creating
306:   Concepts: IS^creating

308: .seealso: ISCreateGeneralWithArray(), ISCreateStride(), ISCreateBlock(), ISAllGather(), ISCreateGeneralNC()
309: @*/
310: PetscErrorCode  ISCreateGeneral(MPI_Comm comm,PetscInt n,const PetscInt idx[],IS *is)
311: {
313:   IS             Nindex;
314:   IS_General     *sub;

318:   if (n < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"length < 0");
320:   *is = PETSC_NULL;
321: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
322:   ISInitializePackage(PETSC_NULL);
323: #endif

325:   PetscHeaderCreate(Nindex,_p_IS,struct _ISOps,IS_COOKIE,IS_GENERAL,"IS",comm,ISDestroy,ISView);
326:   PetscNewLog(Nindex,IS_General,&sub);
327:   Nindex->data   = (void*)sub;
328:   PetscMalloc(n*sizeof(PetscInt),&sub->idx);
329:   PetscLogObjectMemory(Nindex,n*sizeof(PetscInt));
330:   PetscMemcpy(sub->idx,idx,n*sizeof(PetscInt));
331:   sub->n         = n;
332:   sub->allocated = PETSC_TRUE;

334:   *is = Nindex;
335:   ISCreateGeneral_Private(comm,is);

337:   return(0);
338: }

342: /*@C
343:    ISCreateGeneralNC - Creates a data structure for an index set 
344:    containing a list of integers.

346:    Collective on MPI_Comm

348:    Input Parameters:
349: +  comm - the MPI communicator
350: .  n - the length of the index set
351: -  idx - the list of integers

353:    Output Parameter:
354: .  is - the new index set

356:    Notes: This routine does not copy the indices, just keeps the pointer to the
357:    indices. The ISDestroy() will free the space so it must be obtained
358:    with PetscMalloc() and it must not be freed nor modified elsewhere.
359:    Use ISCreateGeneral() if you wish to copy the indices passed into the routine.
360:    Use ISCreateGeneralWithArray() to NOT copy the indices and NOT free the space when
361:    ISDestroy() is called.

363:    When the communicator is not MPI_COMM_SELF, the operations on IS are NOT
364:    conceptually the same as MPI_Group operations. The IS are then
365:    distributed sets of indices and thus certain operations on them are
366:    collective.

368:    Level: beginner

370:   Concepts: index sets^creating
371:   Concepts: IS^creating

373: .seealso: ISCreateGeneral(), ISCreateGeneralWithArray(), ISCreateStride(), ISCreateBlock(), ISAllGather()
374: @*/
375: PetscErrorCode  ISCreateGeneralNC(MPI_Comm comm,PetscInt n,const PetscInt idx[],IS *is)
376: {
378:   IS             Nindex;
379:   IS_General     *sub;

383:   if (n < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"length < 0");
385:   *is = PETSC_NULL;
386: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
387:   ISInitializePackage(PETSC_NULL);
388: #endif

390:   PetscHeaderCreate(Nindex,_p_IS,struct _ISOps,IS_COOKIE,IS_GENERAL,"IS",comm,ISDestroy,ISView);
391:   PetscNewLog(Nindex,IS_General,&sub);
392:   Nindex->data   = (void*)sub;
393:   sub->n         = n;
394:   sub->idx       = (PetscInt*)idx;
395:   sub->allocated = PETSC_TRUE;

397:   *is = Nindex;
398:   ISCreateGeneral_Private(comm,is);

400:   return(0);
401: }

405: /*@C
406:    ISCreateGeneralWithArray - Creates a data structure for an index set 
407:    containing a list of integers.

409:    Collective on MPI_Comm

411:    Input Parameters:
412: +  comm - the MPI communicator
413: .  n - the length of the index set
414: -  idx - the list of integers

416:    Output Parameter:
417: .  is - the new index set

419:    Notes:
420:    Unlike with ISCreateGeneral, the indices are not copied to internally
421:    allocated storage. The user array is not freed by ISDestroy.

423:    When the communicator is not MPI_COMM_SELF, the operations on IS are NOT
424:    conceptually the same as MPI_Group operations. The IS are then
425:    distributed sets of indices and thus certain operations on them are collective.

427:    Level: beginner

429:   Concepts: index sets^creating
430:   Concepts: IS^creating

432: .seealso: ISCreateGeneral(), ISCreateStride(), ISCreateBlock(), ISAllGather()
433: @*/
434: PetscErrorCode  ISCreateGeneralWithArray(MPI_Comm comm,PetscInt n,PetscInt idx[],IS *is)
435: {
437:   IS             Nindex;
438:   IS_General     *sub;

442:   if (n < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"length < 0");
444:   *is = PETSC_NULL;
445: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
446:   ISInitializePackage(PETSC_NULL);
447: #endif

449:   PetscHeaderCreate(Nindex,_p_IS,struct _ISOps,IS_COOKIE,IS_GENERAL,"IS",comm,ISDestroy,ISView);
450:   PetscNewLog(Nindex,IS_General,&sub);
451:   Nindex->data   = (void*)sub;
452:   sub->n         = n;
453:   sub->idx       = idx;
454:   sub->allocated = PETSC_FALSE;

456:   *is = Nindex;
457:   ISCreateGeneral_Private(comm,is);

459:   return(0);
460: }