Actual source code: isltog.c

petsc-3.5.1 2014-08-06
Report Typos and Errors
  2: #include <petsc-private/isimpl.h>    /*I "petscis.h"  I*/
  3: #include <petscsf.h>
  4: #include <petscviewer.h>

  6: PetscClassId IS_LTOGM_CLASSID;


 11: PetscErrorCode ISG2LMapApply(ISLocalToGlobalMapping mapping,PetscInt n,const PetscInt in[],PetscInt out[])
 12: {
 14:   PetscInt       i,start,end;

 17:   if (!mapping->globals) {
 18:     ISGlobalToLocalMappingApply(mapping,IS_GTOLM_MASK,0,0,0,0);
 19:   }
 20:   start = mapping->globalstart;
 21:   end = mapping->globalend;
 22:   for (i=0; i<n; i++) {
 23:     if (in[i] < 0)          out[i] = in[i];
 24:     else if (in[i] < start) out[i] = -1;
 25:     else if (in[i] > end)   out[i] = -1;
 26:     else                    out[i] = mapping->globals[in[i] - start];
 27:   }
 28:   return(0);
 29: }


 34: /*@
 35:     ISLocalToGlobalMappingGetSize - Gets the local size of a local to global mapping

 37:     Not Collective

 39:     Input Parameter:
 40: .   ltog - local to global mapping

 42:     Output Parameter:
 43: .   n - the number of entries in the local mapping, ISLocalToGlobalMappingGetIndices() returns an array of this length

 45:     Level: advanced

 47:     Concepts: mapping^local to global

 49: .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
 50: @*/
 51: PetscErrorCode  ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping mapping,PetscInt *n)
 52: {
 56:   *n = mapping->bs*mapping->n;
 57:   return(0);
 58: }

 62: /*@C
 63:     ISLocalToGlobalMappingView - View a local to global mapping

 65:     Not Collective

 67:     Input Parameters:
 68: +   ltog - local to global mapping
 69: -   viewer - viewer

 71:     Level: advanced

 73:     Concepts: mapping^local to global

 75: .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
 76: @*/
 77: PetscErrorCode  ISLocalToGlobalMappingView(ISLocalToGlobalMapping mapping,PetscViewer viewer)
 78: {
 79:   PetscInt       i;
 80:   PetscMPIInt    rank;
 81:   PetscBool      iascii;

 86:   if (!viewer) {
 87:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mapping),&viewer);
 88:   }

 91:   MPI_Comm_rank(PetscObjectComm((PetscObject)mapping),&rank);
 92:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
 93:   if (iascii) {
 94:     PetscObjectPrintClassNamePrefixType((PetscObject)mapping,viewer);
 95:     PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);
 96:     for (i=0; i<mapping->n; i++) {
 97:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D %D\n",rank,i,mapping->indices[i]);
 98:     }
 99:     PetscViewerFlush(viewer);
100:     PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);
101:   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for ISLocalToGlobalMapping",((PetscObject)viewer)->type_name);
102:   return(0);
103: }

107: /*@
108:     ISLocalToGlobalMappingCreateIS - Creates a mapping between a local (0 to n)
109:     ordering and a global parallel ordering.

111:     Not collective

113:     Input Parameter:
114: .   is - index set containing the global numbers for each local number

116:     Output Parameter:
117: .   mapping - new mapping data structure

119:     Notes: the block size of the IS determines the block size of the mapping
120:     Level: advanced

122:     Concepts: mapping^local to global

124: .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
125: @*/
126: PetscErrorCode  ISLocalToGlobalMappingCreateIS(IS is,ISLocalToGlobalMapping *mapping)
127: {
129:   PetscInt       n,bs;
130:   const PetscInt *indices;
131:   MPI_Comm       comm;
132:   PetscBool      isblock;


138:   PetscObjectGetComm((PetscObject)is,&comm);
139:   ISGetLocalSize(is,&n);
140:   PetscObjectTypeCompare((PetscObject)is,ISBLOCK,&isblock);
141:   ISGetBlockSize(is,&bs);
142:   /*  if (!isblock) { */
143:     ISGetIndices(is,&indices);
144:     ISLocalToGlobalMappingCreate(comm,1,n,indices,PETSC_COPY_VALUES,mapping);
145:     ISRestoreIndices(is,&indices);
146:     /*  } else {
147:     ISBlockGetIndices(is,&indices);
148:     ISLocalToGlobalMappingCreate(comm,bs,n,indices,PETSC_COPY_VALUES,mapping);
149:     ISBlockRestoreIndices(is,&indices);
150:      }*/
151:   return(0);
152: }

156: /*@C
157:     ISLocalToGlobalMappingCreateSF - Creates a mapping between a local (0 to n)
158:     ordering and a global parallel ordering.

160:     Collective

162:     Input Parameter:
163: +   sf - star forest mapping contiguous local indices to (rank, offset)
164: -   start - first global index on this process

166:     Output Parameter:
167: .   mapping - new mapping data structure

169:     Level: advanced

171:     Concepts: mapping^local to global

173: .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingCreateIS()
174: @*/
175: PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF sf,PetscInt start,ISLocalToGlobalMapping *mapping)
176: {
178:   PetscInt       i,maxlocal,nroots,nleaves,*globals,*ltog;
179:   const PetscInt *ilocal;
180:   MPI_Comm       comm;


186:   PetscObjectGetComm((PetscObject)sf,&comm);
187:   PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);
188:   if (ilocal) {
189:     for (i=0,maxlocal=0; i<nleaves; i++) maxlocal = PetscMax(maxlocal,ilocal[i]+1);
190:   }
191:   else maxlocal = nleaves;
192:   PetscMalloc1(nroots,&globals);
193:   PetscMalloc1(maxlocal,&ltog);
194:   for (i=0; i<nroots; i++) globals[i] = start + i;
195:   for (i=0; i<maxlocal; i++) ltog[i] = -1;
196:   PetscSFBcastBegin(sf,MPIU_INT,globals,ltog);
197:   PetscSFBcastEnd(sf,MPIU_INT,globals,ltog);
198:   ISLocalToGlobalMappingCreate(comm,1,maxlocal,ltog,PETSC_OWN_POINTER,mapping);
199:   PetscFree(globals);
200:   return(0);
201: }

205: /*@
206:     ISLocalToGlobalMappingGetBlockSize - Gets the blocksize of the mapping
207:     ordering and a global parallel ordering.

209:     Not Collective

211:     Input Parameters:
212: .   mapping - mapping data structure

214:     Output Parameter:
215: .   bs - the blocksize

217:     Level: advanced

219:     Concepts: mapping^local to global

221: .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS()
222: @*/
223: PetscErrorCode  ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping mapping,PetscInt *bs)
224: {
226:   *bs = mapping->bs;
227:   return(0);
228: }

232: /*@
233:     ISLocalToGlobalMappingCreate - Creates a mapping between a local (0 to n)
234:     ordering and a global parallel ordering.

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

238:     Input Parameters:
239: +   comm - MPI communicator
240: .   bs - the block size
241: .   n - the number of local elements
242: .   indices - the global index for each local element, these do not need to be in increasing order (sorted), these values should not be scaled by the blocksize bs
243: -   mode - see PetscCopyMode

245:     Output Parameter:
246: .   mapping - new mapping data structure

248:     Notes: There is one integer value in indices per block and it represents the actual indices bs*idx + j, where j=0,..,bs-1
249:     Level: advanced

251:     Concepts: mapping^local to global

253: .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS()
254: @*/
255: PetscErrorCode  ISLocalToGlobalMappingCreate(MPI_Comm cm,PetscInt bs,PetscInt n,const PetscInt indices[],PetscCopyMode mode,ISLocalToGlobalMapping *mapping)
256: {
258:   PetscInt       *in;


264:   *mapping = NULL;
265:   ISInitializePackage();

267:   PetscHeaderCreate(*mapping,_p_ISLocalToGlobalMapping,int,IS_LTOGM_CLASSID,"ISLocalToGlobalMapping","Local to global mapping","IS",
268:                            cm,ISLocalToGlobalMappingDestroy,ISLocalToGlobalMappingView);
269:   (*mapping)->n  = n;
270:   (*mapping)->bs = bs;
271:   /*
272:     Do not create the global to local mapping. This is only created if
273:     ISGlobalToLocalMapping() is called
274:   */
275:   (*mapping)->globals = 0;
276:   if (mode == PETSC_COPY_VALUES) {
277:     PetscMalloc1(n,&in);
278:     PetscMemcpy(in,indices,n*sizeof(PetscInt));
279:     (*mapping)->indices = in;
280:     PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt));
281:   } else if (mode == PETSC_OWN_POINTER) {
282:     (*mapping)->indices = (PetscInt*)indices;
283:     PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt));
284:   }
285:   else SETERRQ(cm,PETSC_ERR_SUP,"Cannot currently use PETSC_USE_POINTER");
286:   return(0);
287: }

291: /*@
292:    ISLocalToGlobalMappingDestroy - Destroys a mapping between a local (0 to n)
293:    ordering and a global parallel ordering.

295:    Note Collective

297:    Input Parameters:
298: .  mapping - mapping data structure

300:    Level: advanced

302: .seealso: ISLocalToGlobalMappingCreate()
303: @*/
304: PetscErrorCode  ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping *mapping)
305: {

309:   if (!*mapping) return(0);
311:   if (--((PetscObject)(*mapping))->refct > 0) {*mapping = 0;return(0);}
312:   PetscFree((*mapping)->indices);
313:   PetscFree((*mapping)->globals);
314:   PetscHeaderDestroy(mapping);
315:   *mapping = 0;
316:   return(0);
317: }

321: /*@
322:     ISLocalToGlobalMappingApplyIS - Creates from an IS in the local numbering
323:     a new index set using the global numbering defined in an ISLocalToGlobalMapping
324:     context.

326:     Not collective

328:     Input Parameters:
329: +   mapping - mapping between local and global numbering
330: -   is - index set in local numbering

332:     Output Parameters:
333: .   newis - index set in global numbering

335:     Level: advanced

337:     Concepts: mapping^local to global

339: .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),
340:           ISLocalToGlobalMappingDestroy(), ISGlobalToLocalMappingApply()
341: @*/
342: PetscErrorCode  ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping,IS is,IS *newis)
343: {
345:   PetscInt       n,*idxout;
346:   const PetscInt *idxin;


353:   ISGetLocalSize(is,&n);
354:   ISGetIndices(is,&idxin);
355:   PetscMalloc1(n,&idxout);
356:   ISLocalToGlobalMappingApply(mapping,n,idxin,idxout);
357:   ISRestoreIndices(is,&idxin);
358:   ISCreateGeneral(PetscObjectComm((PetscObject)is),n,idxout,PETSC_OWN_POINTER,newis);
359:   return(0);
360: }

364: /*@
365:    ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering
366:    and converts them to the global numbering.

368:    Not collective

370:    Input Parameters:
371: +  mapping - the local to global mapping context
372: .  N - number of integers
373: -  in - input indices in local numbering

375:    Output Parameter:
376: .  out - indices in global numbering

378:    Notes:
379:    The in and out array parameters may be identical.

381:    Level: advanced

383: .seealso: ISLocalToGlobalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
384:           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
385:           AOPetscToApplication(), ISGlobalToLocalMappingApply()

387:     Concepts: mapping^local to global
388: @*/
389: PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
390: {
391:   PetscInt       i,bs = mapping->bs,Nmax = bs*mapping->n;
392:   const PetscInt *idx = mapping->indices;

395:   if (bs == 1) {
396:     for (i=0; i<N; i++) {
397:       if (in[i] < 0) {
398:         out[i] = in[i];
399:         continue;
400:       }
401:       if (in[i] >= Nmax) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local index %D too large %D (max) at %D",in[i],Nmax-1,i);
402:       out[i] = idx[in[i]];
403:     }
404:   } else {
405:     for (i=0; i<N; i++) {
406:       if (in[i] < 0) {
407:         out[i] = in[i];
408:         continue;
409:       }
410:       if (in[i] >= Nmax) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local index %D too large %D (max) at %D",in[i],Nmax-1,i);
411:       out[i] = idx[in[i]/bs]*bs + (in[i] % bs);
412:     }
413:   }
414:   return(0);
415: }

419: /*@
420:    ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local numbering  and converts them to the global numbering.

422:    Not collective

424:    Input Parameters:
425: +  mapping - the local to global mapping context
426: .  N - number of integers
427: -  in - input indices in local numbering

429:    Output Parameter:
430: .  out - indices in global numbering

432:    Notes:
433:    The in and out array parameters may be identical.

435:    Level: advanced

437: .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
438:           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
439:           AOPetscToApplication(), ISGlobalToLocalMappingApply()

441:     Concepts: mapping^local to global
442: @*/
443: PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
444: {
445:   PetscInt       i,Nmax = mapping->n;
446:   const PetscInt *idx = mapping->indices;

449:   for (i=0; i<N; i++) {
450:     if (in[i] < 0) {
451:       out[i] = in[i];
452:       continue;
453:     }
454:     if (in[i] >= Nmax) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local block index %D too large %D (max) at %D",in[i],Nmax-1,i);
455:     out[i] = idx[in[i]];
456:   }
457:   return(0);
458: }

460: /* -----------------------------------------------------------------------------------------*/

464: /*
465:     Creates the global fields in the ISLocalToGlobalMapping structure
466: */
467: static PetscErrorCode ISGlobalToLocalMappingSetUp_Private(ISLocalToGlobalMapping mapping)
468: {
470:   PetscInt       i,*idx = mapping->indices,n = mapping->n,end,start,*globals;

473:   end   = 0;
474:   start = PETSC_MAX_INT;

476:   for (i=0; i<n; i++) {
477:     if (idx[i] < 0) continue;
478:     if (idx[i] < start) start = idx[i];
479:     if (idx[i] > end)   end   = idx[i];
480:   }
481:   if (start > end) {start = 0; end = -1;}
482:   mapping->globalstart = start;
483:   mapping->globalend   = end;

485:   PetscMalloc1((end-start+2),&globals);
486:   mapping->globals = globals;
487:   for (i=0; i<end-start+1; i++) globals[i] = -1;
488:   for (i=0; i<n; i++) {
489:     if (idx[i] < 0) continue;
490:     globals[idx[i] - start] = i;
491:   }

493:   PetscLogObjectMemory((PetscObject)mapping,(end-start+1)*sizeof(PetscInt));
494:   return(0);
495: }

499: /*@
500:     ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers
501:     specified with a global numbering.

503:     Not collective

505:     Input Parameters:
506: +   mapping - mapping between local and global numbering
507: .   type - IS_GTOLM_MASK - replaces global indices with no local value with -1
508:            IS_GTOLM_DROP - drops the indices with no local value from the output list
509: .   n - number of global indices to map
510: -   idx - global indices to map

512:     Output Parameters:
513: +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
514: -   idxout - local index of each global index, one must pass in an array long enough
515:              to hold all the indices. You can call ISGlobalToLocalMappingApply() with
516:              idxout == NULL to determine the required length (returned in nout)
517:              and then allocate the required space and call ISGlobalToLocalMappingApply()
518:              a second time to set the values.

520:     Notes:
521:     Either nout or idxout may be NULL. idx and idxout may be identical.

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

526:     Level: advanced

528:     Developer Note: The manual page states that idx and idxout may be identical but the calling
529:        sequence declares idx as const so it cannot be the same as idxout.

531:     Concepts: mapping^global to local

533: .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),
534:           ISLocalToGlobalMappingDestroy()
535: @*/
536: PetscErrorCode  ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingType type,
537:                                             PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
538: {
539:   PetscInt       i,*globals,nf = 0,tmp,start,end,bs;

544:   if (!mapping->globals) {
545:     ISGlobalToLocalMappingSetUp_Private(mapping);
546:   }
547:   globals = mapping->globals;
548:   start   = mapping->globalstart;
549:   end     = mapping->globalend;
550:   bs      = mapping->bs;

552:   if (type == IS_GTOLM_MASK) {
553:     if (idxout) {
554:       for (i=0; i<n; i++) {
555:         if (idx[i] < 0) idxout[i] = idx[i];
556:         else if (idx[i] < bs*start) idxout[i] = -1;
557:         else if (idx[i] > bs*end)   idxout[i] = -1;
558:         else                        idxout[i] = bs*globals[idx[i]/bs - start] + (idx[i] % bs);
559:       }
560:     }
561:     if (nout) *nout = n;
562:   } else {
563:     if (idxout) {
564:       for (i=0; i<n; i++) {
565:         if (idx[i] < 0) continue;
566:         if (idx[i] < bs*start) continue;
567:         if (idx[i] > bs*end) continue;
568:         tmp = bs*globals[idx[i]/bs - start] + (idx[i] % bs);
569:         if (tmp < 0) continue;
570:         idxout[nf++] = tmp;
571:       }
572:     } else {
573:       for (i=0; i<n; i++) {
574:         if (idx[i] < 0) continue;
575:         if (idx[i] < bs*start) continue;
576:         if (idx[i] > bs*end) continue;
577:         tmp = bs*globals[idx[i]/bs - start] + (idx[i] % bs);
578:         if (tmp < 0) continue;
579:         nf++;
580:       }
581:     }
582:     if (nout) *nout = nf;
583:   }
584:   return(0);
585: }

589: /*@
590:     ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers
591:     specified with a block global numbering.

593:     Not collective

595:     Input Parameters:
596: +   mapping - mapping between local and global numbering
597: .   type - IS_GTOLM_MASK - replaces global indices with no local value with -1
598:            IS_GTOLM_DROP - drops the indices with no local value from the output list
599: .   n - number of global indices to map
600: -   idx - global indices to map

602:     Output Parameters:
603: +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
604: -   idxout - local index of each global index, one must pass in an array long enough
605:              to hold all the indices. You can call ISGlobalToLocalMappingApplyBlock() with
606:              idxout == NULL to determine the required length (returned in nout)
607:              and then allocate the required space and call ISGlobalToLocalMappingApplyBlock()
608:              a second time to set the values.

610:     Notes:
611:     Either nout or idxout may be NULL. idx and idxout may be identical.

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

616:     Level: advanced

618:     Developer Note: The manual page states that idx and idxout may be identical but the calling
619:        sequence declares idx as const so it cannot be the same as idxout.

621:     Concepts: mapping^global to local

623: .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(),
624:           ISLocalToGlobalMappingDestroy()
625: @*/
626: PetscErrorCode  ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingType type,
627:                                   PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
628: {
629:   PetscInt       i,*globals,nf = 0,tmp,start,end;

634:   if (!mapping->globals) {
635:     ISGlobalToLocalMappingSetUp_Private(mapping);
636:   }
637:   globals = mapping->globals;
638:   start   = mapping->globalstart;
639:   end     = mapping->globalend;

641:   if (type == IS_GTOLM_MASK) {
642:     if (idxout) {
643:       for (i=0; i<n; i++) {
644:         if (idx[i] < 0) idxout[i] = idx[i];
645:         else if (idx[i] < start) idxout[i] = -1;
646:         else if (idx[i] > end)   idxout[i] = -1;
647:         else                     idxout[i] = globals[idx[i] - start];
648:       }
649:     }
650:     if (nout) *nout = n;
651:   } else {
652:     if (idxout) {
653:       for (i=0; i<n; i++) {
654:         if (idx[i] < 0) continue;
655:         if (idx[i] < start) continue;
656:         if (idx[i] > end) continue;
657:         tmp = globals[idx[i] - start];
658:         if (tmp < 0) continue;
659:         idxout[nf++] = tmp;
660:       }
661:     } else {
662:       for (i=0; i<n; i++) {
663:         if (idx[i] < 0) continue;
664:         if (idx[i] < start) continue;
665:         if (idx[i] > end) continue;
666:         tmp = globals[idx[i] - start];
667:         if (tmp < 0) continue;
668:         nf++;
669:       }
670:     }
671:     if (nout) *nout = nf;
672:   }
673:   return(0);
674: }

678: /*@C
679:     ISLocalToGlobalMappingGetBlockInfo - Gets the neighbor information for each processor and
680:      each index shared by more than one processor

682:     Collective on ISLocalToGlobalMapping

684:     Input Parameters:
685: .   mapping - the mapping from local to global indexing

687:     Output Parameter:
688: +   nproc - number of processors that are connected to this one
689: .   proc - neighboring processors
690: .   numproc - number of indices for each subdomain (processor)
691: -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)

693:     Level: advanced

695:     Concepts: mapping^local to global

697:     Fortran Usage:
698: $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
699: $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
700:           PetscInt indices[nproc][numprocmax],ierr)
701:         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
702:         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.


705: .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
706:           ISLocalToGlobalMappingRestoreInfo()
707: @*/
708: PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
709: {
711:   PetscMPIInt    size,rank,tag1,tag2,tag3,*len,*source,imdex;
712:   PetscInt       i,n = mapping->n,Ng,ng,max = 0,*lindices = mapping->indices;
713:   PetscInt       *nprocs,*owner,nsends,*sends,j,*starts,nmax,nrecvs,*recvs,proc;
714:   PetscInt       cnt,scale,*ownedsenders,*nownedsenders,rstart,nowned;
715:   PetscInt       node,nownedm,nt,*sends2,nsends2,*starts2,*lens2,*dest,nrecvs2,*starts3,*recvs2,k,*bprocs,*tmp;
716:   PetscInt       first_procs,first_numprocs,*first_indices;
717:   MPI_Request    *recv_waits,*send_waits;
718:   MPI_Status     recv_status,*send_status,*recv_statuses;
719:   MPI_Comm       comm;
720:   PetscBool      debug = PETSC_FALSE;

724:   PetscObjectGetComm((PetscObject)mapping,&comm);
725:   MPI_Comm_size(comm,&size);
726:   MPI_Comm_rank(comm,&rank);
727:   if (size == 1) {
728:     *nproc         = 0;
729:     *procs         = NULL;
730:     PetscMalloc(sizeof(PetscInt),numprocs);
731:     (*numprocs)[0] = 0;
732:     PetscMalloc(sizeof(PetscInt*),indices);
733:     (*indices)[0]  = NULL;
734:     return(0);
735:   }

737:   PetscOptionsGetBool(NULL,"-islocaltoglobalmappinggetinfo_debug",&debug,NULL);

739:   /*
740:     Notes on ISLocalToGlobalMappingGetBlockInfo

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

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

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

750:     nontrivial locally owned node - node that is not in the interior (i.e. has more than one
751:            local subdomain
752:   */
753:   PetscObjectGetNewTag((PetscObject)mapping,&tag1);
754:   PetscObjectGetNewTag((PetscObject)mapping,&tag2);
755:   PetscObjectGetNewTag((PetscObject)mapping,&tag3);

757:   for (i=0; i<n; i++) {
758:     if (lindices[i] > max) max = lindices[i];
759:   }
760:   MPI_Allreduce(&max,&Ng,1,MPIU_INT,MPI_MAX,comm);
761:   Ng++;
762:   MPI_Comm_size(comm,&size);
763:   MPI_Comm_rank(comm,&rank);
764:   scale  = Ng/size + 1;
765:   ng     = scale; if (rank == size-1) ng = Ng - scale*(size-1); ng = PetscMax(1,ng);
766:   rstart = scale*rank;

768:   /* determine ownership ranges of global indices */
769:   PetscMalloc1(2*size,&nprocs);
770:   PetscMemzero(nprocs,2*size*sizeof(PetscInt));

772:   /* determine owners of each local node  */
773:   PetscMalloc1(n,&owner);
774:   for (i=0; i<n; i++) {
775:     proc             = lindices[i]/scale; /* processor that globally owns this index */
776:     nprocs[2*proc+1] = 1;                 /* processor globally owns at least one of ours */
777:     owner[i]         = proc;
778:     nprocs[2*proc]++;                     /* count of how many that processor globally owns of ours */
779:   }
780:   nsends = 0; for (i=0; i<size; i++) nsends += nprocs[2*i+1];
781:   PetscInfo1(mapping,"Number of global owners for my local data %D\n",nsends);

783:   /* inform other processors of number of messages and max length*/
784:   PetscMaxSum(comm,nprocs,&nmax,&nrecvs);
785:   PetscInfo1(mapping,"Number of local owners for my global data %D\n",nrecvs);

787:   /* post receives for owned rows */
788:   PetscMalloc1((2*nrecvs+1)*(nmax+1),&recvs);
789:   PetscMalloc1((nrecvs+1),&recv_waits);
790:   for (i=0; i<nrecvs; i++) {
791:     MPI_Irecv(recvs+2*nmax*i,2*nmax,MPIU_INT,MPI_ANY_SOURCE,tag1,comm,recv_waits+i);
792:   }

794:   /* pack messages containing lists of local nodes to owners */
795:   PetscMalloc1((2*n+1),&sends);
796:   PetscMalloc1((size+1),&starts);
797:   starts[0] = 0;
798:   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
799:   for (i=0; i<n; i++) {
800:     sends[starts[owner[i]]++] = lindices[i];
801:     sends[starts[owner[i]]++] = i;
802:   }
803:   PetscFree(owner);
804:   starts[0] = 0;
805:   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];

807:   /* send the messages */
808:   PetscMalloc1((nsends+1),&send_waits);
809:   PetscMalloc1((nsends+1),&dest);
810:   cnt = 0;
811:   for (i=0; i<size; i++) {
812:     if (nprocs[2*i]) {
813:       MPI_Isend(sends+starts[i],2*nprocs[2*i],MPIU_INT,i,tag1,comm,send_waits+cnt);
814:       dest[cnt] = i;
815:       cnt++;
816:     }
817:   }
818:   PetscFree(starts);

820:   /* wait on receives */
821:   PetscMalloc1((nrecvs+1),&source);
822:   PetscMalloc1((nrecvs+1),&len);
823:   cnt  = nrecvs;
824:   PetscMalloc1((ng+1),&nownedsenders);
825:   PetscMemzero(nownedsenders,ng*sizeof(PetscInt));
826:   while (cnt) {
827:     MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);
828:     /* unpack receives into our local space */
829:     MPI_Get_count(&recv_status,MPIU_INT,&len[imdex]);
830:     source[imdex] = recv_status.MPI_SOURCE;
831:     len[imdex]    = len[imdex]/2;
832:     /* count how many local owners for each of my global owned indices */
833:     for (i=0; i<len[imdex]; i++) nownedsenders[recvs[2*imdex*nmax+2*i]-rstart]++;
834:     cnt--;
835:   }
836:   PetscFree(recv_waits);

838:   /* count how many globally owned indices are on an edge multiplied by how many processors own them. */
839:   nowned  = 0;
840:   nownedm = 0;
841:   for (i=0; i<ng; i++) {
842:     if (nownedsenders[i] > 1) {nownedm += nownedsenders[i]; nowned++;}
843:   }

845:   /* create single array to contain rank of all local owners of each globally owned index */
846:   PetscMalloc1((nownedm+1),&ownedsenders);
847:   PetscMalloc1((ng+1),&starts);
848:   starts[0] = 0;
849:   for (i=1; i<ng; i++) {
850:     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
851:     else starts[i] = starts[i-1];
852:   }

854:   /* for each nontrival globally owned node list all arriving processors */
855:   for (i=0; i<nrecvs; i++) {
856:     for (j=0; j<len[i]; j++) {
857:       node = recvs[2*i*nmax+2*j]-rstart;
858:       if (nownedsenders[node] > 1) ownedsenders[starts[node]++] = source[i];
859:     }
860:   }

862:   if (debug) { /* -----------------------------------  */
863:     starts[0] = 0;
864:     for (i=1; i<ng; i++) {
865:       if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
866:       else starts[i] = starts[i-1];
867:     }
868:     for (i=0; i<ng; i++) {
869:       if (nownedsenders[i] > 1) {
870:         PetscSynchronizedPrintf(comm,"[%d] global node %D local owner processors: ",rank,i+rstart);
871:         for (j=0; j<nownedsenders[i]; j++) {
872:           PetscSynchronizedPrintf(comm,"%D ",ownedsenders[starts[i]+j]);
873:         }
874:         PetscSynchronizedPrintf(comm,"\n");
875:       }
876:     }
877:     PetscSynchronizedFlush(comm,PETSC_STDOUT);
878:   } /* -----------------------------------  */

880:   /* wait on original sends */
881:   if (nsends) {
882:     PetscMalloc1(nsends,&send_status);
883:     MPI_Waitall(nsends,send_waits,send_status);
884:     PetscFree(send_status);
885:   }
886:   PetscFree(send_waits);
887:   PetscFree(sends);
888:   PetscFree(nprocs);

890:   /* pack messages to send back to local owners */
891:   starts[0] = 0;
892:   for (i=1; i<ng; i++) {
893:     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
894:     else starts[i] = starts[i-1];
895:   }
896:   nsends2 = nrecvs;
897:   PetscMalloc1((nsends2+1),&nprocs); /* length of each message */
898:   for (i=0; i<nrecvs; i++) {
899:     nprocs[i] = 1;
900:     for (j=0; j<len[i]; j++) {
901:       node = recvs[2*i*nmax+2*j]-rstart;
902:       if (nownedsenders[node] > 1) nprocs[i] += 2 + nownedsenders[node];
903:     }
904:   }
905:   nt = 0;
906:   for (i=0; i<nsends2; i++) nt += nprocs[i];

908:   PetscMalloc1((nt+1),&sends2);
909:   PetscMalloc1((nsends2+1),&starts2);

911:   starts2[0] = 0;
912:   for (i=1; i<nsends2; i++) starts2[i] = starts2[i-1] + nprocs[i-1];
913:   /*
914:      Each message is 1 + nprocs[i] long, and consists of
915:        (0) the number of nodes being sent back
916:        (1) the local node number,
917:        (2) the number of processors sharing it,
918:        (3) the processors sharing it
919:   */
920:   for (i=0; i<nsends2; i++) {
921:     cnt = 1;
922:     sends2[starts2[i]] = 0;
923:     for (j=0; j<len[i]; j++) {
924:       node = recvs[2*i*nmax+2*j]-rstart;
925:       if (nownedsenders[node] > 1) {
926:         sends2[starts2[i]]++;
927:         sends2[starts2[i]+cnt++] = recvs[2*i*nmax+2*j+1];
928:         sends2[starts2[i]+cnt++] = nownedsenders[node];
929:         PetscMemcpy(&sends2[starts2[i]+cnt],&ownedsenders[starts[node]],nownedsenders[node]*sizeof(PetscInt));
930:         cnt += nownedsenders[node];
931:       }
932:     }
933:   }

935:   /* receive the message lengths */
936:   nrecvs2 = nsends;
937:   PetscMalloc1((nrecvs2+1),&lens2);
938:   PetscMalloc1((nrecvs2+1),&starts3);
939:   PetscMalloc1((nrecvs2+1),&recv_waits);
940:   for (i=0; i<nrecvs2; i++) {
941:     MPI_Irecv(&lens2[i],1,MPIU_INT,dest[i],tag2,comm,recv_waits+i);
942:   }

944:   /* send the message lengths */
945:   for (i=0; i<nsends2; i++) {
946:     MPI_Send(&nprocs[i],1,MPIU_INT,source[i],tag2,comm);
947:   }

949:   /* wait on receives of lens */
950:   if (nrecvs2) {
951:     PetscMalloc1(nrecvs2,&recv_statuses);
952:     MPI_Waitall(nrecvs2,recv_waits,recv_statuses);
953:     PetscFree(recv_statuses);
954:   }
955:   PetscFree(recv_waits);

957:   starts3[0] = 0;
958:   nt         = 0;
959:   for (i=0; i<nrecvs2-1; i++) {
960:     starts3[i+1] = starts3[i] + lens2[i];
961:     nt          += lens2[i];
962:   }
963:   if (nrecvs2) nt += lens2[nrecvs2-1];

965:   PetscMalloc1((nt+1),&recvs2);
966:   PetscMalloc1((nrecvs2+1),&recv_waits);
967:   for (i=0; i<nrecvs2; i++) {
968:     MPI_Irecv(recvs2+starts3[i],lens2[i],MPIU_INT,dest[i],tag3,comm,recv_waits+i);
969:   }

971:   /* send the messages */
972:   PetscMalloc1((nsends2+1),&send_waits);
973:   for (i=0; i<nsends2; i++) {
974:     MPI_Isend(sends2+starts2[i],nprocs[i],MPIU_INT,source[i],tag3,comm,send_waits+i);
975:   }

977:   /* wait on receives */
978:   if (nrecvs2) {
979:     PetscMalloc1(nrecvs2,&recv_statuses);
980:     MPI_Waitall(nrecvs2,recv_waits,recv_statuses);
981:     PetscFree(recv_statuses);
982:   }
983:   PetscFree(recv_waits);
984:   PetscFree(nprocs);

986:   if (debug) { /* -----------------------------------  */
987:     cnt = 0;
988:     for (i=0; i<nrecvs2; i++) {
989:       nt = recvs2[cnt++];
990:       for (j=0; j<nt; j++) {
991:         PetscSynchronizedPrintf(comm,"[%d] local node %D number of subdomains %D: ",rank,recvs2[cnt],recvs2[cnt+1]);
992:         for (k=0; k<recvs2[cnt+1]; k++) {
993:           PetscSynchronizedPrintf(comm,"%D ",recvs2[cnt+2+k]);
994:         }
995:         cnt += 2 + recvs2[cnt+1];
996:         PetscSynchronizedPrintf(comm,"\n");
997:       }
998:     }
999:     PetscSynchronizedFlush(comm,PETSC_STDOUT);
1000:   } /* -----------------------------------  */

1002:   /* count number subdomains for each local node */
1003:   PetscMalloc1(size,&nprocs);
1004:   PetscMemzero(nprocs,size*sizeof(PetscInt));
1005:   cnt  = 0;
1006:   for (i=0; i<nrecvs2; i++) {
1007:     nt = recvs2[cnt++];
1008:     for (j=0; j<nt; j++) {
1009:       for (k=0; k<recvs2[cnt+1]; k++) nprocs[recvs2[cnt+2+k]]++;
1010:       cnt += 2 + recvs2[cnt+1];
1011:     }
1012:   }
1013:   nt = 0; for (i=0; i<size; i++) nt += (nprocs[i] > 0);
1014:   *nproc    = nt;
1015:   PetscMalloc1((nt+1),procs);
1016:   PetscMalloc1((nt+1),numprocs);
1017:   PetscMalloc1((nt+1),indices);
1018:   for (i=0;i<nt+1;i++) (*indices)[i]=NULL;
1019:   PetscMalloc1(size,&bprocs);
1020:   cnt       = 0;
1021:   for (i=0; i<size; i++) {
1022:     if (nprocs[i] > 0) {
1023:       bprocs[i]        = cnt;
1024:       (*procs)[cnt]    = i;
1025:       (*numprocs)[cnt] = nprocs[i];
1026:       PetscMalloc1(nprocs[i],&(*indices)[cnt]);
1027:       cnt++;
1028:     }
1029:   }

1031:   /* make the list of subdomains for each nontrivial local node */
1032:   PetscMemzero(*numprocs,nt*sizeof(PetscInt));
1033:   cnt  = 0;
1034:   for (i=0; i<nrecvs2; i++) {
1035:     nt = recvs2[cnt++];
1036:     for (j=0; j<nt; j++) {
1037:       for (k=0; k<recvs2[cnt+1]; k++) (*indices)[bprocs[recvs2[cnt+2+k]]][(*numprocs)[bprocs[recvs2[cnt+2+k]]]++] = recvs2[cnt];
1038:       cnt += 2 + recvs2[cnt+1];
1039:     }
1040:   }
1041:   PetscFree(bprocs);
1042:   PetscFree(recvs2);

1044:   /* sort the node indexing by their global numbers */
1045:   nt = *nproc;
1046:   for (i=0; i<nt; i++) {
1047:     PetscMalloc1(((*numprocs)[i]),&tmp);
1048:     for (j=0; j<(*numprocs)[i]; j++) tmp[j] = lindices[(*indices)[i][j]];
1049:     PetscSortIntWithArray((*numprocs)[i],tmp,(*indices)[i]);
1050:     PetscFree(tmp);
1051:   }

1053:   if (debug) { /* -----------------------------------  */
1054:     nt = *nproc;
1055:     for (i=0; i<nt; i++) {
1056:       PetscSynchronizedPrintf(comm,"[%d] subdomain %D number of indices %D: ",rank,(*procs)[i],(*numprocs)[i]);
1057:       for (j=0; j<(*numprocs)[i]; j++) {
1058:         PetscSynchronizedPrintf(comm,"%D ",(*indices)[i][j]);
1059:       }
1060:       PetscSynchronizedPrintf(comm,"\n");
1061:     }
1062:     PetscSynchronizedFlush(comm,PETSC_STDOUT);
1063:   } /* -----------------------------------  */

1065:   /* wait on sends */
1066:   if (nsends2) {
1067:     PetscMalloc1(nsends2,&send_status);
1068:     MPI_Waitall(nsends2,send_waits,send_status);
1069:     PetscFree(send_status);
1070:   }

1072:   PetscFree(starts3);
1073:   PetscFree(dest);
1074:   PetscFree(send_waits);

1076:   PetscFree(nownedsenders);
1077:   PetscFree(ownedsenders);
1078:   PetscFree(starts);
1079:   PetscFree(starts2);
1080:   PetscFree(lens2);

1082:   PetscFree(source);
1083:   PetscFree(len);
1084:   PetscFree(recvs);
1085:   PetscFree(nprocs);
1086:   PetscFree(sends2);

1088:   /* put the information about myself as the first entry in the list */
1089:   first_procs    = (*procs)[0];
1090:   first_numprocs = (*numprocs)[0];
1091:   first_indices  = (*indices)[0];
1092:   for (i=0; i<*nproc; i++) {
1093:     if ((*procs)[i] == rank) {
1094:       (*procs)[0]    = (*procs)[i];
1095:       (*numprocs)[0] = (*numprocs)[i];
1096:       (*indices)[0]  = (*indices)[i];
1097:       (*procs)[i]    = first_procs;
1098:       (*numprocs)[i] = first_numprocs;
1099:       (*indices)[i]  = first_indices;
1100:       break;
1101:     }
1102:   }
1103:   return(0);
1104: }

1108: /*@C
1109:     ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by ISLocalToGlobalMappingGetBlockInfo()

1111:     Collective on ISLocalToGlobalMapping

1113:     Input Parameters:
1114: .   mapping - the mapping from local to global indexing

1116:     Output Parameter:
1117: +   nproc - number of processors that are connected to this one
1118: .   proc - neighboring processors
1119: .   numproc - number of indices for each processor
1120: -   indices - indices of local nodes shared with neighbor (sorted by global numbering)

1122:     Level: advanced

1124: .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
1125:           ISLocalToGlobalMappingGetInfo()
1126: @*/
1127: PetscErrorCode  ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
1128: {
1130:   PetscInt       i;

1133:   PetscFree(*procs);
1134:   PetscFree(*numprocs);
1135:   if (*indices) {
1136:     PetscFree((*indices)[0]);
1137:     for (i=1; i<*nproc; i++) {
1138:       PetscFree((*indices)[i]);
1139:     }
1140:     PetscFree(*indices);
1141:   }
1142:   return(0);
1143: }

1147: /*@C
1148:     ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each processor and
1149:      each index shared by more than one processor

1151:     Collective on ISLocalToGlobalMapping

1153:     Input Parameters:
1154: .   mapping - the mapping from local to global indexing

1156:     Output Parameter:
1157: +   nproc - number of processors that are connected to this one
1158: .   proc - neighboring processors
1159: .   numproc - number of indices for each subdomain (processor)
1160: -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)

1162:     Level: advanced

1164:     Concepts: mapping^local to global

1166:     Fortran Usage:
1167: $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
1168: $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
1169:           PetscInt indices[nproc][numprocmax],ierr)
1170:         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
1171:         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.


1174: .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
1175:           ISLocalToGlobalMappingRestoreInfo()
1176: @*/
1177: PetscErrorCode  ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
1178: {
1180:   PetscInt       **bindices = NULL,bs = mapping->bs,i,j,k;

1184:   ISLocalToGlobalMappingGetBlockInfo(mapping,nproc,procs,numprocs,&bindices);
1185:   PetscCalloc1(*nproc,&*indices);
1186:   for (i=0; i<*nproc; i++) {
1187:     PetscMalloc1(bs*(*numprocs)[i],&(*indices)[i]);
1188:     for (j=0; j<(*numprocs)[i]; j++) {
1189:       for (k=0; k<bs; k++) {
1190:         (*indices)[i][j*bs+k] = bs*bindices[i][j] + k;
1191:       }
1192:     }
1193:     (*numprocs)[i] *= bs;
1194:   }
1195:   if (bindices) {
1196:     PetscFree(bindices[0]);
1197:     for (i=1; i<*nproc; i++) {
1198:       PetscFree(bindices[i]);
1199:     }
1200:     PetscFree(bindices);
1201:   }
1202:   return(0);
1203: }

1207: /*@C
1208:     ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by ISLocalToGlobalMappingGetInfo()

1210:     Collective on ISLocalToGlobalMapping

1212:     Input Parameters:
1213: .   mapping - the mapping from local to global indexing

1215:     Output Parameter:
1216: +   nproc - number of processors that are connected to this one
1217: .   proc - neighboring processors
1218: .   numproc - number of indices for each processor
1219: -   indices - indices of local nodes shared with neighbor (sorted by global numbering)

1221:     Level: advanced

1223: .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
1224:           ISLocalToGlobalMappingGetInfo()
1225: @*/
1226: PetscErrorCode  ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
1227: {

1231:   ISLocalToGlobalMappingRestoreBlockInfo(mapping,nproc,procs,numprocs,indices);
1232:   return(0);
1233: }

1237: /*@C
1238:    ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped

1240:    Not Collective

1242:    Input Arguments:
1243: . ltog - local to global mapping

1245:    Output Arguments:
1246: . array - array of indices, the length of this array may be obtained with ISLocalToGlobalMappingGetSize()

1248:    Level: advanced

1250:    Notes: ISLocalToGlobalMappingGetSize() returns the length the this array

1252: .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreIndices(), ISLocalToGlobalMappingGetBlockIndices(), ISLocalToGlobalMappingRestoreBlockIndices()
1253: @*/
1254: PetscErrorCode  ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
1255: {
1259:   if (ltog->bs == 1) {
1260:     *array = ltog->indices;
1261:   } else {
1262:     PetscInt       *jj,k,i,j,n = ltog->n, bs = ltog->bs;
1263:     const PetscInt *ii;

1266:     PetscMalloc1(bs*n,&jj);
1267:     *array = jj;
1268:     k    = 0;
1269:     ii   = ltog->indices;
1270:     for (i=0; i<n; i++)
1271:       for (j=0; j<bs; j++)
1272:         jj[k++] = bs*ii[i] + j;
1273:   }
1274:   return(0);
1275: }

1279: /*@C
1280:    ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with ISLocalToGlobalMappingRestoreIndices()

1282:    Not Collective

1284:    Input Arguments:
1285: + ltog - local to global mapping
1286: - array - array of indices

1288:    Level: advanced

1290: .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices()
1291: @*/
1292: PetscErrorCode  ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
1293: {
1297:   if (ltog->bs == 1 && *array != ltog->indices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");

1299:   if (ltog->bs > 1) {
1301:     PetscFree(*(void**)array);
1302:   }
1303:   return(0);
1304: }

1308: /*@C
1309:    ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block

1311:    Not Collective

1313:    Input Arguments:
1314: . ltog - local to global mapping

1316:    Output Arguments:
1317: . array - array of indices

1319:    Level: advanced

1321: .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreBlockIndices()
1322: @*/
1323: PetscErrorCode  ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
1324: {
1328:   *array = ltog->indices;
1329:   return(0);
1330: }

1334: /*@C
1335:    ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with ISLocalToGlobalMappingGetBlockIndices()

1337:    Not Collective

1339:    Input Arguments:
1340: + ltog - local to global mapping
1341: - array - array of indices

1343:    Level: advanced

1345: .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices()
1346: @*/
1347: PetscErrorCode  ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
1348: {
1352:   if (*array != ltog->indices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
1353:   *array = NULL;
1354:   return(0);
1355: }

1359: /*@C
1360:    ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings

1362:    Not Collective

1364:    Input Arguments:
1365: + comm - communicator for the new mapping, must contain the communicator of every mapping to concatenate
1366: . n - number of mappings to concatenate
1367: - ltogs - local to global mappings

1369:    Output Arguments:
1370: . ltogcat - new mapping

1372:    Note: this currently always returns a mapping with block size of 1

1374:    Developer Note: If all the input mapping have the same block size we could easily handle that as a special case

1376:    Level: advanced

1378: .seealso: ISLocalToGlobalMappingCreate()
1379: @*/
1380: PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm,PetscInt n,const ISLocalToGlobalMapping ltogs[],ISLocalToGlobalMapping *ltogcat)
1381: {
1382:   PetscInt       i,cnt,m,*idx;

1386:   if (n < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have a non-negative number of mappings, given %D",n);
1390:   for (cnt=0,i=0; i<n; i++) {
1391:     ISLocalToGlobalMappingGetSize(ltogs[i],&m);
1392:     cnt += m;
1393:   }
1394:   PetscMalloc1(cnt,&idx);
1395:   for (cnt=0,i=0; i<n; i++) {
1396:     const PetscInt *subidx;
1397:     ISLocalToGlobalMappingGetSize(ltogs[i],&m);
1398:     ISLocalToGlobalMappingGetIndices(ltogs[i],&subidx);
1399:     PetscMemcpy(&idx[cnt],subidx,m*sizeof(PetscInt));
1400:     ISLocalToGlobalMappingRestoreIndices(ltogs[i],&subidx);
1401:     cnt += m;
1402:   }
1403:   ISLocalToGlobalMappingCreate(comm,1,cnt,idx,PETSC_OWN_POINTER,ltogcat);
1404:   return(0);
1405: }