Actual source code: isltog.c

  1: #define PETSCVEC_DLL

 3:  #include petscvec.h
 4:  #include include/private/isimpl.h

  6: PetscCookie  IS_LTOGM_COOKIE = -1;

 10: /*@C
 11:     ISLocalToGlobalMappingGetSize - Gets the local size of a local to global mapping.

 13:     Not Collective

 15:     Input Parameter:
 16: .   ltog - local to global mapping

 18:     Output Parameter:
 19: .   n - the number of entries in the local mapping

 21:     Level: advanced

 23:     Concepts: mapping^local to global

 25: .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
 26: @*/
 27: PetscErrorCode  ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping mapping,PetscInt *n)
 28: {
 32:   *n = mapping->n;
 33:   return(0);
 34: }

 38: /*@C
 39:     ISLocalToGlobalMappingView - View a local to global mapping

 41:     Not Collective

 43:     Input Parameters:
 44: +   ltog - local to global mapping
 45: -   viewer - viewer

 47:     Level: advanced

 49:     Concepts: mapping^local to global

 51: .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
 52: @*/
 53: PetscErrorCode  ISLocalToGlobalMappingView(ISLocalToGlobalMapping mapping,PetscViewer viewer)
 54: {
 55:   PetscInt        i;
 56:   PetscMPIInt     rank;
 57:   PetscTruth      iascii;
 58:   PetscErrorCode  ierr;

 62:   if (!viewer) {
 63:     PetscViewerASCIIGetStdout(mapping->comm,&viewer);
 64:   }

 67:   MPI_Comm_rank(mapping->comm,&rank);
 68:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
 69:   if (iascii) {
 70:     for (i=0; i<mapping->n; i++) {
 71:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %d %d\n",rank,i,mapping->indices[i]);
 72:     }
 73:     PetscViewerFlush(viewer);
 74:   } else {
 75:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for ISLocalToGlobalMapping",((PetscObject)viewer)->type_name);
 76:   }

 78:   return(0);
 79: }

 83: /*@
 84:     ISLocalToGlobalMappingCreateIS - Creates a mapping between a local (0 to n)
 85:     ordering and a global parallel ordering.

 87:     Not collective

 89:     Input Parameter:
 90: .   is - index set containing the global numbers for each local

 92:     Output Parameter:
 93: .   mapping - new mapping data structure

 95:     Level: advanced

 97:     Concepts: mapping^local to global

 99: .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
100: @*/
101: PetscErrorCode  ISLocalToGlobalMappingCreateIS(IS is,ISLocalToGlobalMapping *mapping)
102: {
104:   PetscInt       n,*indices;
105:   MPI_Comm       comm;


111:   PetscObjectGetComm((PetscObject)is,&comm);
112:   ISGetLocalSize(is,&n);
113:   ISGetIndices(is,&indices);
114:   ISLocalToGlobalMappingCreate(comm,n,indices,mapping);
115:   ISRestoreIndices(is,&indices);

117:   return(0);
118: }


123: /*@
124:     ISLocalToGlobalMappingCreate - Creates a mapping between a local (0 to n)
125:     ordering and a global parallel ordering.

127:     Not Collective, but communicator may have more than one process

129:     Input Parameters:
130: +   comm - MPI communicator
131: .   n - the number of local elements
132: -   indices - the global index for each local element

134:     Output Parameter:
135: .   mapping - new mapping data structure

137:     Level: advanced

139:     Concepts: mapping^local to global

141: .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreateNC()
142: @*/
143: PetscErrorCode  ISLocalToGlobalMappingCreate(MPI_Comm cm,PetscInt n,const PetscInt indices[],ISLocalToGlobalMapping *mapping)
144: {
146:   PetscInt       *in;

151:   PetscMalloc(n*sizeof(PetscInt),&in);
152:   PetscMemcpy(in,indices,n*sizeof(PetscInt));
153:   ISLocalToGlobalMappingCreateNC(cm,n,in,mapping);
154:   return(0);
155: }

159: /*@C
160:     ISLocalToGlobalMappingCreateNC - Creates a mapping between a local (0 to n)
161:     ordering and a global parallel ordering.

163:     Not Collective, but communicator may have more than one process

165:     Input Parameters:
166: +   comm - MPI communicator
167: .   n - the number of local elements
168: -   indices - the global index for each local element

170:     Output Parameter:
171: .   mapping - new mapping data structure

173:     Level: developer

175:     Notes: Does not copy the indices, just keeps the pointer to the indices. The ISLocalToGlobalMappingDestroy()
176:     will free the space so it must be obtained with PetscMalloc() and it must not be freed elsewhere.

178:     Concepts: mapping^local to global

180: .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate()
181: @*/
182: PetscErrorCode  ISLocalToGlobalMappingCreateNC(MPI_Comm cm,PetscInt n,const PetscInt indices[],ISLocalToGlobalMapping *mapping)
183: {

187:   if (n) {
189:   }
191:   *mapping = PETSC_NULL;
192: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
193:   VecInitializePackage(PETSC_NULL);
194: #endif
195:   if (IS_LTOGM_COOKIE == -1) {
196:     PetscLogClassRegister(&IS_LTOGM_COOKIE,"IS L to G Mapping");
197:   }

199:   PetscHeaderCreate(*mapping,_p_ISLocalToGlobalMapping,int,IS_LTOGM_COOKIE,0,"ISLocalToGlobalMapping",
200:                     cm,ISLocalToGlobalMappingDestroy,ISLocalToGlobalMappingView);
201:   PetscLogObjectMemory(*mapping,sizeof(struct _p_ISLocalToGlobalMapping)+n*sizeof(PetscInt));

203:   (*mapping)->n       = n;
204:   (*mapping)->indices = (PetscInt*)indices;

206:   /*
207:       Do not create the global to local mapping. This is only created if 
208:      ISGlobalToLocalMapping() is called 
209:   */
210:   (*mapping)->globals = 0;
211:   return(0);
212: }

216: /*@
217:     ISLocalToGlobalMappingBlock - Creates a blocked index version of an 
218:        ISLocalToGlobalMapping that is appropriate for MatSetLocalToGlobalMappingBlock()
219:        and VecSetLocalToGlobalMappingBlock().

221:     Not Collective, but communicator may have more than one process

223:     Input Parameters:
224: +    inmap - original point-wise mapping
225: -    bs - block size

227:     Output Parameter:
228: .   outmap - block based mapping; the indices are relative to BLOCKS, not individual vector or matrix entries.

230:     Level: advanced

232:     Concepts: mapping^local to global

234: .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingCreateIS()
235: @*/
236: PetscErrorCode  ISLocalToGlobalMappingBlock(ISLocalToGlobalMapping inmap,PetscInt bs,ISLocalToGlobalMapping *outmap)
237: {
239:   PetscInt       *ii,i,n;

242:   if (bs > 1) {
243:     n    = inmap->n/bs;
244:     if (n*bs != inmap->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Pointwise mapping length is not divisible by block size");
245:     PetscMalloc(n*sizeof(PetscInt),&ii);
246:     for (i=0; i<n; i++) {
247:       ii[i] = inmap->indices[bs*i]/bs;
248:     }
249:     ISLocalToGlobalMappingCreate(inmap->comm,n,ii,outmap);
250:     PetscFree(ii);
251:   } else {
252:     PetscObjectReference((PetscObject)inmap);
253:     *outmap = inmap;
254:   }
255:   return(0);
256: }
257: 
260: /*@
261:    ISLocalToGlobalMappingDestroy - Destroys a mapping between a local (0 to n)
262:    ordering and a global parallel ordering.

264:    Note Collective

266:    Input Parameters:
267: .  mapping - mapping data structure

269:    Level: advanced

271: .seealso: ISLocalToGlobalMappingCreate()
272: @*/
273: PetscErrorCode  ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping mapping)
274: {
278:   if (--mapping->refct > 0) return(0);
279:   if (mapping->refct < 0) {
280:     SETERRQ(PETSC_ERR_PLIB,"Mapping already destroyed");
281:   }

283:   PetscFree(mapping->indices);
284:   PetscFree(mapping->globals);
285:   PetscHeaderDestroy(mapping);
286:   return(0);
287: }
288: 
291: /*@
292:     ISLocalToGlobalMappingApplyIS - Creates from an IS in the local numbering
293:     a new index set using the global numbering defined in an ISLocalToGlobalMapping
294:     context.

296:     Not collective

298:     Input Parameters:
299: +   mapping - mapping between local and global numbering
300: -   is - index set in local numbering

302:     Output Parameters:
303: .   newis - index set in global numbering

305:     Level: advanced

307:     Concepts: mapping^local to global

309: .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),
310:           ISLocalToGlobalMappingDestroy(), ISGlobalToLocalMappingApply()
311: @*/
312: PetscErrorCode  ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping,IS is,IS *newis)
313: {
315:   PetscInt       n,i,*idxin,*idxmap,*idxout,Nmax = mapping->n;


322:   ISGetLocalSize(is,&n);
323:   ISGetIndices(is,&idxin);
324:   idxmap = mapping->indices;
325: 
326:   PetscMalloc(n*sizeof(PetscInt),&idxout);
327:   for (i=0; i<n; i++) {
328:     if (idxin[i] >= Nmax) SETERRQ3(PETSC_ERR_ARG_OUTOFRANGE,"Local index %d too large %d (max) at %d",idxin[i],Nmax-1,i);
329:     idxout[i] = idxmap[idxin[i]];
330:   }
331:   ISRestoreIndices(is,&idxin);
332:   ISCreateGeneralNC(PETSC_COMM_SELF,n,idxout,newis);
333:   return(0);
334: }

336: /*MC
337:    ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering
338:    and converts them to the global numbering.

340:    Not collective

342:    Input Parameters:
343: +  mapping - the local to global mapping context
344: .  N - number of integers
345: -  in - input indices in local numbering

347:    Output Parameter:
348: .  out - indices in global numbering

350:    Synopsis:
351:    PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping,int N,int in[],int out[])

353:    Notes: 
354:    The in and out array parameters may be identical.

356:    Level: advanced

358: .seealso: ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(), 
359:           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
360:           AOPetscToApplication(), ISGlobalToLocalMappingApply()

362:     Concepts: mapping^local to global

364: M*/

366: /* -----------------------------------------------------------------------------------------*/

370: /*
371:     Creates the global fields in the ISLocalToGlobalMapping structure
372: */
373: static PetscErrorCode ISGlobalToLocalMappingSetUp_Private(ISLocalToGlobalMapping mapping)
374: {
376:   PetscInt       i,*idx = mapping->indices,n = mapping->n,end,start,*globals;

379:   end   = 0;
380:   start = 100000000;

382:   for (i=0; i<n; i++) {
383:     if (idx[i] < 0) continue;
384:     if (idx[i] < start) start = idx[i];
385:     if (idx[i] > end)   end   = idx[i];
386:   }
387:   if (start > end) {start = 0; end = -1;}
388:   mapping->globalstart = start;
389:   mapping->globalend   = end;

391:   PetscMalloc((end-start+2)*sizeof(PetscInt),&globals);
392:   mapping->globals = globals;
393:   for (i=0; i<end-start+1; i++) {
394:     globals[i] = -1;
395:   }
396:   for (i=0; i<n; i++) {
397:     if (idx[i] < 0) continue;
398:     globals[idx[i] - start] = i;
399:   }

401:   PetscLogObjectMemory(mapping,(end-start+1)*sizeof(PetscInt));
402:   return(0);
403: }

407: /*@
408:     ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers
409:     specified with a global numbering.

411:     Not collective

413:     Input Parameters:
414: +   mapping - mapping between local and global numbering
415: .   type - IS_GTOLM_MASK - replaces global indices with no local value with -1
416:            IS_GTOLM_DROP - drops the indices with no local value from the output list
417: .   n - number of global indices to map
418: -   idx - global indices to map

420:     Output Parameters:
421: +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
422: -   idxout - local index of each global index, one must pass in an array long enough 
423:              to hold all the indices. You can call ISGlobalToLocalMappingApply() with 
424:              idxout == PETSC_NULL to determine the required length (returned in nout)
425:              and then allocate the required space and call ISGlobalToLocalMappingApply()
426:              a second time to set the values.

428:     Notes:
429:     Either nout or idxout may be PETSC_NULL. idx and idxout may be identical.

431:     This is not scalable in memory usage. Each processor requires O(Nglobal) size 
432:     array to compute these.

434:     Level: advanced

436:     Concepts: mapping^global to local

438: .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),
439:           ISLocalToGlobalMappingDestroy()
440: @*/
441: PetscErrorCode  ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingType type,
442:                                   PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
443: {
444:   PetscInt       i,*globals,nf = 0,tmp,start,end;

448:   if (!mapping->globals) {
449:     ISGlobalToLocalMappingSetUp_Private(mapping);
450:   }
451:   globals = mapping->globals;
452:   start   = mapping->globalstart;
453:   end     = mapping->globalend;

455:   if (type == IS_GTOLM_MASK) {
456:     if (idxout) {
457:       for (i=0; i<n; i++) {
458:         if (idx[i] < 0) idxout[i] = idx[i];
459:         else if (idx[i] < start) idxout[i] = -1;
460:         else if (idx[i] > end)   idxout[i] = -1;
461:         else                     idxout[i] = globals[idx[i] - start];
462:       }
463:     }
464:     if (nout) *nout = n;
465:   } else {
466:     if (idxout) {
467:       for (i=0; i<n; i++) {
468:         if (idx[i] < 0) continue;
469:         if (idx[i] < start) continue;
470:         if (idx[i] > end) continue;
471:         tmp = globals[idx[i] - start];
472:         if (tmp < 0) continue;
473:         idxout[nf++] = tmp;
474:       }
475:     } else {
476:       for (i=0; i<n; i++) {
477:         if (idx[i] < 0) continue;
478:         if (idx[i] < start) continue;
479:         if (idx[i] > end) continue;
480:         tmp = globals[idx[i] - start];
481:         if (tmp < 0) continue;
482:         nf++;
483:       }
484:     }
485:     if (nout) *nout = nf;
486:   }

488:   return(0);
489: }

493: /*@C
494:     ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each processor and 
495:      each index shared by more than one processor 

497:     Collective on ISLocalToGlobalMapping

499:     Input Parameters:
500: .   mapping - the mapping from local to global indexing

502:     Output Parameter:
503: +   nproc - number of processors that are connected to this one
504: .   proc - neighboring processors
505: .   numproc - number of indices for each subdomain (processor)
506: -   indices - indices of local nodes shared with neighbor (sorted by global numbering)

508:     Level: advanced

510:     Concepts: mapping^local to global

512:     Fortran Usage: 
513: $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by 
514: $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
515:           PetscInt indices[nproc][numprocmax],ierr)
516:         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and 
517:         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.


520: .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
521:           ISLocalToGlobalMappingRestoreInfo()
522: @*/
523: PetscErrorCode  ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
524: {
526:   PetscMPIInt    size,rank,tag1,tag2,tag3,*len,*source,imdex;
527:   PetscInt       i,n = mapping->n,Ng,ng,max = 0,*lindices = mapping->indices;
528:   PetscInt       *nprocs,*owner,nsends,*sends,j,*starts,nmax,nrecvs,*recvs,proc;
529:   PetscInt       cnt,scale,*ownedsenders,*nownedsenders,rstart,nowned;
530:   PetscInt       node,nownedm,nt,*sends2,nsends2,*starts2,*lens2,*dest,nrecvs2,*starts3,*recvs2,k,*bprocs,*tmp;
531:   PetscInt       first_procs,first_numprocs,*first_indices;
532:   MPI_Request    *recv_waits,*send_waits;
533:   MPI_Status     recv_status,*send_status,*recv_statuses;
534:   MPI_Comm       comm = mapping->comm;
535:   PetscTruth     debug = PETSC_FALSE;

538:   MPI_Comm_size(comm,&size);
539:   MPI_Comm_rank(comm,&rank);
540:   if (size == 1) {
541:     *nproc         = 0;
542:     *procs         = PETSC_NULL;
543:     PetscMalloc(sizeof(PetscInt),numprocs);
544:     (*numprocs)[0] = 0;
545:     PetscMalloc(sizeof(PetscInt*),indices);
546:     (*indices)[0]  = PETSC_NULL;
547:     return(0);
548:   }

550:   PetscOptionsHasName(PETSC_NULL,"-islocaltoglobalmappinggetinfo_debug",&debug);

552:   /*
553:     Notes on ISLocalToGlobalMappingGetInfo

555:     globally owned node - the nodes that have been assigned to this processor in global
556:            numbering, just for this routine.

558:     nontrivial globally owned node - node assigned to this processor that is on a subdomain
559:            boundary (i.e. is has more than one local owner)

561:     locally owned node - node that exists on this processors subdomain

563:     nontrivial locally owned node - node that is not in the interior (i.e. has more than one
564:            local subdomain
565:   */
566:   PetscObjectGetNewTag((PetscObject)mapping,&tag1);
567:   PetscObjectGetNewTag((PetscObject)mapping,&tag2);
568:   PetscObjectGetNewTag((PetscObject)mapping,&tag3);

570:   for (i=0; i<n; i++) {
571:     if (lindices[i] > max) max = lindices[i];
572:   }
573:   MPI_Allreduce(&max,&Ng,1,MPIU_INT,MPI_MAX,comm);
574:   Ng++;
575:   MPI_Comm_size(comm,&size);
576:   MPI_Comm_rank(comm,&rank);
577:   scale  = Ng/size + 1;
578:   ng     = scale; if (rank == size-1) ng = Ng - scale*(size-1); ng = PetscMax(1,ng);
579:   rstart = scale*rank;

581:   /* determine ownership ranges of global indices */
582:   PetscMalloc(2*size*sizeof(PetscInt),&nprocs);
583:   PetscMemzero(nprocs,2*size*sizeof(PetscInt));

585:   /* determine owners of each local node  */
586:   PetscMalloc(n*sizeof(PetscInt),&owner);
587:   for (i=0; i<n; i++) {
588:     proc             = lindices[i]/scale; /* processor that globally owns this index */
589:     nprocs[2*proc+1] = 1;                 /* processor globally owns at least one of ours */
590:     owner[i]         = proc;
591:     nprocs[2*proc]++;                     /* count of how many that processor globally owns of ours */
592:   }
593:   nsends = 0; for (i=0; i<size; i++) nsends += nprocs[2*i+1];
594:   PetscInfo1(0,"Number of global owners for my local data %d\n",nsends);

596:   /* inform other processors of number of messages and max length*/
597:   PetscMaxSum(comm,nprocs,&nmax,&nrecvs);
598:   PetscInfo1(0,"Number of local owners for my global data %d\n",nrecvs);

600:   /* post receives for owned rows */
601:   PetscMalloc((2*nrecvs+1)*(nmax+1)*sizeof(PetscInt),&recvs);
602:   PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);
603:   for (i=0; i<nrecvs; i++) {
604:     MPI_Irecv(recvs+2*nmax*i,2*nmax,MPIU_INT,MPI_ANY_SOURCE,tag1,comm,recv_waits+i);
605:   }

607:   /* pack messages containing lists of local nodes to owners */
608:   PetscMalloc((2*n+1)*sizeof(PetscInt),&sends);
609:   PetscMalloc((size+1)*sizeof(PetscInt),&starts);
610:   starts[0]  = 0;
611:   for (i=1; i<size; i++) { starts[i] = starts[i-1] + 2*nprocs[2*i-2];}
612:   for (i=0; i<n; i++) {
613:     sends[starts[owner[i]]++] = lindices[i];
614:     sends[starts[owner[i]]++] = i;
615:   }
616:   PetscFree(owner);
617:   starts[0]  = 0;
618:   for (i=1; i<size; i++) { starts[i] = starts[i-1] + 2*nprocs[2*i-2];}

620:   /* send the messages */
621:   PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);
622:   PetscMalloc((nsends+1)*sizeof(PetscInt),&dest);
623:   cnt = 0;
624:   for (i=0; i<size; i++) {
625:     if (nprocs[2*i]) {
626:       MPI_Isend(sends+starts[i],2*nprocs[2*i],MPIU_INT,i,tag1,comm,send_waits+cnt);
627:       dest[cnt] = i;
628:       cnt++;
629:     }
630:   }
631:   PetscFree(starts);

633:   /* wait on receives */
634:   PetscMalloc((nrecvs+1)*sizeof(PetscMPIInt),&source);
635:   PetscMalloc((nrecvs+1)*sizeof(PetscMPIInt),&len);
636:   cnt  = nrecvs;
637:   PetscMalloc((ng+1)*sizeof(PetscInt),&nownedsenders);
638:   PetscMemzero(nownedsenders,ng*sizeof(PetscInt));
639:   while (cnt) {
640:     MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);
641:     /* unpack receives into our local space */
642:     MPI_Get_count(&recv_status,MPIU_INT,&len[imdex]);
643:     source[imdex]  = recv_status.MPI_SOURCE;
644:     len[imdex]     = len[imdex]/2;
645:     /* count how many local owners for each of my global owned indices */
646:     for (i=0; i<len[imdex]; i++) nownedsenders[recvs[2*imdex*nmax+2*i]-rstart]++;
647:     cnt--;
648:   }
649:   PetscFree(recv_waits);

651:   /* count how many globally owned indices are on an edge multiplied by how many processors own them. */
652:   nowned  = 0;
653:   nownedm = 0;
654:   for (i=0; i<ng; i++) {
655:     if (nownedsenders[i] > 1) {nownedm += nownedsenders[i]; nowned++;}
656:   }

658:   /* create single array to contain rank of all local owners of each globally owned index */
659:   PetscMalloc((nownedm+1)*sizeof(PetscInt),&ownedsenders);
660:   PetscMalloc((ng+1)*sizeof(PetscInt),&starts);
661:   starts[0] = 0;
662:   for (i=1; i<ng; i++) {
663:     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
664:     else starts[i] = starts[i-1];
665:   }

667:   /* for each nontrival globally owned node list all arriving processors */
668:   for (i=0; i<nrecvs; i++) {
669:     for (j=0; j<len[i]; j++) {
670:       node = recvs[2*i*nmax+2*j]-rstart;
671:       if (nownedsenders[node] > 1) {
672:         ownedsenders[starts[node]++] = source[i];
673:       }
674:     }
675:   }

677:   if (debug) { /* -----------------------------------  */
678:     starts[0]    = 0;
679:     for (i=1; i<ng; i++) {
680:       if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
681:       else starts[i] = starts[i-1];
682:     }
683:     for (i=0; i<ng; i++) {
684:       if (nownedsenders[i] > 1) {
685:         PetscSynchronizedPrintf(comm,"[%d] global node %d local owner processors: ",rank,i+rstart);
686:         for (j=0; j<nownedsenders[i]; j++) {
687:           PetscSynchronizedPrintf(comm,"%d ",ownedsenders[starts[i]+j]);
688:         }
689:         PetscSynchronizedPrintf(comm,"\n");
690:       }
691:     }
692:     PetscSynchronizedFlush(comm);
693:   }/* -----------------------------------  */

695:   /* wait on original sends */
696:   if (nsends) {
697:     PetscMalloc(nsends*sizeof(MPI_Status),&send_status);
698:     MPI_Waitall(nsends,send_waits,send_status);
699:     PetscFree(send_status);
700:   }
701:   PetscFree(send_waits);
702:   PetscFree(sends);
703:   PetscFree(nprocs);

705:   /* pack messages to send back to local owners */
706:   starts[0]    = 0;
707:   for (i=1; i<ng; i++) {
708:     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
709:     else starts[i] = starts[i-1];
710:   }
711:   nsends2 = nrecvs;
712:   PetscMalloc((nsends2+1)*sizeof(PetscInt),&nprocs); /* length of each message */
713:   for (i=0; i<nrecvs; i++) {
714:     nprocs[i] = 1;
715:     for (j=0; j<len[i]; j++) {
716:       node = recvs[2*i*nmax+2*j]-rstart;
717:       if (nownedsenders[node] > 1) {
718:         nprocs[i] += 2 + nownedsenders[node];
719:       }
720:     }
721:   }
722:   nt = 0; for (i=0; i<nsends2; i++) nt += nprocs[i];
723:   PetscMalloc((nt+1)*sizeof(PetscInt),&sends2);
724:   PetscMalloc((nsends2+1)*sizeof(PetscInt),&starts2);
725:   starts2[0] = 0; for (i=1; i<nsends2; i++) starts2[i] = starts2[i-1] + nprocs[i-1];
726:   /*
727:      Each message is 1 + nprocs[i] long, and consists of 
728:        (0) the number of nodes being sent back 
729:        (1) the local node number,
730:        (2) the number of processors sharing it,
731:        (3) the processors sharing it
732:   */
733:   for (i=0; i<nsends2; i++) {
734:     cnt = 1;
735:     sends2[starts2[i]] = 0;
736:     for (j=0; j<len[i]; j++) {
737:       node = recvs[2*i*nmax+2*j]-rstart;
738:       if (nownedsenders[node] > 1) {
739:         sends2[starts2[i]]++;
740:         sends2[starts2[i]+cnt++] = recvs[2*i*nmax+2*j+1];
741:         sends2[starts2[i]+cnt++] = nownedsenders[node];
742:         PetscMemcpy(&sends2[starts2[i]+cnt],&ownedsenders[starts[node]],nownedsenders[node]*sizeof(PetscInt));
743:         cnt += nownedsenders[node];
744:       }
745:     }
746:   }

748:   /* receive the message lengths */
749:   nrecvs2 = nsends;
750:   PetscMalloc((nrecvs2+1)*sizeof(PetscInt),&lens2);
751:   PetscMalloc((nrecvs2+1)*sizeof(PetscInt),&starts3);
752:   PetscMalloc((nrecvs2+1)*sizeof(MPI_Request),&recv_waits);
753:   for (i=0; i<nrecvs2; i++) {
754:     MPI_Irecv(&lens2[i],1,MPIU_INT,dest[i],tag2,comm,recv_waits+i);
755:   }

757:   /* send the message lengths */
758:   for (i=0; i<nsends2; i++) {
759:     MPI_Send(&nprocs[i],1,MPIU_INT,source[i],tag2,comm);
760:   }

762:   /* wait on receives of lens */
763:   if (nrecvs2) {
764:     PetscMalloc(nrecvs2*sizeof(MPI_Status),&recv_statuses);
765:     MPI_Waitall(nrecvs2,recv_waits,recv_statuses);
766:     PetscFree(recv_statuses);
767:   }
768:   PetscFree(recv_waits);

770:   starts3[0] = 0;
771:   nt         = 0;
772:   for (i=0; i<nrecvs2-1; i++) {
773:     starts3[i+1] = starts3[i] + lens2[i];
774:     nt          += lens2[i];
775:   }
776:   nt += lens2[nrecvs2-1];

778:   PetscMalloc((nt+1)*sizeof(PetscInt),&recvs2);
779:   PetscMalloc((nrecvs2+1)*sizeof(MPI_Request),&recv_waits);
780:   for (i=0; i<nrecvs2; i++) {
781:     MPI_Irecv(recvs2+starts3[i],lens2[i],MPIU_INT,dest[i],tag3,comm,recv_waits+i);
782:   }
783: 
784:   /* send the messages */
785:   PetscMalloc((nsends2+1)*sizeof(MPI_Request),&send_waits);
786:   for (i=0; i<nsends2; i++) {
787:     MPI_Isend(sends2+starts2[i],nprocs[i],MPIU_INT,source[i],tag3,comm,send_waits+i);
788:   }

790:   /* wait on receives */
791:   if (nrecvs2) {
792:     PetscMalloc(nrecvs2*sizeof(MPI_Status),&recv_statuses);
793:     MPI_Waitall(nrecvs2,recv_waits,recv_statuses);
794:     PetscFree(recv_statuses);
795:   }
796:   PetscFree(recv_waits);
797:   PetscFree(nprocs);

799:   if (debug) { /* -----------------------------------  */
800:     cnt = 0;
801:     for (i=0; i<nrecvs2; i++) {
802:       nt = recvs2[cnt++];
803:       for (j=0; j<nt; j++) {
804:         PetscSynchronizedPrintf(comm,"[%d] local node %d number of subdomains %d: ",rank,recvs2[cnt],recvs2[cnt+1]);
805:         for (k=0; k<recvs2[cnt+1]; k++) {
806:           PetscSynchronizedPrintf(comm,"%d ",recvs2[cnt+2+k]);
807:         }
808:         cnt += 2 + recvs2[cnt+1];
809:         PetscSynchronizedPrintf(comm,"\n");
810:       }
811:     }
812:     PetscSynchronizedFlush(comm);
813:   } /* -----------------------------------  */

815:   /* count number subdomains for each local node */
816:   PetscMalloc(size*sizeof(PetscInt),&nprocs);
817:   PetscMemzero(nprocs,size*sizeof(PetscInt));
818:   cnt  = 0;
819:   for (i=0; i<nrecvs2; i++) {
820:     nt = recvs2[cnt++];
821:     for (j=0; j<nt; j++) {
822:       for (k=0; k<recvs2[cnt+1]; k++) {
823:         nprocs[recvs2[cnt+2+k]]++;
824:       }
825:       cnt += 2 + recvs2[cnt+1];
826:     }
827:   }
828:   nt = 0; for (i=0; i<size; i++) nt += (nprocs[i] > 0);
829:   *nproc    = nt;
830:   PetscMalloc((nt+1)*sizeof(PetscInt),procs);
831:   PetscMalloc((nt+1)*sizeof(PetscInt),numprocs);
832:   PetscMalloc((nt+1)*sizeof(PetscInt*),indices);
833:   PetscMalloc(size*sizeof(PetscInt),&bprocs);
834:   cnt       = 0;
835:   for (i=0; i<size; i++) {
836:     if (nprocs[i] > 0) {
837:       bprocs[i]        = cnt;
838:       (*procs)[cnt]    = i;
839:       (*numprocs)[cnt] = nprocs[i];
840:       PetscMalloc(nprocs[i]*sizeof(PetscInt),&(*indices)[cnt]);
841:       cnt++;
842:     }
843:   }

845:   /* make the list of subdomains for each nontrivial local node */
846:   PetscMemzero(*numprocs,nt*sizeof(PetscInt));
847:   cnt  = 0;
848:   for (i=0; i<nrecvs2; i++) {
849:     nt = recvs2[cnt++];
850:     for (j=0; j<nt; j++) {
851:       for (k=0; k<recvs2[cnt+1]; k++) {
852:         (*indices)[bprocs[recvs2[cnt+2+k]]][(*numprocs)[bprocs[recvs2[cnt+2+k]]]++] = recvs2[cnt];
853:       }
854:       cnt += 2 + recvs2[cnt+1];
855:     }
856:   }
857:   PetscFree(bprocs);
858:   PetscFree(recvs2);

860:   /* sort the node indexing by their global numbers */
861:   nt = *nproc;
862:   for (i=0; i<nt; i++) {
863:     PetscMalloc(((*numprocs)[i])*sizeof(PetscInt),&tmp);
864:     for (j=0; j<(*numprocs)[i]; j++) {
865:       tmp[j] = lindices[(*indices)[i][j]];
866:     }
867:     PetscSortIntWithArray((*numprocs)[i],tmp,(*indices)[i]);
868:     PetscFree(tmp);
869:   }

871:   if (debug) { /* -----------------------------------  */
872:     nt = *nproc;
873:     for (i=0; i<nt; i++) {
874:       PetscSynchronizedPrintf(comm,"[%d] subdomain %d number of indices %d: ",rank,(*procs)[i],(*numprocs)[i]);
875:       for (j=0; j<(*numprocs)[i]; j++) {
876:         PetscSynchronizedPrintf(comm,"%d ",(*indices)[i][j]);
877:       }
878:       PetscSynchronizedPrintf(comm,"\n");
879:     }
880:     PetscSynchronizedFlush(comm);
881:   } /* -----------------------------------  */

883:   /* wait on sends */
884:   if (nsends2) {
885:     PetscMalloc(nsends2*sizeof(MPI_Status),&send_status);
886:     MPI_Waitall(nsends2,send_waits,send_status);
887:     PetscFree(send_status);
888:   }

890:   PetscFree(starts3);
891:   PetscFree(dest);
892:   PetscFree(send_waits);

894:   PetscFree(nownedsenders);
895:   PetscFree(ownedsenders);
896:   PetscFree(starts);
897:   PetscFree(starts2);
898:   PetscFree(lens2);

900:   PetscFree(source);
901:   PetscFree(len);
902:   PetscFree(recvs);
903:   PetscFree(nprocs);
904:   PetscFree(sends2);

906:   /* put the information about myself as the first entry in the list */
907:   first_procs    = (*procs)[0];
908:   first_numprocs = (*numprocs)[0];
909:   first_indices  = (*indices)[0];
910:   for (i=0; i<*nproc; i++) {
911:     if ((*procs)[i] == rank) {
912:       (*procs)[0]    = (*procs)[i];
913:       (*numprocs)[0] = (*numprocs)[i];
914:       (*indices)[0]  = (*indices)[i];
915:       (*procs)[i]    = first_procs;
916:       (*numprocs)[i] = first_numprocs;
917:       (*indices)[i]  = first_indices;
918:       break;
919:     }
920:   }
921:   return(0);
922: }

926: /*@C
927:     ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by ISLocalToGlobalMappingGetInfo()

929:     Collective on ISLocalToGlobalMapping

931:     Input Parameters:
932: .   mapping - the mapping from local to global indexing

934:     Output Parameter:
935: +   nproc - number of processors that are connected to this one
936: .   proc - neighboring processors
937: .   numproc - number of indices for each processor
938: -   indices - indices of local nodes shared with neighbor (sorted by global numbering)

940:     Level: advanced

942: .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
943:           ISLocalToGlobalMappingGetInfo()
944: @*/
945: PetscErrorCode  ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
946: {
948:   PetscInt       i;

951:   PetscFree(*procs);
952:   PetscFree(*numprocs);
953:   if (*indices) {
954:     PetscFree((*indices)[0]);
955:     for (i=1; i<*nproc; i++) {
956:       PetscFree((*indices)[i]);
957:     }
958:     PetscFree(*indices);
959:   }
960:   return(0);
961: }