Actual source code: iscoloring.c

petsc-master 2019-07-16
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,NULL,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_USE_POINTER,PETSC_IGNORE,&is);
134:   for (i=0; i<iscoloring->n; i++) {
135:     ISView(iscoloring->is[i],viewer);
136:   }
137:   ISColoringRestoreIS(iscoloring,PETSC_USE_POINTER,&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
148: -  mode - if this value is PETSC_OWN_POINTER then the caller owns the pointer and must free the array of IS and each IS in the array

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

154:    Level: advanced

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


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

172: #if defined(PETSC_USE_DEBUG)
173:       for (i=0; i<n; i++) {
174:         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);
175:       }
176: #endif

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

182:       PetscMalloc1(nc,&ii);
183:       PetscMalloc1(n,&ii[0]);
184:       for (i=1; i<nc; i++) ii[i] = ii[i-1] + mcolors[i-1];
185:       PetscArrayzero(mcolors,nc);

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

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

200:       if (mode != PETSC_OWN_POINTER) iscoloring->is = is;
201:       *isis = is;
202:       PetscFree(ii[0]);
203:       PetscFree(ii);
204:       PetscFree(mcolors);
205:     } else {
206:       *isis = iscoloring->is;
207:     }
208:   }
209:   return(0);
210: }

212: /*@C
213:    ISColoringRestoreIS - Restores the index sets extracted from the coloring context

215:    Collective on ISColoring

217:    Input Parameter:
218: +  iscoloring - the coloring context
219: .  mode - who retains ownership of the is
220: -  is - array of index sets

222:    Level: advanced

224: .seealso: ISColoringGetIS(), ISColoringView()
225: @*/
226: PetscErrorCode  ISColoringRestoreIS(ISColoring iscoloring,PetscCopyMode mode,IS *is[])
227: {

231:   /* currently nothing is done here */
232:   return(0);
233: }


236: /*@
237:     ISColoringCreate - Generates an ISColoring context from lists (provided
238:     by each processor) of colors for each node.

240:     Collective

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

249:     Output Parameter:
250: .   iscoloring - the resulting coloring data structure

252:     Options Database Key:
253: .   -is_coloring_view - Activates ISColoringView()

255:    Level: advanced

257:     Notes:
258:     By default sets coloring type to  IS_COLORING_GLOBAL

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

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

272:   if (ncolors != PETSC_DECIDE && ncolors > IS_COLORING_MAX) {
273:     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);
274:     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);
275:   }
276:   PetscNew(iscoloring);
277:   PetscCommDuplicate(comm,&(*iscoloring)->comm,&tag);
278:   comm = (*iscoloring)->comm;

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

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

296:   /* compute the total number of colors */
297:   ncwork = 0;
298:   for (i=0; i<n; i++) {
299:     if (ncwork < colors[i]) ncwork = colors[i];
300:   }
301:   ncwork++;
302:   MPIU_Allreduce(&ncwork,&nc,1,MPIU_INT,MPI_MAX,comm);
303:   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);
304:   (*iscoloring)->n      = nc;
305:   (*iscoloring)->is     = 0;
306:   (*iscoloring)->N      = n;
307:   (*iscoloring)->refct  = 1;
308:   (*iscoloring)->ctype  = IS_COLORING_GLOBAL;
309:   if (mode == PETSC_COPY_VALUES) {
310:     PetscMalloc1(n,&(*iscoloring)->colors);
311:     PetscLogObjectMemory((PetscObject)(*iscoloring),n*sizeof(ISColoringValue));
312:     PetscArraycpy((*iscoloring)->colors,colors,n);
313:     (*iscoloring)->allocated = PETSC_TRUE;
314:   } else if (mode == PETSC_OWN_POINTER) {
315:     (*iscoloring)->colors    = (ISColoringValue*)colors;
316:     (*iscoloring)->allocated = PETSC_TRUE;
317:   } else {
318:     (*iscoloring)->colors    = (ISColoringValue*)colors;
319:     (*iscoloring)->allocated = PETSC_FALSE;
320:   }
321:   ISColoringViewFromOptions(*iscoloring,NULL,"-is_coloring_view");
322:   PetscInfo1(0,"Number of colors %D\n",nc);
323:   return(0);
324: }

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

330:     Collective on IS

332:     Input Parameters
333: +   ito - an IS describes where we will go. Negative target rank will be ignored
334: -   toindx - an IS describes what indices should send. NULL means sending natural numbering

336:     Output Parameter:
337: .   rows - contains new numbers from remote or local

339:    Level: advanced

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

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

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


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

446:     Collective on IS

448:     Input Parameters
449: .   partitioning - a partitioning as generated by MatPartitioningApply()
450:                    or MatPartitioningApplyND()

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

457:    Level: advanced

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

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

473:   /* see if the partitioning comes from nested dissection */
474:   PetscObjectQuery((PetscObject)part,"_petsc_matpartitioning_ndorder",(PetscObject*)&ndorder);
475:   if (ndorder) {
476:     PetscObjectReference((PetscObject)ndorder);
477:     *is  = ndorder;
478:     return(0);
479:   }

481:   PetscObjectGetComm((PetscObject)part,&comm);
482:   /* count the number of partitions, i.e., virtual processors */
483:   ISGetLocalSize(part,&n);
484:   ISGetIndices(part,&indices);
485:   np   = 0;
486:   for (i=0; i<n; i++) np = PetscMax(np,indices[i]);
487:   MPIU_Allreduce(&np,&npt,1,MPIU_INT,MPI_MAX,comm);
488:   np   = npt+1; /* so that it looks like a MPI_Comm_size output */

490:   /*
491:         lsizes - number of elements of each partition on this particular processor
492:         sums - total number of "previous" nodes for any particular partition
493:         starts - global number of first element in each partition on this processor
494:   */
495:   PetscMalloc3(np,&lsizes,np,&starts,np,&sums);
496:   PetscArrayzero(lsizes,np);
497:   for (i=0; i<n; i++) lsizes[indices[i]]++;
498:   MPIU_Allreduce(lsizes,sums,np,MPIU_INT,MPI_SUM,comm);
499:   MPI_Scan(lsizes,starts,np,MPIU_INT,MPI_SUM,comm);
500:   for (i=0; i<np; i++) starts[i] -= lsizes[i];
501:   for (i=1; i<np; i++) {
502:     sums[i]   += sums[i-1];
503:     starts[i] += sums[i-1];
504:   }

506:   /*
507:       For each local index give it the new global number
508:   */
509:   PetscMalloc1(n,&newi);
510:   for (i=0; i<n; i++) newi[i] = starts[indices[i]]++;
511:   PetscFree3(lsizes,starts,sums);

513:   ISRestoreIndices(part,&indices);
514:   ISCreateGeneral(comm,n,newi,PETSC_OWN_POINTER,is);
515:   ISSetPermutation(*is);
516:   return(0);
517: }

519: /*@
520:     ISPartitioningCount - Takes a ISPartitioning and determines the number of
521:     resulting elements on each (partition) process

523:     Collective on IS

525:     Input Parameters:
526: +   partitioning - a partitioning as generated by MatPartitioningApply() or
527:                    MatPartitioningApplyND()
528: -   len - length of the array count, this is the total number of partitions

530:     Output Parameter:
531: .   count - array of length size, to contain the number of elements assigned
532:         to each partition, where size is the number of partitions generated
533:          (see notes below).

535:    Level: advanced

537:     Notes:
538:         By default the number of partitions generated (and thus the length
539:         of count) is the size of the communicator associated with IS,
540:         but it can be set by MatPartitioningSetNParts. The resulting array
541:         of lengths can for instance serve as input of PCBJacobiSetTotalBlocks.
542:         If the partitioning has been obtained by MatPartitioningApplyND(),
543:         the returned count does not include the separators.

545: .seealso: MatPartitioningCreate(), AOCreateBasic(), ISPartitioningToNumbering(),
546:         MatPartitioningSetNParts(), MatPartitioningApply(), MatPartitioningApplyND()

548: @*/
549: PetscErrorCode  ISPartitioningCount(IS part,PetscInt len,PetscInt count[])
550: {
551:   MPI_Comm       comm;
552:   PetscInt       i,n,*lsizes;
553:   const PetscInt *indices;
555:   PetscMPIInt    npp;

558:   PetscObjectGetComm((PetscObject)part,&comm);
559:   if (len == PETSC_DEFAULT) {
560:     PetscMPIInt size;
561:     MPI_Comm_size(comm,&size);
562:     len  = (PetscInt) size;
563:   }

565:   /* count the number of partitions */
566:   ISGetLocalSize(part,&n);
567:   ISGetIndices(part,&indices);
568: #if defined(PETSC_USE_DEBUG)
569:   {
570:     PetscInt np = 0,npt;
571:     for (i=0; i<n; i++) np = PetscMax(np,indices[i]);
572:     MPIU_Allreduce(&np,&npt,1,MPIU_INT,MPI_MAX,comm);
573:     np   = npt+1; /* so that it looks like a MPI_Comm_size output */
574:     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);
575:   }
576: #endif

578:   /*
579:         lsizes - number of elements of each partition on this particular processor
580:         sums - total number of "previous" nodes for any particular partition
581:         starts - global number of first element in each partition on this processor
582:   */
583:   PetscCalloc1(len,&lsizes);
584:   for (i=0; i<n; i++) {
585:     if (indices[i] > -1) lsizes[indices[i]]++;
586:   }
587:   ISRestoreIndices(part,&indices);
588:   PetscMPIIntCast(len,&npp);
589:   MPIU_Allreduce(lsizes,count,npp,MPIU_INT,MPI_SUM,comm);
590:   PetscFree(lsizes);
591:   return(0);
592: }

594: /*@
595:     ISAllGather - Given an index set (IS) on each processor, generates a large
596:     index set (same on each processor) by concatenating together each
597:     processors index set.

599:     Collective on IS

601:     Input Parameter:
602: .   is - the distributed index set

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

607:     Notes:
608:     ISAllGather() is clearly not scalable for large index sets.

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

616:     The communicator for this new IS is PETSC_COMM_SELF

618:     Level: intermediate

620: .seealso: ISCreateGeneral(), ISCreateStride(), ISCreateBlock()
621: @*/
622: PetscErrorCode  ISAllGather(IS is,IS *isout)
623: {
625:   PetscInt       *indices,n,i,N,step,first;
626:   const PetscInt *lindices;
627:   MPI_Comm       comm;
628:   PetscMPIInt    size,*sizes = NULL,*offsets = NULL,nn;
629:   PetscBool      stride;


635:   PetscObjectGetComm((PetscObject)is,&comm);
636:   MPI_Comm_size(comm,&size);
637:   ISGetLocalSize(is,&n);
638:   PetscObjectTypeCompare((PetscObject)is,ISSTRIDE,&stride);
639:   if (size == 1 && stride) { /* should handle parallel ISStride also */
640:     ISStrideGetInfo(is,&first,&step);
641:     ISCreateStride(PETSC_COMM_SELF,n,first,step,isout);
642:   } else {
643:     PetscMalloc2(size,&sizes,size,&offsets);

645:     PetscMPIIntCast(n,&nn);
646:     MPI_Allgather(&nn,1,MPI_INT,sizes,1,MPI_INT,comm);
647:     offsets[0] = 0;
648:     for (i=1; i<size; i++) {
649:       PetscInt s = offsets[i-1] + sizes[i-1];
650:       PetscMPIIntCast(s,&offsets[i]);
651:     }
652:     N = offsets[size-1] + sizes[size-1];

654:     PetscMalloc1(N,&indices);
655:     ISGetIndices(is,&lindices);
656:     MPI_Allgatherv((void*)lindices,nn,MPIU_INT,indices,sizes,offsets,MPIU_INT,comm);
657:     ISRestoreIndices(is,&lindices);
658:     PetscFree2(sizes,offsets);

660:     ISCreateGeneral(PETSC_COMM_SELF,N,indices,PETSC_OWN_POINTER,isout);
661:   }
662:   return(0);
663: }

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

669:     Collective

671:     Input Parameter:
672: +   comm - communicator to share the indices
673: .   n - local size of set
674: -   lindices - local colors

676:     Output Parameter:
677: +   outN - total number of indices
678: -   outindices - all of the colors

680:     Notes:
681:     ISAllGatherColors() is clearly not scalable for large index sets.


684:     Level: intermediate

686: .seealso: ISCreateGeneral(), ISCreateStride(), ISCreateBlock(), ISAllGather()
687: @*/
688: PetscErrorCode  ISAllGatherColors(MPI_Comm comm,PetscInt n,ISColoringValue *lindices,PetscInt *outN,ISColoringValue *outindices[])
689: {
690:   ISColoringValue *indices;
691:   PetscErrorCode  ierr;
692:   PetscInt        i,N;
693:   PetscMPIInt     size,*offsets = NULL,*sizes = NULL, nn = n;

696:   MPI_Comm_size(comm,&size);
697:   PetscMalloc2(size,&sizes,size,&offsets);

699:   MPI_Allgather(&nn,1,MPI_INT,sizes,1,MPI_INT,comm);
700:   offsets[0] = 0;
701:   for (i=1; i<size; i++) offsets[i] = offsets[i-1] + sizes[i-1];
702:   N    = offsets[size-1] + sizes[size-1];
703:   PetscFree2(sizes,offsets);

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

708:   *outindices = indices;
709:   if (outN) *outN = N;
710:   return(0);
711: }

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

717:     Collective on IS

719:     Input Parameter:
720: +   is - the index set
721: .   nmin - the first index desired in the local part of the complement
722: -   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)

724:     Output Parameter:
725: .   isout - the complement

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

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

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

735:     Level: intermediate

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

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

754:   ISGetLocalSize(is,&n);
755:   ISGetIndices(is,&indices);
756: #if defined(PETSC_USE_DEBUG)
757:   for (i=0; i<n; i++) {
758:     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);
759:     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);
760:   }
761: #endif
762:   /* Count number of unique entries */
763:   unique = (n>0);
764:   for (i=0; i<n-1; i++) {
765:     if (indices[i+1] != indices[i]) unique++;
766:   }
767:   PetscMalloc1(nmax-nmin-unique,&nindices);
768:   cnt  = 0;
769:   for (i=nmin,j=0; i<nmax; i++) {
770:     if (j<n && i==indices[j]) do { j++; } while (j<n && i==indices[j]);
771:     else nindices[cnt++] = i;
772:   }
773:   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);
774:   ISCreateGeneral(PetscObjectComm((PetscObject)is),cnt,nindices,PETSC_OWN_POINTER,isout);
775:   ISRestoreIndices(is,&indices);
776:   return(0);
777: }