Actual source code: iscoloring.c

petsc-3.8.4 2018-03-24
Report Typos and Errors

  2:  #include <petsc/private/isimpl.h>
  3:  #include <petscviewer.h>
  4:  #include <petscsf.h>

  6: const char *const ISColoringTypes[] = {"global","ghosted","ISColoringType","IS_COLORING_",0};

  8: PetscErrorCode ISColoringReference(ISColoring coloring)
  9: {
 11:   (coloring)->refct++;
 12:   return(0);
 13: }

 15: PetscErrorCode ISColoringSetType(ISColoring coloring,ISColoringType type)
 16: {
 18:   (coloring)->ctype = type;
 19:   return(0);
 20: }

 22: /*@
 23:    ISColoringDestroy - Destroys a coloring context.

 25:    Collective on ISColoring

 27:    Input Parameter:
 28: .  iscoloring - the coloring context

 30:    Level: advanced

 32: .seealso: ISColoringView(), MatColoring
 33: @*/
 34: PetscErrorCode  ISColoringDestroy(ISColoring *iscoloring)
 35: {
 36:   PetscInt       i;

 40:   if (!*iscoloring) return(0);
 42:   if (--(*iscoloring)->refct > 0) {*iscoloring = 0; return(0);}

 44:   if ((*iscoloring)->is) {
 45:     for (i=0; i<(*iscoloring)->n; i++) {
 46:       ISDestroy(&(*iscoloring)->is[i]);
 47:     }
 48:     PetscFree((*iscoloring)->is);
 49:   }
 50:   if ((*iscoloring)->allocated) {PetscFree((*iscoloring)->colors);}
 51:   PetscCommDestroy(&(*iscoloring)->comm);
 52:   PetscFree((*iscoloring));
 53:   return(0);
 54: }

 56: /*
 57:   ISColoringViewFromOptions - Processes command line options to determine if/how an ISColoring object is to be viewed. 

 59:   Collective on ISColoring

 61:   Input Parameters:
 62: + obj   - the ISColoring object
 63: . prefix - prefix to use for viewing, or NULL to use prefix of 'mat'
 64: - optionname - option to activate viewing

 66:   Level: intermediate

 68:   Developer Note: This cannot use PetscObjectViewFromOptions() because ISColoring is not a PetscObject

 70: */
 71: PetscErrorCode ISColoringViewFromOptions(ISColoring obj,PetscObject bobj,const char optionname[])
 72: {
 73:   PetscErrorCode    ierr;
 74:   PetscViewer       viewer;
 75:   PetscBool         flg;
 76:   PetscViewerFormat format;
 77:   char              *prefix;

 80:   prefix = bobj ? bobj->prefix : NULL;
 81:   PetscOptionsGetViewer(obj->comm,prefix,optionname,&viewer,&format,&flg);
 82:   if (flg) {
 83:     PetscViewerPushFormat(viewer,format);
 84:     ISColoringView(obj,viewer);
 85:     PetscViewerPopFormat(viewer);
 86:     PetscViewerDestroy(&viewer);
 87:   }
 88:   return(0);
 89: }

 91: /*@C
 92:    ISColoringView - Views a coloring context.

 94:    Collective on ISColoring

 96:    Input Parameters:
 97: +  iscoloring - the coloring context
 98: -  viewer - the viewer

100:    Level: advanced

102: .seealso: ISColoringDestroy(), ISColoringGetIS(), MatColoring
103: @*/
104: PetscErrorCode  ISColoringView(ISColoring iscoloring,PetscViewer viewer)
105: {
106:   PetscInt       i;
108:   PetscBool      iascii;
109:   IS             *is;

113:   if (!viewer) {
114:     PetscViewerASCIIGetStdout(iscoloring->comm,&viewer);
115:   }

118:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
119:   if (iascii) {
120:     MPI_Comm    comm;
121:     PetscMPIInt size,rank;

123:     PetscObjectGetComm((PetscObject)viewer,&comm);
124:     MPI_Comm_size(comm,&size);
125:     MPI_Comm_rank(comm,&rank);
126:     PetscViewerASCIIPrintf(viewer,"ISColoring Object: %d MPI processes\n",size);
127:     PetscViewerASCIIPushSynchronized(viewer);
128:     PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of colors %d\n",rank,iscoloring->n);
129:     PetscViewerFlush(viewer);
130:     PetscViewerASCIIPopSynchronized(viewer);
131:   }

133:   ISColoringGetIS(iscoloring,PETSC_IGNORE,&is);
134:   for (i=0; i<iscoloring->n; i++) {
135:     ISView(iscoloring->is[i],viewer);
136:   }
137:   ISColoringRestoreIS(iscoloring,&is);
138:   return(0);
139: }

141: /*@C
142:    ISColoringGetIS - Extracts index sets from the coloring context

144:    Collective on ISColoring

146:    Input Parameter:
147: .  iscoloring - the coloring context

149:    Output Parameters:
150: +  nn - number of index sets in the coloring context
151: -  is - array of index sets

153:    Level: advanced

155: .seealso: ISColoringRestoreIS(), ISColoringView()
156: @*/
157: PetscErrorCode  ISColoringGetIS(ISColoring iscoloring,PetscInt *nn,IS *isis[])
158: {


164:   if (nn) *nn = iscoloring->n;
165:   if (isis) {
166:     if (!iscoloring->is) {
167:       PetscInt        *mcolors,**ii,nc = iscoloring->n,i,base, n = iscoloring->N;
168:       ISColoringValue *colors = iscoloring->colors;
169:       IS              *is;

171: #if defined(PETSC_USE_DEBUG)
172:       for (i=0; i<n; i++) {
173:         if (((PetscInt)colors[i]) >= nc) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Coloring is our of range index %d value %d number colors %d",(int)i,(int)colors[i],(int)nc);
174:       }
175: #endif

177:       /* generate the lists of nodes for each color */
178:       PetscCalloc1(nc,&mcolors);
179:       for (i=0; i<n; i++) mcolors[colors[i]]++;

181:       PetscMalloc1(nc,&ii);
182:       PetscMalloc1(n,&ii[0]);
183:       for (i=1; i<nc; i++) ii[i] = ii[i-1] + mcolors[i-1];
184:       PetscMemzero(mcolors,nc*sizeof(PetscInt));

186:       if (iscoloring->ctype == IS_COLORING_GLOBAL) {
187:         MPI_Scan(&iscoloring->N,&base,1,MPIU_INT,MPI_SUM,iscoloring->comm);
188:         base -= iscoloring->N;
189:         for (i=0; i<n; i++) ii[colors[i]][mcolors[colors[i]]++] = i + base; /* global idx */
190:       } else if (iscoloring->ctype == IS_COLORING_LOCAL) {
191:         for (i=0; i<n; i++) ii[colors[i]][mcolors[colors[i]]++] = i;   /* local idx */
192:       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this ISColoringType type");

194:       PetscMalloc1(nc,&is);
195:       for (i=0; i<nc; i++) {
196:         ISCreateGeneral(iscoloring->comm,mcolors[i],ii[i],PETSC_COPY_VALUES,is+i);
197:       }

199:       iscoloring->is = is;
200:       PetscFree(ii[0]);
201:       PetscFree(ii);
202:       PetscFree(mcolors);
203:     }
204:     *isis = iscoloring->is;
205:   }
206:   return(0);
207: }

209: /*@C
210:    ISColoringRestoreIS - Restores the index sets extracted from the coloring context

212:    Collective on ISColoring

214:    Input Parameter:
215: +  iscoloring - the coloring context
216: -  is - array of index sets

218:    Level: advanced

220: .seealso: ISColoringGetIS(), ISColoringView()
221: @*/
222: PetscErrorCode  ISColoringRestoreIS(ISColoring iscoloring,IS *is[])
223: {

227:   /* currently nothing is done here */
228:   return(0);
229: }


232: /*@
233:     ISColoringCreate - Generates an ISColoring context from lists (provided
234:     by each processor) of colors for each node.

236:     Collective on MPI_Comm

238:     Input Parameters:
239: +   comm - communicator for the processors creating the coloring
240: .   ncolors - max color value
241: .   n - number of nodes on this processor
242: .   colors - array containing the colors for this processor, color numbers begin at 0.
243: -   mode - see PetscCopyMode for meaning of this flag.

245:     Output Parameter:
246: .   iscoloring - the resulting coloring data structure

248:     Options Database Key:
249: .   -is_coloring_view - Activates ISColoringView()

251:    Level: advanced

253:     Notes: By default sets coloring type to  IS_COLORING_GLOBAL

255: .seealso: MatColoringCreate(), ISColoringView(), ISColoringDestroy(), ISColoringSetType()

257: @*/
258: PetscErrorCode  ISColoringCreate(MPI_Comm comm,PetscInt ncolors,PetscInt n,const ISColoringValue colors[],PetscCopyMode mode,ISColoring *iscoloring)
259: {
261:   PetscMPIInt    size,rank,tag;
262:   PetscInt       base,top,i;
263:   PetscInt       nc,ncwork;
264:   MPI_Status     status;

267:   if (ncolors != PETSC_DECIDE && ncolors > IS_COLORING_MAX) {
268:     if (ncolors > 65535) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Max color value exeeds 65535 limit. This number is unrealistic. Perhaps a bug in code?\nCurrent max: %d user rewuested: %d",IS_COLORING_MAX,ncolors);
269:     else                 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Max color value exeeds limit. Perhaps reconfigure PETSc with --with-is-color-value-type=short?\n Current max: %d user rewuested: %d",IS_COLORING_MAX,ncolors);
270:   }
271:   PetscNew(iscoloring);
272:   PetscCommDuplicate(comm,&(*iscoloring)->comm,&tag);
273:   comm = (*iscoloring)->comm;

275:   /* compute the number of the first node on my processor */
276:   MPI_Comm_size(comm,&size);

278:   /* should use MPI_Scan() */
279:   MPI_Comm_rank(comm,&rank);
280:   if (!rank) {
281:     base = 0;
282:     top  = n;
283:   } else {
284:     MPI_Recv(&base,1,MPIU_INT,rank-1,tag,comm,&status);
285:     top  = base+n;
286:   }
287:   if (rank < size-1) {
288:     MPI_Send(&top,1,MPIU_INT,rank+1,tag,comm);
289:   }

291:   /* compute the total number of colors */
292:   ncwork = 0;
293:   for (i=0; i<n; i++) {
294:     if (ncwork < colors[i]) ncwork = colors[i];
295:   }
296:   ncwork++;
297:   MPIU_Allreduce(&ncwork,&nc,1,MPIU_INT,MPI_MAX,comm);
298:   if (nc > ncolors) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of colors passed in %D is less then the actual number of colors in array %D",ncolors,nc);
299:   (*iscoloring)->n      = nc;
300:   (*iscoloring)->is     = 0;
301:   (*iscoloring)->N      = n;
302:   (*iscoloring)->refct  = 1;
303:   (*iscoloring)->ctype  = IS_COLORING_GLOBAL;
304:   if (mode == PETSC_COPY_VALUES) {
305:     PetscMalloc1(n,&(*iscoloring)->colors);
306:     PetscLogObjectMemory((PetscObject)(*iscoloring),n*sizeof(ISColoringValue));
307:     PetscMemcpy((*iscoloring)->colors,colors,n*sizeof(ISColoringValue));
308:     (*iscoloring)->allocated = PETSC_TRUE;
309:   } else if (mode == PETSC_OWN_POINTER) {
310:     (*iscoloring)->colors    = (ISColoringValue*)colors;
311:     (*iscoloring)->allocated = PETSC_TRUE;
312:   } else {
313:     (*iscoloring)->colors    = (ISColoringValue*)colors;
314:     (*iscoloring)->allocated = PETSC_FALSE;
315:   }
316:   ISColoringViewFromOptions(*iscoloring,NULL,"-is_coloring_view");
317:   PetscInfo1(0,"Number of colors %D\n",nc);
318:   return(0);
319: }

321: /*@
322:     ISBuildTwoSided - Takes an IS that describes where we will go. Generates an IS that contains new numbers from remote or local
323:     on the IS.

325:     Collective on IS

327:     Input Parameters
328: .   to - an IS describes where we will go. Negative target rank will be ignored
329: .   toindx - an IS describes what indices should send. NULL means sending natural numbering

331:     Output Parameter:
332: .   rows - contains new numbers from remote or local

334:    Level: advanced

336: .seealso: MatPartitioningCreate(), ISPartitioningToNumbering(), ISPartitioningCount()

338: @*/
339: PetscErrorCode  ISBuildTwoSided(IS ito,IS toindx, IS *rows)
340: {
341:    const PetscInt       *ito_indices,*toindx_indices;
342:    PetscInt             *send_indices,rstart,*recv_indices,nrecvs,nsends;
343:    PetscInt             *tosizes,*fromsizes,i,j,*tosizes_tmp,*tooffsets_tmp,ito_ln;
344:    PetscMPIInt          *toranks,*fromranks,size,target_rank,*fromperm_newtoold,nto,nfrom;
345:    PetscLayout           isrmap;
346:    MPI_Comm              comm;
347:    PetscSF               sf;
348:    PetscSFNode          *iremote;
349:    PetscErrorCode        ierr;

352:    PetscObjectGetComm((PetscObject)ito,&comm);
353:    MPI_Comm_size(comm,&size);
354:    ISGetLocalSize(ito,&ito_ln);
355:    /* why we do not have ISGetLayout? */
356:    isrmap = ito->map;
357:    PetscLayoutGetRange(isrmap,&rstart,NULL);
358:    ISGetIndices(ito,&ito_indices);
359:    PetscCalloc2(size,&tosizes_tmp,size+1,&tooffsets_tmp);
360:    for(i=0; i<ito_ln; i++){
361:      if(ito_indices[i]<0) continue;
362: #if defined(PETSC_USE_DEBUG)
363:      if(ito_indices[i]>=size) SETERRQ2(comm,PETSC_ERR_ARG_OUTOFRANGE,"target rank %d is larger than communicator size %d ",ito_indices[i],size);
364: #endif
365:      tosizes_tmp[ito_indices[i]]++;
366:    }
367:    nto = 0;
368:    for(i=0; i<size; i++){
369:          tooffsets_tmp[i+1] = tooffsets_tmp[i]+tosizes_tmp[i];
370:      if(tosizes_tmp[i]>0) nto++;
371:     }
372:    PetscCalloc2(nto,&toranks,2*nto,&tosizes);
373:    nto = 0;
374:    for(i=0; i<size; i++){
375:      if(tosizes_tmp[i]>0){
376:         toranks[nto]      = i;
377:         tosizes[2*nto]    = tosizes_tmp[i];/* size */
378:         tosizes[2*nto+1]  = tooffsets_tmp[i];/* offset */
379:         nto++;
380:      }
381:    }
382:    nsends = tooffsets_tmp[size];
383:    PetscCalloc1(nsends,&send_indices);
384:    if(toindx){
385:          ISGetIndices(toindx,&toindx_indices);
386:    }
387:    for(i=0; i<ito_ln; i++){
388:          if(ito_indices[i]<0) continue;
389:          target_rank = ito_indices[i];
390:          send_indices[tooffsets_tmp[target_rank]] = toindx? toindx_indices[i]:(i+rstart);
391:          tooffsets_tmp[target_rank]++;
392:    }
393:    if(toindx){
394:             ISRestoreIndices(toindx,&toindx_indices);
395:    }
396:    ISRestoreIndices(ito,&ito_indices);
397:    PetscFree2(tosizes_tmp,tooffsets_tmp);
398:    PetscCommBuildTwoSided(comm,2,MPIU_INT,nto,toranks,tosizes,&nfrom,&fromranks,&fromsizes);
399:    PetscFree2(toranks,tosizes);
400:    PetscCalloc1(nfrom,&fromperm_newtoold);
401:    for(i=0; i<nfrom; i++){
402:          fromperm_newtoold[i] = i;
403:    }
404:    PetscSortMPIIntWithArray(nfrom,fromranks,fromperm_newtoold);
405:    nrecvs   = 0;
406:    for(i=0; i<nfrom; i++){
407:          nrecvs += fromsizes[i*2];
408:    }
409:    PetscCalloc1(nrecvs,&recv_indices);
410:    PetscCalloc1(nrecvs,&iremote);
411:    nrecvs = 0;
412:    for(i=0; i<nfrom; i++){
413:      for(j=0; j<fromsizes[2*fromperm_newtoold[i]]; j++){
414:        iremote[nrecvs].rank    = fromranks[i];
415:        iremote[nrecvs++].index = fromsizes[2*fromperm_newtoold[i]+1]+j;
416:      }
417:    }
418:    PetscSFCreate(comm,&sf);
419:    PetscSFSetGraph(sf,nsends,nrecvs,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);
420:    PetscSFSetType(sf,PETSCSFBASIC);
421:    /* how to put a prefix ? */
422:    PetscSFSetFromOptions(sf);
423:    PetscSFBcastBegin(sf,MPIU_INT,send_indices,recv_indices);
424:    PetscSFBcastEnd(sf,MPIU_INT,send_indices,recv_indices);
425:    PetscSFDestroy(&sf);
426:    PetscFree(fromranks);
427:    PetscFree(fromsizes);
428:    PetscFree(fromperm_newtoold);
429:    PetscFree(send_indices);
430:    if(rows){
431:          PetscSortInt(nrecvs,recv_indices);
432:      ISCreateGeneral(comm, nrecvs,recv_indices,PETSC_OWN_POINTER,rows);
433:    }else{
434:          PetscFree(recv_indices);
435:    }
436:    return(0);
437: }


440: /*@
441:     ISPartitioningToNumbering - Takes an ISPartitioning and on each processor
442:     generates an IS that contains a new global node number for each index based
443:     on the partitioing.

445:     Collective on IS

447:     Input Parameters
448: .   partitioning - a partitioning as generated by MatPartitioningApply()

450:     Output Parameter:
451: .   is - on each processor the index set that defines the global numbers
452:          (in the new numbering) for all the nodes currently (before the partitioning)
453:          on that processor

455:    Level: advanced

457: .seealso: MatPartitioningCreate(), AOCreateBasic(), ISPartitioningCount()

459: @*/
460: PetscErrorCode  ISPartitioningToNumbering(IS part,IS *is)
461: {
462:   MPI_Comm       comm;
463:   PetscInt       i,np,npt,n,*starts = NULL,*sums = NULL,*lsizes = NULL,*newi = NULL;
464:   const PetscInt *indices = NULL;

468:   PetscObjectGetComm((PetscObject)part,&comm);

470:   /* count the number of partitions, i.e., virtual processors */
471:   ISGetLocalSize(part,&n);
472:   ISGetIndices(part,&indices);
473:   np   = 0;
474:   for (i=0; i<n; i++) np = PetscMax(np,indices[i]);
475:   MPIU_Allreduce(&np,&npt,1,MPIU_INT,MPI_MAX,comm);
476:   np   = npt+1; /* so that it looks like a MPI_Comm_size output */

478:   /*
479:         lsizes - number of elements of each partition on this particular processor
480:         sums - total number of "previous" nodes for any particular partition
481:         starts - global number of first element in each partition on this processor
482:   */
483:   PetscMalloc3(np,&lsizes,np,&starts,np,&sums);
484:   PetscMemzero(lsizes,np*sizeof(PetscInt));
485:   for (i=0; i<n; i++) lsizes[indices[i]]++;
486:   MPIU_Allreduce(lsizes,sums,np,MPIU_INT,MPI_SUM,comm);
487:   MPI_Scan(lsizes,starts,np,MPIU_INT,MPI_SUM,comm);
488:   for (i=0; i<np; i++) starts[i] -= lsizes[i];
489:   for (i=1; i<np; i++) {
490:     sums[i]   += sums[i-1];
491:     starts[i] += sums[i-1];
492:   }

494:   /*
495:       For each local index give it the new global number
496:   */
497:   PetscMalloc1(n,&newi);
498:   for (i=0; i<n; i++) newi[i] = starts[indices[i]]++;
499:   PetscFree3(lsizes,starts,sums);

501:   ISRestoreIndices(part,&indices);
502:   ISCreateGeneral(comm,n,newi,PETSC_OWN_POINTER,is);
503:   ISSetPermutation(*is);
504:   return(0);
505: }

507: /*@
508:     ISPartitioningCount - Takes a ISPartitioning and determines the number of
509:     resulting elements on each (partition) process

511:     Collective on IS

513:     Input Parameters:
514: +   partitioning - a partitioning as generated by MatPartitioningApply()
515: -   len - length of the array count, this is the total number of partitions

517:     Output Parameter:
518: .   count - array of length size, to contain the number of elements assigned
519:         to each partition, where size is the number of partitions generated
520:          (see notes below).

522:    Level: advanced

524:     Notes:
525:         By default the number of partitions generated (and thus the length
526:         of count) is the size of the communicator associated with IS,
527:         but it can be set by MatPartitioningSetNParts. The resulting array
528:         of lengths can for instance serve as input of PCBJacobiSetTotalBlocks.


531: .seealso: MatPartitioningCreate(), AOCreateBasic(), ISPartitioningToNumbering(),
532:         MatPartitioningSetNParts(), MatPartitioningApply()

534: @*/
535: PetscErrorCode  ISPartitioningCount(IS part,PetscInt len,PetscInt count[])
536: {
537:   MPI_Comm       comm;
538:   PetscInt       i,n,*lsizes;
539:   const PetscInt *indices;
541:   PetscMPIInt    npp;

544:   PetscObjectGetComm((PetscObject)part,&comm);
545:   if (len == PETSC_DEFAULT) {
546:     PetscMPIInt size;
547:     MPI_Comm_size(comm,&size);
548:     len  = (PetscInt) size;
549:   }

551:   /* count the number of partitions */
552:   ISGetLocalSize(part,&n);
553:   ISGetIndices(part,&indices);
554: #if defined(PETSC_USE_DEBUG)
555:   {
556:     PetscInt np = 0,npt;
557:     for (i=0; i<n; i++) np = PetscMax(np,indices[i]);
558:     MPIU_Allreduce(&np,&npt,1,MPIU_INT,MPI_MAX,comm);
559:     np   = npt+1; /* so that it looks like a MPI_Comm_size output */
560:     if (np > len) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Length of count array %D is less than number of partitions %D",len,np);
561:   }
562: #endif

564:   /*
565:         lsizes - number of elements of each partition on this particular processor
566:         sums - total number of "previous" nodes for any particular partition
567:         starts - global number of first element in each partition on this processor
568:   */
569:   PetscCalloc1(len,&lsizes);
570:   for (i=0; i<n; i++) lsizes[indices[i]]++;
571:   ISRestoreIndices(part,&indices);
572:   PetscMPIIntCast(len,&npp);
573:   MPIU_Allreduce(lsizes,count,npp,MPIU_INT,MPI_SUM,comm);
574:   PetscFree(lsizes);
575:   return(0);
576: }

578: /*@
579:     ISAllGather - Given an index set (IS) on each processor, generates a large
580:     index set (same on each processor) by concatenating together each
581:     processors index set.

583:     Collective on IS

585:     Input Parameter:
586: .   is - the distributed index set

588:     Output Parameter:
589: .   isout - the concatenated index set (same on all processors)

591:     Notes:
592:     ISAllGather() is clearly not scalable for large index sets.

594:     The IS created on each processor must be created with a common
595:     communicator (e.g., PETSC_COMM_WORLD). If the index sets were created
596:     with PETSC_COMM_SELF, this routine will not work as expected, since
597:     each process will generate its own new IS that consists only of
598:     itself.

600:     The communicator for this new IS is PETSC_COMM_SELF

602:     Level: intermediate

604:     Concepts: gather^index sets
605:     Concepts: index sets^gathering to all processors
606:     Concepts: IS^gathering to all processors

608: .seealso: ISCreateGeneral(), ISCreateStride(), ISCreateBlock()
609: @*/
610: PetscErrorCode  ISAllGather(IS is,IS *isout)
611: {
613:   PetscInt       *indices,n,i,N,step,first;
614:   const PetscInt *lindices;
615:   MPI_Comm       comm;
616:   PetscMPIInt    size,*sizes = NULL,*offsets = NULL,nn;
617:   PetscBool      stride;


623:   PetscObjectGetComm((PetscObject)is,&comm);
624:   MPI_Comm_size(comm,&size);
625:   ISGetLocalSize(is,&n);
626:   PetscObjectTypeCompare((PetscObject)is,ISSTRIDE,&stride);
627:   if (size == 1 && stride) { /* should handle parallel ISStride also */
628:     ISStrideGetInfo(is,&first,&step);
629:     ISCreateStride(PETSC_COMM_SELF,n,first,step,isout);
630:   } else {
631:     PetscMalloc2(size,&sizes,size,&offsets);

633:     PetscMPIIntCast(n,&nn);
634:     MPI_Allgather(&nn,1,MPI_INT,sizes,1,MPI_INT,comm);
635:     offsets[0] = 0;
636:     for (i=1; i<size; i++) {
637:       PetscInt s = offsets[i-1] + sizes[i-1];
638:       PetscMPIIntCast(s,&offsets[i]);
639:     }
640:     N = offsets[size-1] + sizes[size-1];

642:     PetscMalloc1(N,&indices);
643:     ISGetIndices(is,&lindices);
644:     MPI_Allgatherv((void*)lindices,nn,MPIU_INT,indices,sizes,offsets,MPIU_INT,comm);
645:     ISRestoreIndices(is,&lindices);
646:     PetscFree2(sizes,offsets);

648:     ISCreateGeneral(PETSC_COMM_SELF,N,indices,PETSC_OWN_POINTER,isout);
649:   }
650:   return(0);
651: }

653: /*@C
654:     ISAllGatherColors - Given a a set of colors on each processor, generates a large
655:     set (same on each processor) by concatenating together each processors colors

657:     Collective on MPI_Comm

659:     Input Parameter:
660: +   comm - communicator to share the indices
661: .   n - local size of set
662: -   lindices - local colors

664:     Output Parameter:
665: +   outN - total number of indices
666: -   outindices - all of the colors

668:     Notes:
669:     ISAllGatherColors() is clearly not scalable for large index sets.


672:     Level: intermediate

674:     Concepts: gather^index sets
675:     Concepts: index sets^gathering to all processors
676:     Concepts: IS^gathering to all processors

678: .seealso: ISCreateGeneral(), ISCreateStride(), ISCreateBlock(), ISAllGather()
679: @*/
680: PetscErrorCode  ISAllGatherColors(MPI_Comm comm,PetscInt n,ISColoringValue *lindices,PetscInt *outN,ISColoringValue *outindices[])
681: {
682:   ISColoringValue *indices;
683:   PetscErrorCode  ierr;
684:   PetscInt        i,N;
685:   PetscMPIInt     size,*offsets = NULL,*sizes = NULL, nn = n;

688:   MPI_Comm_size(comm,&size);
689:   PetscMalloc2(size,&sizes,size,&offsets);

691:   MPI_Allgather(&nn,1,MPI_INT,sizes,1,MPI_INT,comm);
692:   offsets[0] = 0;
693:   for (i=1; i<size; i++) offsets[i] = offsets[i-1] + sizes[i-1];
694:   N    = offsets[size-1] + sizes[size-1];
695:   PetscFree2(sizes,offsets);

697:   PetscMalloc1(N+1,&indices);
698:   MPI_Allgatherv(lindices,(PetscMPIInt)n,MPIU_COLORING_VALUE,indices,sizes,offsets,MPIU_COLORING_VALUE,comm);

700:   *outindices = indices;
701:   if (outN) *outN = N;
702:   return(0);
703: }

705: /*@
706:     ISComplement - Given an index set (IS) generates the complement index set. That is all
707:        all indices that are NOT in the given set.

709:     Collective on IS

711:     Input Parameter:
712: +   is - the index set
713: .   nmin - the first index desired in the local part of the complement
714: -   nmax - the largest index desired in the local part of the complement (note that all indices in is must be greater or equal to nmin and less than nmax)

716:     Output Parameter:
717: .   isout - the complement

719:     Notes:  The communicator for this new IS is the same as for the input IS

721:       For a parallel IS, this will generate the local part of the complement on each process

723:       To generate the entire complement (on each process) of a parallel IS, first call ISAllGather() and then
724:     call this routine.

726:     Level: intermediate

728:     Concepts: gather^index sets
729:     Concepts: index sets^gathering to all processors
730:     Concepts: IS^gathering to all processors

732: .seealso: ISCreateGeneral(), ISCreateStride(), ISCreateBlock(), ISAllGather()
733: @*/
734: PetscErrorCode  ISComplement(IS is,PetscInt nmin,PetscInt nmax,IS *isout)
735: {
737:   const PetscInt *indices;
738:   PetscInt       n,i,j,unique,cnt,*nindices;
739:   PetscBool      sorted;

744:   if (nmin < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nmin %D cannot be negative",nmin);
745:   if (nmin > nmax) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nmin %D cannot be greater than nmax %D",nmin,nmax);
746:   ISSorted(is,&sorted);
747:   if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Index set must be sorted");

749:   ISGetLocalSize(is,&n);
750:   ISGetIndices(is,&indices);
751: #if defined(PETSC_USE_DEBUG)
752:   for (i=0; i<n; i++) {
753:     if (indices[i] <  nmin) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index %D's value %D is smaller than minimum given %D",i,indices[i],nmin);
754:     if (indices[i] >= nmax) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index %D's value %D is larger than maximum given %D",i,indices[i],nmax);
755:   }
756: #endif
757:   /* Count number of unique entries */
758:   unique = (n>0);
759:   for (i=0; i<n-1; i++) {
760:     if (indices[i+1] != indices[i]) unique++;
761:   }
762:   PetscMalloc1(nmax-nmin-unique,&nindices);
763:   cnt  = 0;
764:   for (i=nmin,j=0; i<nmax; i++) {
765:     if (j<n && i==indices[j]) do { j++; } while (j<n && i==indices[j]);
766:     else nindices[cnt++] = i;
767:   }
768:   if (cnt != nmax-nmin-unique) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of entries found in complement %D does not match expected %D",cnt,nmax-nmin-unique);
769:   ISCreateGeneral(PetscObjectComm((PetscObject)is),cnt,nindices,PETSC_OWN_POINTER,isout);
770:   ISRestoreIndices(is,&indices);
771:   return(0);
772: }