xref: /petsc/src/vec/is/utils/isltog.c (revision 7e99dc12794bd48581d0b21da67a7e4cbbce6ea1)
12362add9SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/isimpl.h>    /*I "petscis.h"  I*/
30c312b8eSJed Brown #include <petscsf.h>
4665c2dedSJed Brown #include <petscviewer.h>
52362add9SBarry Smith 
67087cfbeSBarry Smith PetscClassId IS_LTOGM_CLASSID;
7268a049cSStefano Zampini static PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping,PetscInt*,PetscInt**,PetscInt**,PetscInt***);
88e58c17dSMatthew Knepley 
96658fb44Sstefano_zampini /*@
106658fb44Sstefano_zampini     ISLocalToGlobalMappingDuplicate - Duplicates the local to global mapping object
116658fb44Sstefano_zampini 
126658fb44Sstefano_zampini     Not Collective
136658fb44Sstefano_zampini 
146658fb44Sstefano_zampini     Input Parameter:
156658fb44Sstefano_zampini .   ltog - local to global mapping
166658fb44Sstefano_zampini 
176658fb44Sstefano_zampini     Output Parameter:
186658fb44Sstefano_zampini .   nltog - the duplicated local to global mapping
196658fb44Sstefano_zampini 
206658fb44Sstefano_zampini     Level: advanced
216658fb44Sstefano_zampini 
226658fb44Sstefano_zampini     Concepts: mapping^local to global
236658fb44Sstefano_zampini 
246658fb44Sstefano_zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
256658fb44Sstefano_zampini @*/
266658fb44Sstefano_zampini PetscErrorCode  ISLocalToGlobalMappingDuplicate(ISLocalToGlobalMapping ltog,ISLocalToGlobalMapping* nltog)
276658fb44Sstefano_zampini {
286658fb44Sstefano_zampini   PetscErrorCode ierr;
296658fb44Sstefano_zampini 
306658fb44Sstefano_zampini   PetscFunctionBegin;
316658fb44Sstefano_zampini   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
326658fb44Sstefano_zampini   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)ltog),ltog->bs,ltog->n,ltog->indices,PETSC_COPY_VALUES,nltog);CHKERRQ(ierr);
336658fb44Sstefano_zampini   PetscFunctionReturn(0);
346658fb44Sstefano_zampini }
356658fb44Sstefano_zampini 
36565245c5SBarry Smith /*@
37107e9a97SBarry Smith     ISLocalToGlobalMappingGetSize - Gets the local size of a local to global mapping
383b9aefa3SBarry Smith 
393b9aefa3SBarry Smith     Not Collective
403b9aefa3SBarry Smith 
413b9aefa3SBarry Smith     Input Parameter:
423b9aefa3SBarry Smith .   ltog - local to global mapping
433b9aefa3SBarry Smith 
443b9aefa3SBarry Smith     Output Parameter:
45107e9a97SBarry Smith .   n - the number of entries in the local mapping, ISLocalToGlobalMappingGetIndices() returns an array of this length
463b9aefa3SBarry Smith 
473b9aefa3SBarry Smith     Level: advanced
483b9aefa3SBarry Smith 
49273d9f13SBarry Smith     Concepts: mapping^local to global
503b9aefa3SBarry Smith 
513b9aefa3SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
523b9aefa3SBarry Smith @*/
537087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping mapping,PetscInt *n)
543b9aefa3SBarry Smith {
553b9aefa3SBarry Smith   PetscFunctionBegin;
560700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
574482741eSBarry Smith   PetscValidIntPointer(n,2);
58107e9a97SBarry Smith   *n = mapping->bs*mapping->n;
593b9aefa3SBarry Smith   PetscFunctionReturn(0);
603b9aefa3SBarry Smith }
613b9aefa3SBarry Smith 
625a5d4f66SBarry Smith /*@C
635a5d4f66SBarry Smith     ISLocalToGlobalMappingView - View a local to global mapping
645a5d4f66SBarry Smith 
65b9cd556bSLois Curfman McInnes     Not Collective
66b9cd556bSLois Curfman McInnes 
675a5d4f66SBarry Smith     Input Parameters:
683b9aefa3SBarry Smith +   ltog - local to global mapping
693b9aefa3SBarry Smith -   viewer - viewer
705a5d4f66SBarry Smith 
71a997ad1aSLois Curfman McInnes     Level: advanced
72a997ad1aSLois Curfman McInnes 
73273d9f13SBarry Smith     Concepts: mapping^local to global
745a5d4f66SBarry Smith 
755a5d4f66SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
765a5d4f66SBarry Smith @*/
777087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingView(ISLocalToGlobalMapping mapping,PetscViewer viewer)
785a5d4f66SBarry Smith {
7932dcc486SBarry Smith   PetscInt       i;
8032dcc486SBarry Smith   PetscMPIInt    rank;
81ace3abfcSBarry Smith   PetscBool      iascii;
826849ba73SBarry Smith   PetscErrorCode ierr;
835a5d4f66SBarry Smith 
845a5d4f66SBarry Smith   PetscFunctionBegin;
850700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
863050cee2SBarry Smith   if (!viewer) {
87ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mapping),&viewer);CHKERRQ(ierr);
883050cee2SBarry Smith   }
890700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
905a5d4f66SBarry Smith 
91ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mapping),&rank);CHKERRQ(ierr);
92251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
9332077d6dSBarry Smith   if (iascii) {
9498c3331eSBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)mapping,viewer);CHKERRQ(ierr);
951575c14dSBarry Smith     ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
965a5d4f66SBarry Smith     for (i=0; i<mapping->n; i++) {
977904a332SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D %D\n",rank,i,mapping->indices[i]);CHKERRQ(ierr);
986831982aSBarry Smith     }
99b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1001575c14dSBarry Smith     ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
1011575c14dSBarry Smith   }
1025a5d4f66SBarry Smith   PetscFunctionReturn(0);
1035a5d4f66SBarry Smith }
1045a5d4f66SBarry Smith 
1051f428162SBarry Smith /*@
1062bdab257SBarry Smith     ISLocalToGlobalMappingCreateIS - Creates a mapping between a local (0 to n)
1072bdab257SBarry Smith     ordering and a global parallel ordering.
1082bdab257SBarry Smith 
1090f5bd95cSBarry Smith     Not collective
110b9cd556bSLois Curfman McInnes 
111a997ad1aSLois Curfman McInnes     Input Parameter:
1128c03b21aSDmitry Karpeev .   is - index set containing the global numbers for each local number
1132bdab257SBarry Smith 
114a997ad1aSLois Curfman McInnes     Output Parameter:
1152bdab257SBarry Smith .   mapping - new mapping data structure
1162bdab257SBarry Smith 
117f0413b6fSBarry Smith     Notes: the block size of the IS determines the block size of the mapping
118a997ad1aSLois Curfman McInnes     Level: advanced
119a997ad1aSLois Curfman McInnes 
120273d9f13SBarry Smith     Concepts: mapping^local to global
1212bdab257SBarry Smith 
122*7e99dc12SLawrence Mitchell .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetFromOptions()
1232bdab257SBarry Smith @*/
1247087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingCreateIS(IS is,ISLocalToGlobalMapping *mapping)
1252bdab257SBarry Smith {
1266849ba73SBarry Smith   PetscErrorCode ierr;
1273bbf0e92SBarry Smith   PetscInt       n,bs;
1285d0c19d7SBarry Smith   const PetscInt *indices;
1292bdab257SBarry Smith   MPI_Comm       comm;
1303bbf0e92SBarry Smith   PetscBool      isblock;
1313a40ed3dSBarry Smith 
1323a40ed3dSBarry Smith   PetscFunctionBegin;
1330700a824SBarry Smith   PetscValidHeaderSpecific(is,IS_CLASSID,1);
1344482741eSBarry Smith   PetscValidPointer(mapping,2);
1352bdab257SBarry Smith 
1362bdab257SBarry Smith   ierr = PetscObjectGetComm((PetscObject)is,&comm);CHKERRQ(ierr);
1373b9aefa3SBarry Smith   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
1383bbf0e92SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)is,ISBLOCK,&isblock);CHKERRQ(ierr);
1396006e8d2SBarry Smith   if (!isblock) {
140f0413b6fSBarry Smith     ierr = ISGetIndices(is,&indices);CHKERRQ(ierr);
141f0413b6fSBarry Smith     ierr = ISLocalToGlobalMappingCreate(comm,1,n,indices,PETSC_COPY_VALUES,mapping);CHKERRQ(ierr);
1422bdab257SBarry Smith     ierr = ISRestoreIndices(is,&indices);CHKERRQ(ierr);
1436006e8d2SBarry Smith   } else {
1446006e8d2SBarry Smith     ierr = ISGetBlockSize(is,&bs);CHKERRQ(ierr);
145f0413b6fSBarry Smith     ierr = ISBlockGetIndices(is,&indices);CHKERRQ(ierr);
14628bc9809SBarry Smith     ierr = ISLocalToGlobalMappingCreate(comm,bs,n/bs,indices,PETSC_COPY_VALUES,mapping);CHKERRQ(ierr);
147f0413b6fSBarry Smith     ierr = ISBlockRestoreIndices(is,&indices);CHKERRQ(ierr);
1486006e8d2SBarry Smith   }
1493a40ed3dSBarry Smith   PetscFunctionReturn(0);
1502bdab257SBarry Smith }
1515a5d4f66SBarry Smith 
152a4d96a55SJed Brown /*@C
153a4d96a55SJed Brown     ISLocalToGlobalMappingCreateSF - Creates a mapping between a local (0 to n)
154a4d96a55SJed Brown     ordering and a global parallel ordering.
155a4d96a55SJed Brown 
156a4d96a55SJed Brown     Collective
157a4d96a55SJed Brown 
158a4d96a55SJed Brown     Input Parameter:
159a4d96a55SJed Brown +   sf - star forest mapping contiguous local indices to (rank, offset)
160a4d96a55SJed Brown -   start - first global index on this process
161a4d96a55SJed Brown 
162a4d96a55SJed Brown     Output Parameter:
163a4d96a55SJed Brown .   mapping - new mapping data structure
164a4d96a55SJed Brown 
165a4d96a55SJed Brown     Level: advanced
166a4d96a55SJed Brown 
167a4d96a55SJed Brown     Concepts: mapping^local to global
168a4d96a55SJed Brown 
169*7e99dc12SLawrence Mitchell .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingSetFromOptions()
170a4d96a55SJed Brown @*/
171a4d96a55SJed Brown PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF sf,PetscInt start,ISLocalToGlobalMapping *mapping)
172a4d96a55SJed Brown {
173a4d96a55SJed Brown   PetscErrorCode ierr;
174a4d96a55SJed Brown   PetscInt       i,maxlocal,nroots,nleaves,*globals,*ltog;
175a4d96a55SJed Brown   const PetscInt *ilocal;
176a4d96a55SJed Brown   MPI_Comm       comm;
177a4d96a55SJed Brown 
178a4d96a55SJed Brown   PetscFunctionBegin;
179a4d96a55SJed Brown   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
180a4d96a55SJed Brown   PetscValidPointer(mapping,3);
181a4d96a55SJed Brown 
182a4d96a55SJed Brown   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
1830298fd71SBarry Smith   ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr);
184f6e5521dSKarl Rupp   if (ilocal) {
185f6e5521dSKarl Rupp     for (i=0,maxlocal=0; i<nleaves; i++) maxlocal = PetscMax(maxlocal,ilocal[i]+1);
186f6e5521dSKarl Rupp   }
187a4d96a55SJed Brown   else maxlocal = nleaves;
188785e854fSJed Brown   ierr = PetscMalloc1(nroots,&globals);CHKERRQ(ierr);
189785e854fSJed Brown   ierr = PetscMalloc1(maxlocal,&ltog);CHKERRQ(ierr);
190a4d96a55SJed Brown   for (i=0; i<nroots; i++) globals[i] = start + i;
191a4d96a55SJed Brown   for (i=0; i<maxlocal; i++) ltog[i] = -1;
192a4d96a55SJed Brown   ierr = PetscSFBcastBegin(sf,MPIU_INT,globals,ltog);CHKERRQ(ierr);
193a4d96a55SJed Brown   ierr = PetscSFBcastEnd(sf,MPIU_INT,globals,ltog);CHKERRQ(ierr);
194f0413b6fSBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,1,maxlocal,ltog,PETSC_OWN_POINTER,mapping);CHKERRQ(ierr);
195a4d96a55SJed Brown   ierr = PetscFree(globals);CHKERRQ(ierr);
196a4d96a55SJed Brown   PetscFunctionReturn(0);
197a4d96a55SJed Brown }
198b46b645bSBarry Smith 
19963fa5c83Sstefano_zampini /*@
20063fa5c83Sstefano_zampini     ISLocalToGlobalMappingSetBlockSize - Sets the blocksize of the mapping
20163fa5c83Sstefano_zampini 
20263fa5c83Sstefano_zampini     Not collective
20363fa5c83Sstefano_zampini 
20463fa5c83Sstefano_zampini     Input Parameters:
20563fa5c83Sstefano_zampini .   mapping - mapping data structure
20663fa5c83Sstefano_zampini .   bs - the blocksize
20763fa5c83Sstefano_zampini 
20863fa5c83Sstefano_zampini     Level: advanced
20963fa5c83Sstefano_zampini 
21063fa5c83Sstefano_zampini     Concepts: mapping^local to global
21163fa5c83Sstefano_zampini 
21263fa5c83Sstefano_zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS()
21363fa5c83Sstefano_zampini @*/
21463fa5c83Sstefano_zampini PetscErrorCode  ISLocalToGlobalMappingSetBlockSize(ISLocalToGlobalMapping mapping,PetscInt bs)
21563fa5c83Sstefano_zampini {
216a59f3c4dSstefano_zampini   PetscInt       *nid;
217a59f3c4dSstefano_zampini   const PetscInt *oid;
218a59f3c4dSstefano_zampini   PetscInt       i,cn,on,obs,nn;
21963fa5c83Sstefano_zampini   PetscErrorCode ierr;
22063fa5c83Sstefano_zampini 
22163fa5c83Sstefano_zampini   PetscFunctionBegin;
22263fa5c83Sstefano_zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
22363fa5c83Sstefano_zampini   if (bs < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid block size %D",bs);
22463fa5c83Sstefano_zampini   if (bs == mapping->bs) PetscFunctionReturn(0);
22563fa5c83Sstefano_zampini   on  = mapping->n;
22663fa5c83Sstefano_zampini   obs = mapping->bs;
22763fa5c83Sstefano_zampini   oid = mapping->indices;
22863fa5c83Sstefano_zampini   nn  = (on*obs)/bs;
22963fa5c83Sstefano_zampini   if ((on*obs)%bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Block size %D is inconsistent with block size %D and number of block indices %D",bs,obs,on);
230a59f3c4dSstefano_zampini 
23163fa5c83Sstefano_zampini   ierr = PetscMalloc1(nn,&nid);CHKERRQ(ierr);
232a59f3c4dSstefano_zampini   ierr = ISLocalToGlobalMappingGetIndices(mapping,&oid);CHKERRQ(ierr);
233a59f3c4dSstefano_zampini   for (i=0;i<nn;i++) {
234a59f3c4dSstefano_zampini     PetscInt j;
235a59f3c4dSstefano_zampini     for (j=0,cn=0;j<bs-1;j++) {
236a59f3c4dSstefano_zampini       if (oid[i*bs+j] < 0) { cn++; continue; }
237a59f3c4dSstefano_zampini       if (oid[i*bs+j] != oid[i*bs+j+1]-1) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Block sizes %D and %D are incompatible with the block indices: non consecutive indices %D %D",bs,obs,oid[i*bs+j],oid[i*bs+j+1]);
238a59f3c4dSstefano_zampini     }
239a59f3c4dSstefano_zampini     if (oid[i*bs+j] < 0) cn++;
2408b7cb0e6Sstefano_zampini     if (cn) {
241a59f3c4dSstefano_zampini       if (cn != bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Block sizes %D and %D are incompatible with the block indices: invalid number of negative entries in block %D",bs,obs,cn);
242a59f3c4dSstefano_zampini       nid[i] = -1;
2438b7cb0e6Sstefano_zampini     } else {
244a59f3c4dSstefano_zampini       nid[i] = oid[i*bs]/bs;
24563fa5c83Sstefano_zampini     }
24663fa5c83Sstefano_zampini   }
247a59f3c4dSstefano_zampini   ierr = ISLocalToGlobalMappingRestoreIndices(mapping,&oid);CHKERRQ(ierr);
248a59f3c4dSstefano_zampini 
24963fa5c83Sstefano_zampini   mapping->n           = nn;
25063fa5c83Sstefano_zampini   mapping->bs          = bs;
25163fa5c83Sstefano_zampini   ierr                 = PetscFree(mapping->indices);CHKERRQ(ierr);
25263fa5c83Sstefano_zampini   mapping->indices     = nid;
253c9345713Sstefano_zampini   ierr                 = PetscFree(mapping->globals);CHKERRQ(ierr);
254c9345713Sstefano_zampini   mapping->globalstart = 0;
255c9345713Sstefano_zampini   mapping->globalend   = 0;
25663fa5c83Sstefano_zampini   PetscFunctionReturn(0);
25763fa5c83Sstefano_zampini }
25863fa5c83Sstefano_zampini 
25945b6f7e9SBarry Smith /*@
26045b6f7e9SBarry Smith     ISLocalToGlobalMappingGetBlockSize - Gets the blocksize of the mapping
26145b6f7e9SBarry Smith     ordering and a global parallel ordering.
26245b6f7e9SBarry Smith 
26345b6f7e9SBarry Smith     Not Collective
26445b6f7e9SBarry Smith 
26545b6f7e9SBarry Smith     Input Parameters:
26645b6f7e9SBarry Smith .   mapping - mapping data structure
26745b6f7e9SBarry Smith 
26845b6f7e9SBarry Smith     Output Parameter:
26945b6f7e9SBarry Smith .   bs - the blocksize
27045b6f7e9SBarry Smith 
27145b6f7e9SBarry Smith     Level: advanced
27245b6f7e9SBarry Smith 
27345b6f7e9SBarry Smith     Concepts: mapping^local to global
27445b6f7e9SBarry Smith 
27545b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS()
27645b6f7e9SBarry Smith @*/
27745b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping mapping,PetscInt *bs)
27845b6f7e9SBarry Smith {
27945b6f7e9SBarry Smith   PetscFunctionBegin;
280cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
28145b6f7e9SBarry Smith   *bs = mapping->bs;
28245b6f7e9SBarry Smith   PetscFunctionReturn(0);
28345b6f7e9SBarry Smith }
28445b6f7e9SBarry Smith 
285ba5bb76aSSatish Balay /*@
28690f02eecSBarry Smith     ISLocalToGlobalMappingCreate - Creates a mapping between a local (0 to n)
28790f02eecSBarry Smith     ordering and a global parallel ordering.
2882362add9SBarry Smith 
28989d82c54SBarry Smith     Not Collective, but communicator may have more than one process
290b9cd556bSLois Curfman McInnes 
2912362add9SBarry Smith     Input Parameters:
29289d82c54SBarry Smith +   comm - MPI communicator
293f0413b6fSBarry Smith .   bs - the block size
29428bc9809SBarry Smith .   n - the number of local elements divided by the block size, or equivalently the number of block indices
29528bc9809SBarry Smith .   indices - the global index for each local element, these do not need to be in increasing order (sorted), these values should not be scaled (i.e. multiplied) by the blocksize bs
296d5ad8652SBarry Smith -   mode - see PetscCopyMode
2972362add9SBarry Smith 
298a997ad1aSLois Curfman McInnes     Output Parameter:
29990f02eecSBarry Smith .   mapping - new mapping data structure
3002362add9SBarry Smith 
301f0413b6fSBarry Smith     Notes: There is one integer value in indices per block and it represents the actual indices bs*idx + j, where j=0,..,bs-1
302a997ad1aSLois Curfman McInnes     Level: advanced
303a997ad1aSLois Curfman McInnes 
304273d9f13SBarry Smith     Concepts: mapping^local to global
3052362add9SBarry Smith 
306*7e99dc12SLawrence Mitchell .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingSetFromOptions()
3072362add9SBarry Smith @*/
30860c7cefcSBarry Smith PetscErrorCode  ISLocalToGlobalMappingCreate(MPI_Comm comm,PetscInt bs,PetscInt n,const PetscInt indices[],PetscCopyMode mode,ISLocalToGlobalMapping *mapping)
3092362add9SBarry Smith {
3106849ba73SBarry Smith   PetscErrorCode ierr;
31132dcc486SBarry Smith   PetscInt       *in;
312b46b645bSBarry Smith 
313b46b645bSBarry Smith   PetscFunctionBegin;
31473911063SBarry Smith   if (n) PetscValidIntPointer(indices,3);
3154482741eSBarry Smith   PetscValidPointer(mapping,4);
316b46b645bSBarry Smith 
3170298fd71SBarry Smith   *mapping = NULL;
318607a6623SBarry Smith   ierr = ISInitializePackage();CHKERRQ(ierr);
3192362add9SBarry Smith 
32073107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(*mapping,IS_LTOGM_CLASSID,"ISLocalToGlobalMapping","Local to global mapping","IS",
32160c7cefcSBarry Smith                            comm,ISLocalToGlobalMappingDestroy,ISLocalToGlobalMappingView);CHKERRQ(ierr);
322d4bb536fSBarry Smith   (*mapping)->n             = n;
323f0413b6fSBarry Smith   (*mapping)->bs            = bs;
324268a049cSStefano Zampini   (*mapping)->info_cached   = PETSC_FALSE;
325268a049cSStefano Zampini   (*mapping)->info_free     = PETSC_FALSE;
326268a049cSStefano Zampini   (*mapping)->info_procs    = NULL;
327268a049cSStefano Zampini   (*mapping)->info_numprocs = NULL;
328268a049cSStefano Zampini   (*mapping)->info_indices  = NULL;
329d4bb536fSBarry Smith   /*
330d4bb536fSBarry Smith     Do not create the global to local mapping. This is only created if
331d4bb536fSBarry Smith     ISGlobalToLocalMapping() is called
332d4bb536fSBarry Smith   */
333d4bb536fSBarry Smith   (*mapping)->globals        = 0;
334cd5a51dfSLawrence Mitchell   (*mapping)->globalht       = 0;
335*7e99dc12SLawrence Mitchell   (*mapping)->use_hash_table = PETSC_FALSE;
336d5ad8652SBarry Smith   if (mode == PETSC_COPY_VALUES) {
337785e854fSJed Brown     ierr = PetscMalloc1(n,&in);CHKERRQ(ierr);
338d5ad8652SBarry Smith     ierr = PetscMemcpy(in,indices,n*sizeof(PetscInt));CHKERRQ(ierr);
339d5ad8652SBarry Smith     (*mapping)->indices = in;
3406389a1a1SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt));CHKERRQ(ierr);
3416389a1a1SBarry Smith   } else if (mode == PETSC_OWN_POINTER) {
3426389a1a1SBarry Smith     (*mapping)->indices = (PetscInt*)indices;
3436389a1a1SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt));CHKERRQ(ierr);
3446389a1a1SBarry Smith   }
34560c7cefcSBarry Smith   else SETERRQ(comm,PETSC_ERR_SUP,"Cannot currently use PETSC_USE_POINTER");
346d96308ebSBarry Smith   ierr = PetscStrallocpy("basic",&((PetscObject)*mapping)->type_name);CHKERRQ(ierr);
3473a40ed3dSBarry Smith   PetscFunctionReturn(0);
3482362add9SBarry Smith }
3492362add9SBarry Smith 
35090f02eecSBarry Smith /*@
351*7e99dc12SLawrence Mitchell    ISLocalToGlobalMappingSetFromOptions - Set mapping options from the options database.
352*7e99dc12SLawrence Mitchell 
353*7e99dc12SLawrence Mitchell    Not collective
354*7e99dc12SLawrence Mitchell 
355*7e99dc12SLawrence Mitchell    Input Parameters:
356*7e99dc12SLawrence Mitchell .  mapping - mapping data structure
357*7e99dc12SLawrence Mitchell 
358*7e99dc12SLawrence Mitchell    Options Database:
359*7e99dc12SLawrence Mitchell .   -ltogm_use_hash_table true,false see ISLocalToGlobalMappingSetUseHashTable()
360*7e99dc12SLawrence Mitchell 
361*7e99dc12SLawrence Mitchell    Level: advanced
362*7e99dc12SLawrence Mitchell 
363*7e99dc12SLawrence Mitchell .seealso: ISLocalToGlobalMappingSetUseHashTable()
364*7e99dc12SLawrence Mitchell @*/
365*7e99dc12SLawrence Mitchell 
366*7e99dc12SLawrence Mitchell PetscErrorCode ISLocalToGlobalMappingSetFromOptions(ISLocalToGlobalMapping mapping)
367*7e99dc12SLawrence Mitchell {
368*7e99dc12SLawrence Mitchell   PetscErrorCode ierr;
369*7e99dc12SLawrence Mitchell   PetscBool      flg;
370*7e99dc12SLawrence Mitchell 
371*7e99dc12SLawrence Mitchell   PetscFunctionBegin;
372*7e99dc12SLawrence Mitchell   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
373*7e99dc12SLawrence Mitchell 
374*7e99dc12SLawrence Mitchell   ierr = PetscObjectOptionsBegin((PetscObject)mapping);CHKERRQ(ierr);
375*7e99dc12SLawrence Mitchell 
376*7e99dc12SLawrence Mitchell   ierr = ISLocalToGlobalMappingGetUseHashTable(mapping,&flg);CHKERRQ(ierr);
377*7e99dc12SLawrence Mitchell   ierr = PetscOptionsBool("-ltogm_use_hash_table", "use hash table (memory scalable for global to local in ISLtoGMappings","ISLocalToGlobalMappingSetUseHashTable",flg,&flg,NULL);CHKERRQ(ierr);
378*7e99dc12SLawrence Mitchell   ierr = ISLocalToGlobalMappingSetUseHashTable(mapping,flg);CHKERRQ(ierr);
379*7e99dc12SLawrence Mitchell   ierr = PetscOptionsEnd();CHKERRQ(ierr);
380*7e99dc12SLawrence Mitchell   PetscFunctionReturn(0);
381*7e99dc12SLawrence Mitchell }
382*7e99dc12SLawrence Mitchell 
383*7e99dc12SLawrence Mitchell /*@
38490f02eecSBarry Smith    ISLocalToGlobalMappingDestroy - Destroys a mapping between a local (0 to n)
38590f02eecSBarry Smith    ordering and a global parallel ordering.
38690f02eecSBarry Smith 
3870f5bd95cSBarry Smith    Note Collective
388b9cd556bSLois Curfman McInnes 
38990f02eecSBarry Smith    Input Parameters:
39090f02eecSBarry Smith .  mapping - mapping data structure
39190f02eecSBarry Smith 
392a997ad1aSLois Curfman McInnes    Level: advanced
393a997ad1aSLois Curfman McInnes 
3943acfe500SLois Curfman McInnes .seealso: ISLocalToGlobalMappingCreate()
39590f02eecSBarry Smith @*/
3966bf464f9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping *mapping)
39790f02eecSBarry Smith {
398dfbe8321SBarry Smith   PetscErrorCode ierr;
3995fd66863SKarl Rupp 
4003a40ed3dSBarry Smith   PetscFunctionBegin;
4016bf464f9SBarry Smith   if (!*mapping) PetscFunctionReturn(0);
4026bf464f9SBarry Smith   PetscValidHeaderSpecific((*mapping),IS_LTOGM_CLASSID,1);
403997056adSBarry Smith   if (--((PetscObject)(*mapping))->refct > 0) {*mapping = 0;PetscFunctionReturn(0);}
4046bf464f9SBarry Smith   ierr = PetscFree((*mapping)->indices);CHKERRQ(ierr);
4056bf464f9SBarry Smith   ierr = PetscFree((*mapping)->globals);CHKERRQ(ierr);
406268a049cSStefano Zampini   ierr = PetscFree((*mapping)->info_procs);CHKERRQ(ierr);
407268a049cSStefano Zampini   ierr = PetscFree((*mapping)->info_numprocs);CHKERRQ(ierr);
408cd5a51dfSLawrence Mitchell   PetscHashIDestroy((*mapping)->globalht);
409268a049cSStefano Zampini   if ((*mapping)->info_indices) {
410268a049cSStefano Zampini     PetscInt i;
411268a049cSStefano Zampini 
412268a049cSStefano Zampini     ierr = PetscFree(((*mapping)->info_indices)[0]);CHKERRQ(ierr);
413268a049cSStefano Zampini     for (i=1; i<(*mapping)->info_nproc; i++) {
414268a049cSStefano Zampini       ierr = PetscFree(((*mapping)->info_indices)[i]);CHKERRQ(ierr);
415268a049cSStefano Zampini     }
416268a049cSStefano Zampini     ierr = PetscFree((*mapping)->info_indices);CHKERRQ(ierr);
417268a049cSStefano Zampini   }
418d38fa0fbSBarry Smith   ierr     = PetscHeaderDestroy(mapping);CHKERRQ(ierr);
419992144d0SBarry Smith   *mapping = 0;
4203a40ed3dSBarry Smith   PetscFunctionReturn(0);
42190f02eecSBarry Smith }
42290f02eecSBarry Smith 
42390f02eecSBarry Smith /*@
4243acfe500SLois Curfman McInnes     ISLocalToGlobalMappingApplyIS - Creates from an IS in the local numbering
4253acfe500SLois Curfman McInnes     a new index set using the global numbering defined in an ISLocalToGlobalMapping
4263acfe500SLois Curfman McInnes     context.
42790f02eecSBarry Smith 
428b9cd556bSLois Curfman McInnes     Not collective
429b9cd556bSLois Curfman McInnes 
43090f02eecSBarry Smith     Input Parameters:
431b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
432b9cd556bSLois Curfman McInnes -   is - index set in local numbering
43390f02eecSBarry Smith 
43490f02eecSBarry Smith     Output Parameters:
43590f02eecSBarry Smith .   newis - index set in global numbering
43690f02eecSBarry Smith 
437a997ad1aSLois Curfman McInnes     Level: advanced
438a997ad1aSLois Curfman McInnes 
439273d9f13SBarry Smith     Concepts: mapping^local to global
4403acfe500SLois Curfman McInnes 
44190f02eecSBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),
442d4bb536fSBarry Smith           ISLocalToGlobalMappingDestroy(), ISGlobalToLocalMappingApply()
44390f02eecSBarry Smith @*/
4447087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping,IS is,IS *newis)
44590f02eecSBarry Smith {
4466849ba73SBarry Smith   PetscErrorCode ierr;
447e24637baSBarry Smith   PetscInt       n,*idxout;
4485d0c19d7SBarry Smith   const PetscInt *idxin;
4493a40ed3dSBarry Smith 
4503a40ed3dSBarry Smith   PetscFunctionBegin;
4510700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
4520700a824SBarry Smith   PetscValidHeaderSpecific(is,IS_CLASSID,2);
4534482741eSBarry Smith   PetscValidPointer(newis,3);
45490f02eecSBarry Smith 
4553b9aefa3SBarry Smith   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
45690f02eecSBarry Smith   ierr = ISGetIndices(is,&idxin);CHKERRQ(ierr);
457785e854fSJed Brown   ierr = PetscMalloc1(n,&idxout);CHKERRQ(ierr);
458e24637baSBarry Smith   ierr = ISLocalToGlobalMappingApply(mapping,n,idxin,idxout);CHKERRQ(ierr);
4593b9aefa3SBarry Smith   ierr = ISRestoreIndices(is,&idxin);CHKERRQ(ierr);
460543f3098SMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),n,idxout,PETSC_OWN_POINTER,newis);CHKERRQ(ierr);
4613a40ed3dSBarry Smith   PetscFunctionReturn(0);
46290f02eecSBarry Smith }
46390f02eecSBarry Smith 
464b89cb25eSSatish Balay /*@
4653acfe500SLois Curfman McInnes    ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering
4663acfe500SLois Curfman McInnes    and converts them to the global numbering.
46790f02eecSBarry Smith 
468b9cd556bSLois Curfman McInnes    Not collective
469b9cd556bSLois Curfman McInnes 
470bb25748dSBarry Smith    Input Parameters:
471b9cd556bSLois Curfman McInnes +  mapping - the local to global mapping context
472bb25748dSBarry Smith .  N - number of integers
473b9cd556bSLois Curfman McInnes -  in - input indices in local numbering
474bb25748dSBarry Smith 
475bb25748dSBarry Smith    Output Parameter:
476bb25748dSBarry Smith .  out - indices in global numbering
477bb25748dSBarry Smith 
478b9cd556bSLois Curfman McInnes    Notes:
479b9cd556bSLois Curfman McInnes    The in and out array parameters may be identical.
480d4bb536fSBarry Smith 
481a997ad1aSLois Curfman McInnes    Level: advanced
482a997ad1aSLois Curfman McInnes 
48345b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
4840752156aSBarry Smith           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
485d4bb536fSBarry Smith           AOPetscToApplication(), ISGlobalToLocalMappingApply()
486bb25748dSBarry Smith 
487273d9f13SBarry Smith     Concepts: mapping^local to global
488afcb2eb5SJed Brown @*/
489afcb2eb5SJed Brown PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
490afcb2eb5SJed Brown {
491cbc1caf0SMatthew G. Knepley   PetscInt i,bs,Nmax;
49245b6f7e9SBarry Smith 
49345b6f7e9SBarry Smith   PetscFunctionBegin;
494cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
495cbc1caf0SMatthew G. Knepley   bs   = mapping->bs;
496cbc1caf0SMatthew G. Knepley   Nmax = bs*mapping->n;
49745b6f7e9SBarry Smith   if (bs == 1) {
498cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
49945b6f7e9SBarry Smith     for (i=0; i<N; i++) {
50045b6f7e9SBarry Smith       if (in[i] < 0) {
50145b6f7e9SBarry Smith         out[i] = in[i];
50245b6f7e9SBarry Smith         continue;
50345b6f7e9SBarry Smith       }
504e24637baSBarry Smith       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);
50545b6f7e9SBarry Smith       out[i] = idx[in[i]];
50645b6f7e9SBarry Smith     }
50745b6f7e9SBarry Smith   } else {
508cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
50945b6f7e9SBarry Smith     for (i=0; i<N; i++) {
51045b6f7e9SBarry Smith       if (in[i] < 0) {
51145b6f7e9SBarry Smith         out[i] = in[i];
51245b6f7e9SBarry Smith         continue;
51345b6f7e9SBarry Smith       }
514e24637baSBarry Smith       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);
51545b6f7e9SBarry Smith       out[i] = idx[in[i]/bs]*bs + (in[i] % bs);
51645b6f7e9SBarry Smith     }
51745b6f7e9SBarry Smith   }
51845b6f7e9SBarry Smith   PetscFunctionReturn(0);
51945b6f7e9SBarry Smith }
52045b6f7e9SBarry Smith 
52145b6f7e9SBarry Smith /*@
5226006e8d2SBarry Smith    ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local block numbering and converts them to the global block numbering
52345b6f7e9SBarry Smith 
52445b6f7e9SBarry Smith    Not collective
52545b6f7e9SBarry Smith 
52645b6f7e9SBarry Smith    Input Parameters:
52745b6f7e9SBarry Smith +  mapping - the local to global mapping context
52845b6f7e9SBarry Smith .  N - number of integers
5296006e8d2SBarry Smith -  in - input indices in local block numbering
53045b6f7e9SBarry Smith 
53145b6f7e9SBarry Smith    Output Parameter:
5326006e8d2SBarry Smith .  out - indices in global block numbering
53345b6f7e9SBarry Smith 
53445b6f7e9SBarry Smith    Notes:
53545b6f7e9SBarry Smith    The in and out array parameters may be identical.
53645b6f7e9SBarry Smith 
5376006e8d2SBarry Smith    Example:
5386006e8d2SBarry Smith      If the index values are {0,1,6,7} set with a call to ISLocalToGlobalMappingCreate(PETSC_COMM_SELF,2,2,{0,3}) then the mapping applied to 0
5396006e8d2SBarry Smith      (the first block) would produce 0 and the mapping applied to 1 (the second block) would produce 3.
5406006e8d2SBarry Smith 
54145b6f7e9SBarry Smith    Level: advanced
54245b6f7e9SBarry Smith 
54345b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
54445b6f7e9SBarry Smith           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
54545b6f7e9SBarry Smith           AOPetscToApplication(), ISGlobalToLocalMappingApply()
54645b6f7e9SBarry Smith 
54745b6f7e9SBarry Smith     Concepts: mapping^local to global
54845b6f7e9SBarry Smith @*/
54945b6f7e9SBarry Smith PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
55045b6f7e9SBarry Smith {
551cbc1caf0SMatthew G. Knepley 
552cbc1caf0SMatthew G. Knepley   PetscFunctionBegin;
553cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
554cbc1caf0SMatthew G. Knepley   {
555afcb2eb5SJed Brown     PetscInt i,Nmax = mapping->n;
556afcb2eb5SJed Brown     const PetscInt *idx = mapping->indices;
557d4bb536fSBarry Smith 
558afcb2eb5SJed Brown     for (i=0; i<N; i++) {
559afcb2eb5SJed Brown       if (in[i] < 0) {
560afcb2eb5SJed Brown         out[i] = in[i];
561afcb2eb5SJed Brown         continue;
562afcb2eb5SJed Brown       }
563e24637baSBarry Smith       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);
564afcb2eb5SJed Brown       out[i] = idx[in[i]];
565afcb2eb5SJed Brown     }
566cbc1caf0SMatthew G. Knepley   }
567afcb2eb5SJed Brown   PetscFunctionReturn(0);
568afcb2eb5SJed Brown }
569d4bb536fSBarry Smith 
570d4bb536fSBarry Smith /* -----------------------------------------------------------------------------------------*/
571d4bb536fSBarry Smith 
572d4bb536fSBarry Smith /*
573d4bb536fSBarry Smith     Creates the global fields in the ISLocalToGlobalMapping structure
574d4bb536fSBarry Smith */
5756849ba73SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp_Private(ISLocalToGlobalMapping mapping)
576d4bb536fSBarry Smith {
5776849ba73SBarry Smith   PetscErrorCode ierr;
57832dcc486SBarry Smith   PetscInt       i,*idx = mapping->indices,n = mapping->n,end,start,*globals;
579*7e99dc12SLawrence Mitchell   PetscBool      use_hash_table;
580d4bb536fSBarry Smith 
5813a40ed3dSBarry Smith   PetscFunctionBegin;
582*7e99dc12SLawrence Mitchell   ierr = ISLocalToGlobalMappingGetUseHashTable(mapping,&use_hash_table);CHKERRQ(ierr);
583*7e99dc12SLawrence Mitchell   if (use_hash_table && mapping->globalht) PetscFunctionReturn(0);
584cd5a51dfSLawrence Mitchell   else if (mapping->globals) PetscFunctionReturn(0);
585d4bb536fSBarry Smith   end   = 0;
586ec268f7cSJed Brown   start = PETSC_MAX_INT;
587d4bb536fSBarry Smith 
588d4bb536fSBarry Smith   for (i=0; i<n; i++) {
589d4bb536fSBarry Smith     if (idx[i] < 0) continue;
590d4bb536fSBarry Smith     if (idx[i] < start) start = idx[i];
591d4bb536fSBarry Smith     if (idx[i] > end)   end   = idx[i];
592d4bb536fSBarry Smith   }
593d4bb536fSBarry Smith   if (start > end) {start = 0; end = -1;}
594d4bb536fSBarry Smith   mapping->globalstart = start;
595d4bb536fSBarry Smith   mapping->globalend   = end;
596d4bb536fSBarry Smith 
597cd5a51dfSLawrence Mitchell   if (mapping->use_hash_table) {
598cd5a51dfSLawrence Mitchell     PetscHashICreate(mapping->globalht);
599cd5a51dfSLawrence Mitchell     for (i=0; i<n; i++ ) {
600cd5a51dfSLawrence Mitchell       if (idx[i] < 0) continue;
601cd5a51dfSLawrence Mitchell       PetscHashIAdd(mapping->globalht, idx[i], i);
602cd5a51dfSLawrence Mitchell     }
603cd5a51dfSLawrence Mitchell     ierr = PetscLogObjectMemory((PetscObject)mapping,2*n*sizeof(PetscInt));CHKERRQ(ierr);
604cd5a51dfSLawrence Mitchell   } else {
605854ce69bSBarry Smith     ierr             = PetscMalloc1(end-start+2,&globals);CHKERRQ(ierr);
606b0a32e0cSBarry Smith     mapping->globals = globals;
607f6e5521dSKarl Rupp     for (i=0; i<end-start+1; i++) globals[i] = -1;
608d4bb536fSBarry Smith     for (i=0; i<n; i++) {
609d4bb536fSBarry Smith       if (idx[i] < 0) continue;
610d4bb536fSBarry Smith       globals[idx[i] - start] = i;
611d4bb536fSBarry Smith     }
612d4bb536fSBarry Smith 
6133bb1ff40SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)mapping,(end-start+1)*sizeof(PetscInt));CHKERRQ(ierr);
614cd5a51dfSLawrence Mitchell   }
6153a40ed3dSBarry Smith   PetscFunctionReturn(0);
616d4bb536fSBarry Smith }
617d4bb536fSBarry Smith 
618d4bb536fSBarry Smith /*@
619*7e99dc12SLawrence Mitchell     ISLocalToGlobalMappingSetUseHashTable - Should global to local lookups use a hash table?
620*7e99dc12SLawrence Mitchell 
621*7e99dc12SLawrence Mitchell     This is memory-scalable, but slightly slower, implementation than the default flat array.
622*7e99dc12SLawrence Mitchell 
623*7e99dc12SLawrence Mitchell     Not collective
624*7e99dc12SLawrence Mitchell 
625*7e99dc12SLawrence Mitchell     Input parameters:
626*7e99dc12SLawrence Mitchell +   mapping - mapping between local and global numbering
627*7e99dc12SLawrence Mitchell -   flg - use a hash table?
628*7e99dc12SLawrence Mitchell 
629*7e99dc12SLawrence Mitchell     Level: advanced
630*7e99dc12SLawrence Mitchell 
631*7e99dc12SLawrence Mitchell .seealso: ISGlobalToLocalMappingApply(), ISGlobalToLocalMappingApplyBlock(), ISLocalToGlobalMappingSetFromOptions(), ISLocalToGlobalMappingGetUseHashTable()
632*7e99dc12SLawrence Mitchell 
633*7e99dc12SLawrence Mitchell */
634*7e99dc12SLawrence Mitchell PetscErrorCode ISLocalToGlobalMappingSetUseHashTable(ISLocalToGlobalMapping mapping, PetscBool flg)
635*7e99dc12SLawrence Mitchell {
636*7e99dc12SLawrence Mitchell   PetscFunctionBegin;
637*7e99dc12SLawrence Mitchell 
638*7e99dc12SLawrence Mitchell   mapping->use_hash_table = flg;
639*7e99dc12SLawrence Mitchell 
640*7e99dc12SLawrence Mitchell   PetscFunctionReturn(0);
641*7e99dc12SLawrence Mitchell }
642*7e99dc12SLawrence Mitchell 
643*7e99dc12SLawrence Mitchell /*@
644*7e99dc12SLawrence Mitchell     ISLocalToGlobalMappingGetUseHashTable - Do global to local lookups use a hash table?
645*7e99dc12SLawrence Mitchell 
646*7e99dc12SLawrence Mitchell     This is memory-scalable, but slightly slower, implementation than the default flat array.
647*7e99dc12SLawrence Mitchell 
648*7e99dc12SLawrence Mitchell     Not collective
649*7e99dc12SLawrence Mitchell 
650*7e99dc12SLawrence Mitchell     Input parameter:
651*7e99dc12SLawrence Mitchell +   mapping - mapping between local and global numbering
652*7e99dc12SLawrence Mitchell 
653*7e99dc12SLawrence Mitchell     Output parameter:
654*7e99dc12SLawrence Mitchell +   flg - does the mapping use a hash table?
655*7e99dc12SLawrence Mitchell 
656*7e99dc12SLawrence Mitchell     Level: advanced
657*7e99dc12SLawrence Mitchell 
658*7e99dc12SLawrence Mitchell .seealso: ISGlobalToLocalMappingApply(), ISGlobalToLocalMappingApplyBlock(), ISLocalToGlobalMappingSetFromOptions(), ISLocalToGlobalMappingSetUseHashTable()
659*7e99dc12SLawrence Mitchell */
660*7e99dc12SLawrence Mitchell PetscErrorCode ISLocalToGlobalMappingGetUseHashTable(ISLocalToGlobalMapping mapping, PetscBool *flg)
661*7e99dc12SLawrence Mitchell {
662*7e99dc12SLawrence Mitchell   PetscFunctionBegin;
663*7e99dc12SLawrence Mitchell 
664*7e99dc12SLawrence Mitchell   *flg = mapping->use_hash_table;
665*7e99dc12SLawrence Mitchell 
666*7e99dc12SLawrence Mitchell   PetscFunctionReturn(0);
667*7e99dc12SLawrence Mitchell }
668*7e99dc12SLawrence Mitchell 
669*7e99dc12SLawrence Mitchell /*@
670a997ad1aSLois Curfman McInnes     ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers
671a997ad1aSLois Curfman McInnes     specified with a global numbering.
672d4bb536fSBarry Smith 
673b9cd556bSLois Curfman McInnes     Not collective
674b9cd556bSLois Curfman McInnes 
675d4bb536fSBarry Smith     Input Parameters:
676b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
677d4bb536fSBarry Smith .   type - IS_GTOLM_MASK - replaces global indices with no local value with -1
678d4bb536fSBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
679d4bb536fSBarry Smith .   n - number of global indices to map
680b9cd556bSLois Curfman McInnes -   idx - global indices to map
681d4bb536fSBarry Smith 
682d4bb536fSBarry Smith     Output Parameters:
683b9cd556bSLois Curfman McInnes +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
684b9cd556bSLois Curfman McInnes -   idxout - local index of each global index, one must pass in an array long enough
685e182c471SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApply() with
6860298fd71SBarry Smith              idxout == NULL to determine the required length (returned in nout)
687e182c471SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApply()
688e182c471SBarry Smith              a second time to set the values.
689d4bb536fSBarry Smith 
690b9cd556bSLois Curfman McInnes     Notes:
6910298fd71SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
692d4bb536fSBarry Smith 
693*7e99dc12SLawrence Mitchell     Unless the option -ltogm_use_hash_table is used, this is not scalable in memory usage. Each processor requires
694*7e99dc12SLawrence Mitchell     O(Nglobal) size array to compute these.
6950f5bd95cSBarry Smith 
696a997ad1aSLois Curfman McInnes     Level: advanced
697a997ad1aSLois Curfman McInnes 
69832fd6b96SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
69932fd6b96SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
70032fd6b96SBarry Smith 
701273d9f13SBarry Smith     Concepts: mapping^global to local
702d4bb536fSBarry Smith 
7039d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),
704*7e99dc12SLawrence Mitchell           ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingGetUseHashTable(), ISLocalToGlobalMappingSetUseHashTable()
705d4bb536fSBarry Smith @*/
7067087cfbeSBarry Smith PetscErrorCode  ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingType type,
70732dcc486SBarry Smith                                             PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
708d4bb536fSBarry Smith {
709cd5a51dfSLawrence Mitchell   PetscInt       i,*globals,nf = 0,tmp,start,end,bs,local;
710*7e99dc12SLawrence Mitchell   PetscBool      use_hash_table;
7119d90f715SBarry Smith   PetscErrorCode ierr;
7129d90f715SBarry Smith 
7139d90f715SBarry Smith   PetscFunctionBegin;
7149d90f715SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
7159d90f715SBarry Smith   ierr = ISGlobalToLocalMappingSetUp_Private(mapping);CHKERRQ(ierr);
716*7e99dc12SLawrence Mitchell   ierr = ISLocalToGlobalMappingGetUseHashTable(mapping,&use_hash_table);CHKERRQ(ierr);
7179d90f715SBarry Smith   globals = mapping->globals;
7189d90f715SBarry Smith   start   = mapping->globalstart;
7199d90f715SBarry Smith   end     = mapping->globalend;
7209d90f715SBarry Smith   bs      = mapping->bs;
7219d90f715SBarry Smith 
722cd5a51dfSLawrence Mitchell #define GTOL(g, local) do {                        \
723cd5a51dfSLawrence Mitchell     if (use_hash_table) {                          \
724cd5a51dfSLawrence Mitchell       PetscHashIMap(mapping->globalht,g/bs,local); \
725cd5a51dfSLawrence Mitchell     } else {                                       \
726cd5a51dfSLawrence Mitchell       local = globals[g/bs - start];               \
727cd5a51dfSLawrence Mitchell     }                                              \
728cd5a51dfSLawrence Mitchell     local = bs*local + (g % bs);                   \
729cd5a51dfSLawrence Mitchell   } while (0)
730cd5a51dfSLawrence Mitchell 
7319d90f715SBarry Smith   if (type == IS_GTOLM_MASK) {
7329d90f715SBarry Smith     if (idxout) {
7339d90f715SBarry Smith       for (i=0; i<n; i++) {
7349d90f715SBarry Smith         if (idx[i] < 0)                   idxout[i] = idx[i];
7359d90f715SBarry Smith         else if (idx[i] < bs*start)       idxout[i] = -1;
736663bb84eSBarry Smith         else if (idx[i] > bs*(end+1)-1)   idxout[i] = -1;
737cd5a51dfSLawrence Mitchell         else                              GTOL(idx[i], local);
7389d90f715SBarry Smith       }
7399d90f715SBarry Smith     }
7409d90f715SBarry Smith     if (nout) *nout = n;
7419d90f715SBarry Smith   } else {
7429d90f715SBarry Smith     if (idxout) {
7439d90f715SBarry Smith       for (i=0; i<n; i++) {
7449d90f715SBarry Smith         if (idx[i] < 0) continue;
7459d90f715SBarry Smith         if (idx[i] < bs*start) continue;
746663bb84eSBarry Smith         if (idx[i] > bs*(end+1)-1) continue;
747cd5a51dfSLawrence Mitchell         GTOL(idx[i], tmp);
7489d90f715SBarry Smith         if (tmp < 0) continue;
7499d90f715SBarry Smith         idxout[nf++] = tmp;
7509d90f715SBarry Smith       }
7519d90f715SBarry Smith     } else {
7529d90f715SBarry Smith       for (i=0; i<n; i++) {
7539d90f715SBarry Smith         if (idx[i] < 0) continue;
7549d90f715SBarry Smith         if (idx[i] < bs*start) continue;
755663bb84eSBarry Smith         if (idx[i] > bs*(end+1)-1) continue;
756cd5a51dfSLawrence Mitchell         GTOL(idx[i], tmp);
7579d90f715SBarry Smith         if (tmp < 0) continue;
7589d90f715SBarry Smith         nf++;
7599d90f715SBarry Smith       }
7609d90f715SBarry Smith     }
7619d90f715SBarry Smith     if (nout) *nout = nf;
7629d90f715SBarry Smith   }
763cd5a51dfSLawrence Mitchell #undef GTOL
7649d90f715SBarry Smith   PetscFunctionReturn(0);
7659d90f715SBarry Smith }
7669d90f715SBarry Smith 
767d4fe737eSStefano Zampini /*@
768d4fe737eSStefano Zampini     ISGlobalToLocalMappingApplyIS - Creates from an IS in the global numbering
769d4fe737eSStefano Zampini     a new index set using the local numbering defined in an ISLocalToGlobalMapping
770d4fe737eSStefano Zampini     context.
771d4fe737eSStefano Zampini 
772d4fe737eSStefano Zampini     Not collective
773d4fe737eSStefano Zampini 
774d4fe737eSStefano Zampini     Input Parameters:
775d4fe737eSStefano Zampini +   mapping - mapping between local and global numbering
776d4fe737eSStefano Zampini -   is - index set in global numbering
777d4fe737eSStefano Zampini 
778d4fe737eSStefano Zampini     Output Parameters:
779d4fe737eSStefano Zampini .   newis - index set in local numbering
780d4fe737eSStefano Zampini 
781d4fe737eSStefano Zampini     Level: advanced
782d4fe737eSStefano Zampini 
783d4fe737eSStefano Zampini     Concepts: mapping^local to global
784d4fe737eSStefano Zampini 
785d4fe737eSStefano Zampini .seealso: ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(),
786d4fe737eSStefano Zampini           ISLocalToGlobalMappingDestroy()
787d4fe737eSStefano Zampini @*/
788d4fe737eSStefano Zampini PetscErrorCode  ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingType type, IS is,IS *newis)
789d4fe737eSStefano Zampini {
790d4fe737eSStefano Zampini   PetscErrorCode ierr;
791d4fe737eSStefano Zampini   PetscInt       n,nout,*idxout;
792d4fe737eSStefano Zampini   const PetscInt *idxin;
793d4fe737eSStefano Zampini 
794d4fe737eSStefano Zampini   PetscFunctionBegin;
795d4fe737eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
796d4fe737eSStefano Zampini   PetscValidHeaderSpecific(is,IS_CLASSID,3);
797d4fe737eSStefano Zampini   PetscValidPointer(newis,4);
798d4fe737eSStefano Zampini 
799d4fe737eSStefano Zampini   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
800d4fe737eSStefano Zampini   ierr = ISGetIndices(is,&idxin);CHKERRQ(ierr);
801d4fe737eSStefano Zampini   if (type == IS_GTOLM_MASK) {
802d4fe737eSStefano Zampini     ierr = PetscMalloc1(n,&idxout);CHKERRQ(ierr);
803d4fe737eSStefano Zampini   } else {
804d4fe737eSStefano Zampini     ierr = ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,NULL);CHKERRQ(ierr);
805d4fe737eSStefano Zampini     ierr = PetscMalloc1(nout,&idxout);CHKERRQ(ierr);
806d4fe737eSStefano Zampini   }
807d4fe737eSStefano Zampini   ierr = ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,idxout);CHKERRQ(ierr);
808d4fe737eSStefano Zampini   ierr = ISRestoreIndices(is,&idxin);CHKERRQ(ierr);
809d4fe737eSStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),nout,idxout,PETSC_OWN_POINTER,newis);CHKERRQ(ierr);
810d4fe737eSStefano Zampini   PetscFunctionReturn(0);
811d4fe737eSStefano Zampini }
812d4fe737eSStefano Zampini 
8139d90f715SBarry Smith /*@
8149d90f715SBarry Smith     ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers
8159d90f715SBarry Smith     specified with a block global numbering.
8169d90f715SBarry Smith 
8179d90f715SBarry Smith     Not collective
8189d90f715SBarry Smith 
8199d90f715SBarry Smith     Input Parameters:
8209d90f715SBarry Smith +   mapping - mapping between local and global numbering
8219d90f715SBarry Smith .   type - IS_GTOLM_MASK - replaces global indices with no local value with -1
8229d90f715SBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
8239d90f715SBarry Smith .   n - number of global indices to map
8249d90f715SBarry Smith -   idx - global indices to map
8259d90f715SBarry Smith 
8269d90f715SBarry Smith     Output Parameters:
8279d90f715SBarry Smith +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
8289d90f715SBarry Smith -   idxout - local index of each global index, one must pass in an array long enough
8299d90f715SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApplyBlock() with
8309d90f715SBarry Smith              idxout == NULL to determine the required length (returned in nout)
8319d90f715SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApplyBlock()
8329d90f715SBarry Smith              a second time to set the values.
8339d90f715SBarry Smith 
8349d90f715SBarry Smith     Notes:
8359d90f715SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
8369d90f715SBarry Smith 
837*7e99dc12SLawrence Mitchell     Unless the option -ltogm_use_hash_table is used, this is not scalable in memory usage. Each processor requires
838*7e99dc12SLawrence Mitchell     O(Nglobal) size array to compute these.
8399d90f715SBarry Smith 
8409d90f715SBarry Smith     Level: advanced
8419d90f715SBarry Smith 
8429d90f715SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
8439d90f715SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
8449d90f715SBarry Smith 
8459d90f715SBarry Smith     Concepts: mapping^global to local
8469d90f715SBarry Smith 
8479d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(),
848*7e99dc12SLawrence Mitchell           ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingGetUseHashTable(), ISLocalToGlobalMappingSetUseHashTable()
8499d90f715SBarry Smith @*/
8509d90f715SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingType type,
8519d90f715SBarry Smith                                   PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
8529d90f715SBarry Smith {
85332dcc486SBarry Smith   PetscInt       i,*globals,nf = 0,tmp,start,end;
854*7e99dc12SLawrence Mitchell   PetscBool      use_hash_table;
8556849ba73SBarry Smith   PetscErrorCode ierr;
856d4bb536fSBarry Smith 
8573a40ed3dSBarry Smith   PetscFunctionBegin;
8580700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
859d4bb536fSBarry Smith   ierr = ISGlobalToLocalMappingSetUp_Private(mapping);CHKERRQ(ierr);
860*7e99dc12SLawrence Mitchell   ierr = ISLocalToGlobalMappingGetUseHashTable(mapping,&use_hash_table);CHKERRQ(ierr);
861d4bb536fSBarry Smith   globals = mapping->globals;
862d4bb536fSBarry Smith   start   = mapping->globalstart;
863d4bb536fSBarry Smith   end     = mapping->globalend;
864d4bb536fSBarry Smith 
865cd5a51dfSLawrence Mitchell #define GTOL(g, local) do {                     \
866cd5a51dfSLawrence Mitchell     if (use_hash_table) {                       \
867cd5a51dfSLawrence Mitchell       PetscHashIMap(mapping->globalht,g,local); \
868cd5a51dfSLawrence Mitchell     } else {                                    \
869cd5a51dfSLawrence Mitchell       local = globals[g - start];               \
870cd5a51dfSLawrence Mitchell     }                                           \
871cd5a51dfSLawrence Mitchell   } while (0)
872cd5a51dfSLawrence Mitchell 
873d4bb536fSBarry Smith   if (type == IS_GTOLM_MASK) {
874d4bb536fSBarry Smith     if (idxout) {
875d4bb536fSBarry Smith       for (i=0; i<n; i++) {
876d4bb536fSBarry Smith         if (idx[i] < 0) idxout[i] = idx[i];
877d4bb536fSBarry Smith         else if (idx[i] < start) idxout[i] = -1;
878d4bb536fSBarry Smith         else if (idx[i] > end)   idxout[i] = -1;
879cd5a51dfSLawrence Mitchell         else                     GTOL(idx[i], idxout[i]);
880d4bb536fSBarry Smith       }
881d4bb536fSBarry Smith     }
882d4bb536fSBarry Smith     if (nout) *nout = n;
883d4bb536fSBarry Smith   } else {
884d4bb536fSBarry Smith     if (idxout) {
885d4bb536fSBarry Smith       for (i=0; i<n; i++) {
886d4bb536fSBarry Smith         if (idx[i] < 0) continue;
887d4bb536fSBarry Smith         if (idx[i] < start) continue;
888d4bb536fSBarry Smith         if (idx[i] > end) continue;
889cd5a51dfSLawrence Mitchell         GTOL(idx[i], tmp);
890d4bb536fSBarry Smith         if (tmp < 0) continue;
891d4bb536fSBarry Smith         idxout[nf++] = tmp;
892d4bb536fSBarry Smith       }
893d4bb536fSBarry Smith     } else {
894d4bb536fSBarry Smith       for (i=0; i<n; i++) {
895d4bb536fSBarry Smith         if (idx[i] < 0) continue;
896d4bb536fSBarry Smith         if (idx[i] < start) continue;
897d4bb536fSBarry Smith         if (idx[i] > end) continue;
898cd5a51dfSLawrence Mitchell         GTOL(idx[i], tmp);
899d4bb536fSBarry Smith         if (tmp < 0) continue;
900d4bb536fSBarry Smith         nf++;
901d4bb536fSBarry Smith       }
902d4bb536fSBarry Smith     }
903d4bb536fSBarry Smith     if (nout) *nout = nf;
904d4bb536fSBarry Smith   }
905cd5a51dfSLawrence Mitchell #undef GTOL
9063a40ed3dSBarry Smith   PetscFunctionReturn(0);
907d4bb536fSBarry Smith }
90890f02eecSBarry Smith 
90989d82c54SBarry Smith /*@C
9106a818285SBarry Smith     ISLocalToGlobalMappingGetBlockInfo - Gets the neighbor information for each processor and
91189d82c54SBarry Smith      each index shared by more than one processor
91289d82c54SBarry Smith 
91389d82c54SBarry Smith     Collective on ISLocalToGlobalMapping
91489d82c54SBarry Smith 
91589d82c54SBarry Smith     Input Parameters:
91689d82c54SBarry Smith .   mapping - the mapping from local to global indexing
91789d82c54SBarry Smith 
91889d82c54SBarry Smith     Output Parameter:
91989d82c54SBarry Smith +   nproc - number of processors that are connected to this one
92089d82c54SBarry Smith .   proc - neighboring processors
92107b52d57SBarry Smith .   numproc - number of indices for each subdomain (processor)
9223463a7baSJed Brown -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
92389d82c54SBarry Smith 
92489d82c54SBarry Smith     Level: advanced
92589d82c54SBarry Smith 
926273d9f13SBarry Smith     Concepts: mapping^local to global
92789d82c54SBarry Smith 
9282cfcea29SBarry Smith     Fortran Usage:
9292cfcea29SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
9302cfcea29SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
9312cfcea29SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
9322cfcea29SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
9332cfcea29SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
9342cfcea29SBarry Smith 
9352cfcea29SBarry Smith 
93607b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
93707b52d57SBarry Smith           ISLocalToGlobalMappingRestoreInfo()
93889d82c54SBarry Smith @*/
9396a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
94089d82c54SBarry Smith {
9416849ba73SBarry Smith   PetscErrorCode ierr;
942268a049cSStefano Zampini 
943268a049cSStefano Zampini   PetscFunctionBegin;
944268a049cSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
945268a049cSStefano Zampini   if (mapping->info_cached) {
946268a049cSStefano Zampini     *nproc = mapping->info_nproc;
947268a049cSStefano Zampini     *procs = mapping->info_procs;
948268a049cSStefano Zampini     *numprocs = mapping->info_numprocs;
949268a049cSStefano Zampini     *indices = mapping->info_indices;
950268a049cSStefano Zampini   } else {
951268a049cSStefano Zampini     ierr = ISLocalToGlobalMappingGetBlockInfo_Private(mapping,nproc,procs,numprocs,indices);CHKERRQ(ierr);
952268a049cSStefano Zampini   }
953268a049cSStefano Zampini   PetscFunctionReturn(0);
954268a049cSStefano Zampini }
955268a049cSStefano Zampini 
956268a049cSStefano Zampini static PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
957268a049cSStefano Zampini {
958268a049cSStefano Zampini   PetscErrorCode ierr;
95997f1f81fSBarry Smith   PetscMPIInt    size,rank,tag1,tag2,tag3,*len,*source,imdex;
96032dcc486SBarry Smith   PetscInt       i,n = mapping->n,Ng,ng,max = 0,*lindices = mapping->indices;
96132dcc486SBarry Smith   PetscInt       *nprocs,*owner,nsends,*sends,j,*starts,nmax,nrecvs,*recvs,proc;
96297f1f81fSBarry Smith   PetscInt       cnt,scale,*ownedsenders,*nownedsenders,rstart,nowned;
96332dcc486SBarry Smith   PetscInt       node,nownedm,nt,*sends2,nsends2,*starts2,*lens2,*dest,nrecvs2,*starts3,*recvs2,k,*bprocs,*tmp;
96432dcc486SBarry Smith   PetscInt       first_procs,first_numprocs,*first_indices;
96589d82c54SBarry Smith   MPI_Request    *recv_waits,*send_waits;
96630dcb7c9SBarry Smith   MPI_Status     recv_status,*send_status,*recv_statuses;
967ce94432eSBarry Smith   MPI_Comm       comm;
968ace3abfcSBarry Smith   PetscBool      debug = PETSC_FALSE;
96989d82c54SBarry Smith 
97089d82c54SBarry Smith   PetscFunctionBegin;
971ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)mapping,&comm);CHKERRQ(ierr);
97224cf384cSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
97324cf384cSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
97424cf384cSBarry Smith   if (size == 1) {
97524cf384cSBarry Smith     *nproc         = 0;
9760298fd71SBarry Smith     *procs         = NULL;
97795dccacaSBarry Smith     ierr           = PetscNew(numprocs);CHKERRQ(ierr);
9781e2105dcSBarry Smith     (*numprocs)[0] = 0;
97995dccacaSBarry Smith     ierr           = PetscNew(indices);CHKERRQ(ierr);
9800298fd71SBarry Smith     (*indices)[0]  = NULL;
981268a049cSStefano Zampini     /* save info for reuse */
982268a049cSStefano Zampini     mapping->info_nproc = *nproc;
983268a049cSStefano Zampini     mapping->info_procs = *procs;
984268a049cSStefano Zampini     mapping->info_numprocs = *numprocs;
985268a049cSStefano Zampini     mapping->info_indices = *indices;
986268a049cSStefano Zampini     mapping->info_cached = PETSC_TRUE;
98724cf384cSBarry Smith     PetscFunctionReturn(0);
98824cf384cSBarry Smith   }
98924cf384cSBarry Smith 
990c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)mapping)->options,NULL,"-islocaltoglobalmappinggetinfo_debug",&debug,NULL);CHKERRQ(ierr);
99107b52d57SBarry Smith 
9923677ff5aSBarry Smith   /*
9936a818285SBarry Smith     Notes on ISLocalToGlobalMappingGetBlockInfo
9943677ff5aSBarry Smith 
9953677ff5aSBarry Smith     globally owned node - the nodes that have been assigned to this processor in global
9963677ff5aSBarry Smith            numbering, just for this routine.
9973677ff5aSBarry Smith 
9983677ff5aSBarry Smith     nontrivial globally owned node - node assigned to this processor that is on a subdomain
9993677ff5aSBarry Smith            boundary (i.e. is has more than one local owner)
10003677ff5aSBarry Smith 
10013677ff5aSBarry Smith     locally owned node - node that exists on this processors subdomain
10023677ff5aSBarry Smith 
10033677ff5aSBarry Smith     nontrivial locally owned node - node that is not in the interior (i.e. has more than one
10043677ff5aSBarry Smith            local subdomain
10053677ff5aSBarry Smith   */
100624cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag1);CHKERRQ(ierr);
100724cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag2);CHKERRQ(ierr);
100824cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag3);CHKERRQ(ierr);
100989d82c54SBarry Smith 
101089d82c54SBarry Smith   for (i=0; i<n; i++) {
101189d82c54SBarry Smith     if (lindices[i] > max) max = lindices[i];
101289d82c54SBarry Smith   }
1013b2566f29SBarry Smith   ierr   = MPIU_Allreduce(&max,&Ng,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr);
101478058e43SBarry Smith   Ng++;
101589d82c54SBarry Smith   ierr   = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
101689d82c54SBarry Smith   ierr   = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1017bc8ff85bSBarry Smith   scale  = Ng/size + 1;
1018a2e34c3dSBarry Smith   ng     = scale; if (rank == size-1) ng = Ng - scale*(size-1); ng = PetscMax(1,ng);
1019caba0dd0SBarry Smith   rstart = scale*rank;
102089d82c54SBarry Smith 
102189d82c54SBarry Smith   /* determine ownership ranges of global indices */
1022785e854fSJed Brown   ierr = PetscMalloc1(2*size,&nprocs);CHKERRQ(ierr);
102332dcc486SBarry Smith   ierr = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr);
102489d82c54SBarry Smith 
102589d82c54SBarry Smith   /* determine owners of each local node  */
1026785e854fSJed Brown   ierr = PetscMalloc1(n,&owner);CHKERRQ(ierr);
102789d82c54SBarry Smith   for (i=0; i<n; i++) {
10283677ff5aSBarry Smith     proc             = lindices[i]/scale; /* processor that globally owns this index */
102927c402fcSBarry Smith     nprocs[2*proc+1] = 1;                 /* processor globally owns at least one of ours */
10303677ff5aSBarry Smith     owner[i]         = proc;
103127c402fcSBarry Smith     nprocs[2*proc]++;                     /* count of how many that processor globally owns of ours */
103289d82c54SBarry Smith   }
103327c402fcSBarry Smith   nsends = 0; for (i=0; i<size; i++) nsends += nprocs[2*i+1];
10347904a332SBarry Smith   ierr = PetscInfo1(mapping,"Number of global owners for my local data %D\n",nsends);CHKERRQ(ierr);
103589d82c54SBarry Smith 
103689d82c54SBarry Smith   /* inform other processors of number of messages and max length*/
103727c402fcSBarry Smith   ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr);
10387904a332SBarry Smith   ierr = PetscInfo1(mapping,"Number of local owners for my global data %D\n",nrecvs);CHKERRQ(ierr);
103989d82c54SBarry Smith 
104089d82c54SBarry Smith   /* post receives for owned rows */
1041785e854fSJed Brown   ierr = PetscMalloc1((2*nrecvs+1)*(nmax+1),&recvs);CHKERRQ(ierr);
1042854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&recv_waits);CHKERRQ(ierr);
104389d82c54SBarry Smith   for (i=0; i<nrecvs; i++) {
104432dcc486SBarry Smith     ierr = MPI_Irecv(recvs+2*nmax*i,2*nmax,MPIU_INT,MPI_ANY_SOURCE,tag1,comm,recv_waits+i);CHKERRQ(ierr);
104589d82c54SBarry Smith   }
104689d82c54SBarry Smith 
104789d82c54SBarry Smith   /* pack messages containing lists of local nodes to owners */
1048854ce69bSBarry Smith   ierr      = PetscMalloc1(2*n+1,&sends);CHKERRQ(ierr);
1049854ce69bSBarry Smith   ierr      = PetscMalloc1(size+1,&starts);CHKERRQ(ierr);
105089d82c54SBarry Smith   starts[0] = 0;
1051f6e5521dSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
105289d82c54SBarry Smith   for (i=0; i<n; i++) {
105389d82c54SBarry Smith     sends[starts[owner[i]]++] = lindices[i];
105430dcb7c9SBarry Smith     sends[starts[owner[i]]++] = i;
105589d82c54SBarry Smith   }
105689d82c54SBarry Smith   ierr = PetscFree(owner);CHKERRQ(ierr);
105789d82c54SBarry Smith   starts[0] = 0;
1058f6e5521dSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
105989d82c54SBarry Smith 
106089d82c54SBarry Smith   /* send the messages */
1061854ce69bSBarry Smith   ierr = PetscMalloc1(nsends+1,&send_waits);CHKERRQ(ierr);
1062854ce69bSBarry Smith   ierr = PetscMalloc1(nsends+1,&dest);CHKERRQ(ierr);
106389d82c54SBarry Smith   cnt = 0;
106489d82c54SBarry Smith   for (i=0; i<size; i++) {
106527c402fcSBarry Smith     if (nprocs[2*i]) {
106632dcc486SBarry Smith       ierr      = MPI_Isend(sends+starts[i],2*nprocs[2*i],MPIU_INT,i,tag1,comm,send_waits+cnt);CHKERRQ(ierr);
106730dcb7c9SBarry Smith       dest[cnt] = i;
106889d82c54SBarry Smith       cnt++;
106989d82c54SBarry Smith     }
107089d82c54SBarry Smith   }
107189d82c54SBarry Smith   ierr = PetscFree(starts);CHKERRQ(ierr);
107289d82c54SBarry Smith 
107389d82c54SBarry Smith   /* wait on receives */
1074854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&source);CHKERRQ(ierr);
1075854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&len);CHKERRQ(ierr);
107689d82c54SBarry Smith   cnt  = nrecvs;
1077854ce69bSBarry Smith   ierr = PetscMalloc1(ng+1,&nownedsenders);CHKERRQ(ierr);
107832dcc486SBarry Smith   ierr = PetscMemzero(nownedsenders,ng*sizeof(PetscInt));CHKERRQ(ierr);
107989d82c54SBarry Smith   while (cnt) {
108089d82c54SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
108189d82c54SBarry Smith     /* unpack receives into our local space */
108232dcc486SBarry Smith     ierr          = MPI_Get_count(&recv_status,MPIU_INT,&len[imdex]);CHKERRQ(ierr);
108389d82c54SBarry Smith     source[imdex] = recv_status.MPI_SOURCE;
108430dcb7c9SBarry Smith     len[imdex]    = len[imdex]/2;
1085caba0dd0SBarry Smith     /* count how many local owners for each of my global owned indices */
108630dcb7c9SBarry Smith     for (i=0; i<len[imdex]; i++) nownedsenders[recvs[2*imdex*nmax+2*i]-rstart]++;
108789d82c54SBarry Smith     cnt--;
108889d82c54SBarry Smith   }
108989d82c54SBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
109089d82c54SBarry Smith 
109130dcb7c9SBarry Smith   /* count how many globally owned indices are on an edge multiplied by how many processors own them. */
1092bc8ff85bSBarry Smith   nowned  = 0;
1093bc8ff85bSBarry Smith   nownedm = 0;
1094bc8ff85bSBarry Smith   for (i=0; i<ng; i++) {
1095bc8ff85bSBarry Smith     if (nownedsenders[i] > 1) {nownedm += nownedsenders[i]; nowned++;}
1096bc8ff85bSBarry Smith   }
1097bc8ff85bSBarry Smith 
1098bc8ff85bSBarry Smith   /* create single array to contain rank of all local owners of each globally owned index */
1099854ce69bSBarry Smith   ierr      = PetscMalloc1(nownedm+1,&ownedsenders);CHKERRQ(ierr);
1100854ce69bSBarry Smith   ierr      = PetscMalloc1(ng+1,&starts);CHKERRQ(ierr);
1101bc8ff85bSBarry Smith   starts[0] = 0;
1102bc8ff85bSBarry Smith   for (i=1; i<ng; i++) {
1103bc8ff85bSBarry Smith     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
1104bc8ff85bSBarry Smith     else starts[i] = starts[i-1];
1105bc8ff85bSBarry Smith   }
1106bc8ff85bSBarry Smith 
110730dcb7c9SBarry Smith   /* for each nontrival globally owned node list all arriving processors */
1108bc8ff85bSBarry Smith   for (i=0; i<nrecvs; i++) {
1109bc8ff85bSBarry Smith     for (j=0; j<len[i]; j++) {
111030dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
1111f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) ownedsenders[starts[node]++] = source[i];
1112bc8ff85bSBarry Smith     }
1113bc8ff85bSBarry Smith   }
1114bc8ff85bSBarry Smith 
111507b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
111630dcb7c9SBarry Smith     starts[0] = 0;
111730dcb7c9SBarry Smith     for (i=1; i<ng; i++) {
111830dcb7c9SBarry Smith       if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
111930dcb7c9SBarry Smith       else starts[i] = starts[i-1];
112030dcb7c9SBarry Smith     }
112130dcb7c9SBarry Smith     for (i=0; i<ng; i++) {
112230dcb7c9SBarry Smith       if (nownedsenders[i] > 1) {
11237904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"[%d] global node %D local owner processors: ",rank,i+rstart);CHKERRQ(ierr);
112430dcb7c9SBarry Smith         for (j=0; j<nownedsenders[i]; j++) {
11257904a332SBarry Smith           ierr = PetscSynchronizedPrintf(comm,"%D ",ownedsenders[starts[i]+j]);CHKERRQ(ierr);
112630dcb7c9SBarry Smith         }
112730dcb7c9SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
112830dcb7c9SBarry Smith       }
112930dcb7c9SBarry Smith     }
11300ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
113107b52d57SBarry Smith   } /* -----------------------------------  */
113230dcb7c9SBarry Smith 
11333677ff5aSBarry Smith   /* wait on original sends */
11343a96401aSBarry Smith   if (nsends) {
1135785e854fSJed Brown     ierr = PetscMalloc1(nsends,&send_status);CHKERRQ(ierr);
11363a96401aSBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
11373a96401aSBarry Smith     ierr = PetscFree(send_status);CHKERRQ(ierr);
11383a96401aSBarry Smith   }
113989d82c54SBarry Smith   ierr = PetscFree(send_waits);CHKERRQ(ierr);
11403a96401aSBarry Smith   ierr = PetscFree(sends);CHKERRQ(ierr);
11413677ff5aSBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
11423677ff5aSBarry Smith 
11433677ff5aSBarry Smith   /* pack messages to send back to local owners */
114430dcb7c9SBarry Smith   starts[0] = 0;
114530dcb7c9SBarry Smith   for (i=1; i<ng; i++) {
114630dcb7c9SBarry Smith     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
114730dcb7c9SBarry Smith     else starts[i] = starts[i-1];
114830dcb7c9SBarry Smith   }
114930dcb7c9SBarry Smith   nsends2 = nrecvs;
1150854ce69bSBarry Smith   ierr    = PetscMalloc1(nsends2+1,&nprocs);CHKERRQ(ierr); /* length of each message */
115130dcb7c9SBarry Smith   for (i=0; i<nrecvs; i++) {
115230dcb7c9SBarry Smith     nprocs[i] = 1;
115330dcb7c9SBarry Smith     for (j=0; j<len[i]; j++) {
115430dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
1155f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) nprocs[i] += 2 + nownedsenders[node];
115630dcb7c9SBarry Smith     }
115730dcb7c9SBarry Smith   }
1158f6e5521dSKarl Rupp   nt = 0;
1159f6e5521dSKarl Rupp   for (i=0; i<nsends2; i++) nt += nprocs[i];
1160f6e5521dSKarl Rupp 
1161854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,&sends2);CHKERRQ(ierr);
1162854ce69bSBarry Smith   ierr = PetscMalloc1(nsends2+1,&starts2);CHKERRQ(ierr);
1163f6e5521dSKarl Rupp 
1164f6e5521dSKarl Rupp   starts2[0] = 0;
1165f6e5521dSKarl Rupp   for (i=1; i<nsends2; i++) starts2[i] = starts2[i-1] + nprocs[i-1];
116630dcb7c9SBarry Smith   /*
116730dcb7c9SBarry Smith      Each message is 1 + nprocs[i] long, and consists of
116830dcb7c9SBarry Smith        (0) the number of nodes being sent back
116930dcb7c9SBarry Smith        (1) the local node number,
117030dcb7c9SBarry Smith        (2) the number of processors sharing it,
117130dcb7c9SBarry Smith        (3) the processors sharing it
117230dcb7c9SBarry Smith   */
117330dcb7c9SBarry Smith   for (i=0; i<nsends2; i++) {
117430dcb7c9SBarry Smith     cnt = 1;
117530dcb7c9SBarry Smith     sends2[starts2[i]] = 0;
117630dcb7c9SBarry Smith     for (j=0; j<len[i]; j++) {
117730dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
117830dcb7c9SBarry Smith       if (nownedsenders[node] > 1) {
117930dcb7c9SBarry Smith         sends2[starts2[i]]++;
118030dcb7c9SBarry Smith         sends2[starts2[i]+cnt++] = recvs[2*i*nmax+2*j+1];
118130dcb7c9SBarry Smith         sends2[starts2[i]+cnt++] = nownedsenders[node];
118232dcc486SBarry Smith         ierr = PetscMemcpy(&sends2[starts2[i]+cnt],&ownedsenders[starts[node]],nownedsenders[node]*sizeof(PetscInt));CHKERRQ(ierr);
118330dcb7c9SBarry Smith         cnt += nownedsenders[node];
118430dcb7c9SBarry Smith       }
118530dcb7c9SBarry Smith     }
118630dcb7c9SBarry Smith   }
118730dcb7c9SBarry Smith 
118830dcb7c9SBarry Smith   /* receive the message lengths */
118930dcb7c9SBarry Smith   nrecvs2 = nsends;
1190854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&lens2);CHKERRQ(ierr);
1191854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&starts3);CHKERRQ(ierr);
1192854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&recv_waits);CHKERRQ(ierr);
119330dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
1194d44834fbSBarry Smith     ierr = MPI_Irecv(&lens2[i],1,MPIU_INT,dest[i],tag2,comm,recv_waits+i);CHKERRQ(ierr);
119530dcb7c9SBarry Smith   }
1196d44834fbSBarry Smith 
11978a8e0b3aSBarry Smith   /* send the message lengths */
11988a8e0b3aSBarry Smith   for (i=0; i<nsends2; i++) {
11998a8e0b3aSBarry Smith     ierr = MPI_Send(&nprocs[i],1,MPIU_INT,source[i],tag2,comm);CHKERRQ(ierr);
12008a8e0b3aSBarry Smith   }
12018a8e0b3aSBarry Smith 
1202d44834fbSBarry Smith   /* wait on receives of lens */
12030c468ba9SBarry Smith   if (nrecvs2) {
1204785e854fSJed Brown     ierr = PetscMalloc1(nrecvs2,&recv_statuses);CHKERRQ(ierr);
1205d44834fbSBarry Smith     ierr = MPI_Waitall(nrecvs2,recv_waits,recv_statuses);CHKERRQ(ierr);
1206d44834fbSBarry Smith     ierr = PetscFree(recv_statuses);CHKERRQ(ierr);
12070c468ba9SBarry Smith   }
1208a2ea699eSBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
1209d44834fbSBarry Smith 
121030dcb7c9SBarry Smith   starts3[0] = 0;
1211d44834fbSBarry Smith   nt         = 0;
121230dcb7c9SBarry Smith   for (i=0; i<nrecvs2-1; i++) {
121330dcb7c9SBarry Smith     starts3[i+1] = starts3[i] + lens2[i];
1214d44834fbSBarry Smith     nt          += lens2[i];
121530dcb7c9SBarry Smith   }
121676466f69SStefano Zampini   if (nrecvs2) nt += lens2[nrecvs2-1];
1217d44834fbSBarry Smith 
1218854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,&recvs2);CHKERRQ(ierr);
1219854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs2+1,&recv_waits);CHKERRQ(ierr);
122052b72c4aSBarry Smith   for (i=0; i<nrecvs2; i++) {
122132dcc486SBarry Smith     ierr = MPI_Irecv(recvs2+starts3[i],lens2[i],MPIU_INT,dest[i],tag3,comm,recv_waits+i);CHKERRQ(ierr);
122230dcb7c9SBarry Smith   }
122330dcb7c9SBarry Smith 
122430dcb7c9SBarry Smith   /* send the messages */
1225854ce69bSBarry Smith   ierr = PetscMalloc1(nsends2+1,&send_waits);CHKERRQ(ierr);
122630dcb7c9SBarry Smith   for (i=0; i<nsends2; i++) {
122732dcc486SBarry Smith     ierr = MPI_Isend(sends2+starts2[i],nprocs[i],MPIU_INT,source[i],tag3,comm,send_waits+i);CHKERRQ(ierr);
122830dcb7c9SBarry Smith   }
122930dcb7c9SBarry Smith 
123030dcb7c9SBarry Smith   /* wait on receives */
12310c468ba9SBarry Smith   if (nrecvs2) {
1232785e854fSJed Brown     ierr = PetscMalloc1(nrecvs2,&recv_statuses);CHKERRQ(ierr);
123330dcb7c9SBarry Smith     ierr = MPI_Waitall(nrecvs2,recv_waits,recv_statuses);CHKERRQ(ierr);
123430dcb7c9SBarry Smith     ierr = PetscFree(recv_statuses);CHKERRQ(ierr);
12350c468ba9SBarry Smith   }
123630dcb7c9SBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
123730dcb7c9SBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
123830dcb7c9SBarry Smith 
123907b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
124030dcb7c9SBarry Smith     cnt = 0;
124130dcb7c9SBarry Smith     for (i=0; i<nrecvs2; i++) {
124230dcb7c9SBarry Smith       nt = recvs2[cnt++];
124330dcb7c9SBarry Smith       for (j=0; j<nt; j++) {
12447904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"[%d] local node %D number of subdomains %D: ",rank,recvs2[cnt],recvs2[cnt+1]);CHKERRQ(ierr);
124530dcb7c9SBarry Smith         for (k=0; k<recvs2[cnt+1]; k++) {
12467904a332SBarry Smith           ierr = PetscSynchronizedPrintf(comm,"%D ",recvs2[cnt+2+k]);CHKERRQ(ierr);
124730dcb7c9SBarry Smith         }
124830dcb7c9SBarry Smith         cnt += 2 + recvs2[cnt+1];
124930dcb7c9SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
125030dcb7c9SBarry Smith       }
125130dcb7c9SBarry Smith     }
12520ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
125307b52d57SBarry Smith   } /* -----------------------------------  */
125430dcb7c9SBarry Smith 
125530dcb7c9SBarry Smith   /* count number subdomains for each local node */
1256785e854fSJed Brown   ierr = PetscMalloc1(size,&nprocs);CHKERRQ(ierr);
125732dcc486SBarry Smith   ierr = PetscMemzero(nprocs,size*sizeof(PetscInt));CHKERRQ(ierr);
125830dcb7c9SBarry Smith   cnt  = 0;
125930dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
126030dcb7c9SBarry Smith     nt = recvs2[cnt++];
126130dcb7c9SBarry Smith     for (j=0; j<nt; j++) {
1262f6e5521dSKarl Rupp       for (k=0; k<recvs2[cnt+1]; k++) nprocs[recvs2[cnt+2+k]]++;
126330dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt+1];
126430dcb7c9SBarry Smith     }
126530dcb7c9SBarry Smith   }
126630dcb7c9SBarry Smith   nt = 0; for (i=0; i<size; i++) nt += (nprocs[i] > 0);
126730dcb7c9SBarry Smith   *nproc    = nt;
1268854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,procs);CHKERRQ(ierr);
1269854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,numprocs);CHKERRQ(ierr);
1270854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,indices);CHKERRQ(ierr);
12710298fd71SBarry Smith   for (i=0;i<nt+1;i++) (*indices)[i]=NULL;
1272785e854fSJed Brown   ierr = PetscMalloc1(size,&bprocs);CHKERRQ(ierr);
127330dcb7c9SBarry Smith   cnt  = 0;
127430dcb7c9SBarry Smith   for (i=0; i<size; i++) {
127530dcb7c9SBarry Smith     if (nprocs[i] > 0) {
127630dcb7c9SBarry Smith       bprocs[i]        = cnt;
127730dcb7c9SBarry Smith       (*procs)[cnt]    = i;
127830dcb7c9SBarry Smith       (*numprocs)[cnt] = nprocs[i];
1279785e854fSJed Brown       ierr             = PetscMalloc1(nprocs[i],&(*indices)[cnt]);CHKERRQ(ierr);
128030dcb7c9SBarry Smith       cnt++;
128130dcb7c9SBarry Smith     }
128230dcb7c9SBarry Smith   }
128330dcb7c9SBarry Smith 
128430dcb7c9SBarry Smith   /* make the list of subdomains for each nontrivial local node */
128532dcc486SBarry Smith   ierr = PetscMemzero(*numprocs,nt*sizeof(PetscInt));CHKERRQ(ierr);
128630dcb7c9SBarry Smith   cnt  = 0;
128730dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
128830dcb7c9SBarry Smith     nt = recvs2[cnt++];
128930dcb7c9SBarry Smith     for (j=0; j<nt; j++) {
1290f6e5521dSKarl Rupp       for (k=0; k<recvs2[cnt+1]; k++) (*indices)[bprocs[recvs2[cnt+2+k]]][(*numprocs)[bprocs[recvs2[cnt+2+k]]]++] = recvs2[cnt];
129130dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt+1];
129230dcb7c9SBarry Smith     }
129330dcb7c9SBarry Smith   }
129430dcb7c9SBarry Smith   ierr = PetscFree(bprocs);CHKERRQ(ierr);
129507b52d57SBarry Smith   ierr = PetscFree(recvs2);CHKERRQ(ierr);
129630dcb7c9SBarry Smith 
129707b52d57SBarry Smith   /* sort the node indexing by their global numbers */
129807b52d57SBarry Smith   nt = *nproc;
129907b52d57SBarry Smith   for (i=0; i<nt; i++) {
1300854ce69bSBarry Smith     ierr = PetscMalloc1((*numprocs)[i],&tmp);CHKERRQ(ierr);
1301f6e5521dSKarl Rupp     for (j=0; j<(*numprocs)[i]; j++) tmp[j] = lindices[(*indices)[i][j]];
130207b52d57SBarry Smith     ierr = PetscSortIntWithArray((*numprocs)[i],tmp,(*indices)[i]);CHKERRQ(ierr);
130307b52d57SBarry Smith     ierr = PetscFree(tmp);CHKERRQ(ierr);
130407b52d57SBarry Smith   }
130507b52d57SBarry Smith 
130607b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
130730dcb7c9SBarry Smith     nt = *nproc;
130830dcb7c9SBarry Smith     for (i=0; i<nt; i++) {
13097904a332SBarry Smith       ierr = PetscSynchronizedPrintf(comm,"[%d] subdomain %D number of indices %D: ",rank,(*procs)[i],(*numprocs)[i]);CHKERRQ(ierr);
131030dcb7c9SBarry Smith       for (j=0; j<(*numprocs)[i]; j++) {
13117904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"%D ",(*indices)[i][j]);CHKERRQ(ierr);
131230dcb7c9SBarry Smith       }
131330dcb7c9SBarry Smith       ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
131430dcb7c9SBarry Smith     }
13150ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
131607b52d57SBarry Smith   } /* -----------------------------------  */
131730dcb7c9SBarry Smith 
131830dcb7c9SBarry Smith   /* wait on sends */
131930dcb7c9SBarry Smith   if (nsends2) {
1320785e854fSJed Brown     ierr = PetscMalloc1(nsends2,&send_status);CHKERRQ(ierr);
132130dcb7c9SBarry Smith     ierr = MPI_Waitall(nsends2,send_waits,send_status);CHKERRQ(ierr);
132230dcb7c9SBarry Smith     ierr = PetscFree(send_status);CHKERRQ(ierr);
132330dcb7c9SBarry Smith   }
132430dcb7c9SBarry Smith 
132530dcb7c9SBarry Smith   ierr = PetscFree(starts3);CHKERRQ(ierr);
132630dcb7c9SBarry Smith   ierr = PetscFree(dest);CHKERRQ(ierr);
132730dcb7c9SBarry Smith   ierr = PetscFree(send_waits);CHKERRQ(ierr);
13283677ff5aSBarry Smith 
1329bc8ff85bSBarry Smith   ierr = PetscFree(nownedsenders);CHKERRQ(ierr);
1330bc8ff85bSBarry Smith   ierr = PetscFree(ownedsenders);CHKERRQ(ierr);
1331bc8ff85bSBarry Smith   ierr = PetscFree(starts);CHKERRQ(ierr);
133230dcb7c9SBarry Smith   ierr = PetscFree(starts2);CHKERRQ(ierr);
133330dcb7c9SBarry Smith   ierr = PetscFree(lens2);CHKERRQ(ierr);
133489d82c54SBarry Smith 
133589d82c54SBarry Smith   ierr = PetscFree(source);CHKERRQ(ierr);
133697f1f81fSBarry Smith   ierr = PetscFree(len);CHKERRQ(ierr);
133789d82c54SBarry Smith   ierr = PetscFree(recvs);CHKERRQ(ierr);
13383a96401aSBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
133930dcb7c9SBarry Smith   ierr = PetscFree(sends2);CHKERRQ(ierr);
134024cf384cSBarry Smith 
134124cf384cSBarry Smith   /* put the information about myself as the first entry in the list */
134224cf384cSBarry Smith   first_procs    = (*procs)[0];
134324cf384cSBarry Smith   first_numprocs = (*numprocs)[0];
134424cf384cSBarry Smith   first_indices  = (*indices)[0];
134524cf384cSBarry Smith   for (i=0; i<*nproc; i++) {
134624cf384cSBarry Smith     if ((*procs)[i] == rank) {
134724cf384cSBarry Smith       (*procs)[0]    = (*procs)[i];
134824cf384cSBarry Smith       (*numprocs)[0] = (*numprocs)[i];
134924cf384cSBarry Smith       (*indices)[0]  = (*indices)[i];
135024cf384cSBarry Smith       (*procs)[i]    = first_procs;
135124cf384cSBarry Smith       (*numprocs)[i] = first_numprocs;
135224cf384cSBarry Smith       (*indices)[i]  = first_indices;
135324cf384cSBarry Smith       break;
135424cf384cSBarry Smith     }
135524cf384cSBarry Smith   }
1356268a049cSStefano Zampini 
1357268a049cSStefano Zampini   /* save info for reuse */
1358268a049cSStefano Zampini   mapping->info_nproc = *nproc;
1359268a049cSStefano Zampini   mapping->info_procs = *procs;
1360268a049cSStefano Zampini   mapping->info_numprocs = *numprocs;
1361268a049cSStefano Zampini   mapping->info_indices = *indices;
1362268a049cSStefano Zampini   mapping->info_cached = PETSC_TRUE;
136389d82c54SBarry Smith   PetscFunctionReturn(0);
136489d82c54SBarry Smith }
136589d82c54SBarry Smith 
13666a818285SBarry Smith /*@C
13676a818285SBarry Smith     ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by ISLocalToGlobalMappingGetBlockInfo()
13686a818285SBarry Smith 
13696a818285SBarry Smith     Collective on ISLocalToGlobalMapping
13706a818285SBarry Smith 
13716a818285SBarry Smith     Input Parameters:
13726a818285SBarry Smith .   mapping - the mapping from local to global indexing
13736a818285SBarry Smith 
13746a818285SBarry Smith     Output Parameter:
13756a818285SBarry Smith +   nproc - number of processors that are connected to this one
13766a818285SBarry Smith .   proc - neighboring processors
13776a818285SBarry Smith .   numproc - number of indices for each processor
13786a818285SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
13796a818285SBarry Smith 
13806a818285SBarry Smith     Level: advanced
13816a818285SBarry Smith 
13826a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
13836a818285SBarry Smith           ISLocalToGlobalMappingGetInfo()
13846a818285SBarry Smith @*/
13856a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
13866a818285SBarry Smith {
13876a818285SBarry Smith   PetscErrorCode ierr;
13886a818285SBarry Smith 
13896a818285SBarry Smith   PetscFunctionBegin;
1390cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1391268a049cSStefano Zampini   if (mapping->info_free) {
13926a818285SBarry Smith     ierr = PetscFree(*numprocs);CHKERRQ(ierr);
13936a818285SBarry Smith     if (*indices) {
1394268a049cSStefano Zampini       PetscInt i;
1395268a049cSStefano Zampini 
13966a818285SBarry Smith       ierr = PetscFree((*indices)[0]);CHKERRQ(ierr);
13976a818285SBarry Smith       for (i=1; i<*nproc; i++) {
13986a818285SBarry Smith         ierr = PetscFree((*indices)[i]);CHKERRQ(ierr);
13996a818285SBarry Smith       }
14006a818285SBarry Smith       ierr = PetscFree(*indices);CHKERRQ(ierr);
14016a818285SBarry Smith     }
1402268a049cSStefano Zampini   }
1403268a049cSStefano Zampini   *nproc = 0;
1404268a049cSStefano Zampini   *procs = NULL;
1405268a049cSStefano Zampini   *numprocs = NULL;
1406268a049cSStefano Zampini   *indices = NULL;
14076a818285SBarry Smith   PetscFunctionReturn(0);
14086a818285SBarry Smith }
14096a818285SBarry Smith 
14106a818285SBarry Smith /*@C
14116a818285SBarry Smith     ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each processor and
14126a818285SBarry Smith      each index shared by more than one processor
14136a818285SBarry Smith 
14146a818285SBarry Smith     Collective on ISLocalToGlobalMapping
14156a818285SBarry Smith 
14166a818285SBarry Smith     Input Parameters:
14176a818285SBarry Smith .   mapping - the mapping from local to global indexing
14186a818285SBarry Smith 
14196a818285SBarry Smith     Output Parameter:
14206a818285SBarry Smith +   nproc - number of processors that are connected to this one
14216a818285SBarry Smith .   proc - neighboring processors
14226a818285SBarry Smith .   numproc - number of indices for each subdomain (processor)
14236a818285SBarry Smith -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
14246a818285SBarry Smith 
14256a818285SBarry Smith     Level: advanced
14266a818285SBarry Smith 
14276a818285SBarry Smith     Concepts: mapping^local to global
14286a818285SBarry Smith 
14296a818285SBarry Smith     Fortran Usage:
14306a818285SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
14316a818285SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
14326a818285SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
14336a818285SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
14346a818285SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
14356a818285SBarry Smith 
14366a818285SBarry Smith 
14376a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
14386a818285SBarry Smith           ISLocalToGlobalMappingRestoreInfo()
14396a818285SBarry Smith @*/
14406a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
14416a818285SBarry Smith {
14426a818285SBarry Smith   PetscErrorCode ierr;
1443268a049cSStefano Zampini   PetscInt       **bindices = NULL,*bnumprocs = NULL,bs = mapping->bs,i,j,k;
14446a818285SBarry Smith 
14456a818285SBarry Smith   PetscFunctionBegin;
14466a818285SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1447268a049cSStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockInfo(mapping,nproc,procs,&bnumprocs,&bindices);CHKERRQ(ierr);
1448268a049cSStefano Zampini   if (bs > 1) { /* we need to expand the cached info */
1449732f65e3SBarry Smith     ierr = PetscCalloc1(*nproc,&*indices);CHKERRQ(ierr);
1450268a049cSStefano Zampini     ierr = PetscCalloc1(*nproc,&*numprocs);CHKERRQ(ierr);
14516a818285SBarry Smith     for (i=0; i<*nproc; i++) {
1452268a049cSStefano Zampini       ierr = PetscMalloc1(bs*bnumprocs[i],&(*indices)[i]);CHKERRQ(ierr);
1453268a049cSStefano Zampini       for (j=0; j<bnumprocs[i]; j++) {
14546a818285SBarry Smith         for (k=0; k<bs; k++) {
14556a818285SBarry Smith           (*indices)[i][j*bs+k] = bs*bindices[i][j] + k;
14566a818285SBarry Smith         }
14576a818285SBarry Smith       }
1458268a049cSStefano Zampini       (*numprocs)[i] = bnumprocs[i]*bs;
14596a818285SBarry Smith     }
1460268a049cSStefano Zampini     mapping->info_free = PETSC_TRUE;
1461268a049cSStefano Zampini   } else {
1462268a049cSStefano Zampini     *numprocs = bnumprocs;
1463268a049cSStefano Zampini     *indices  = bindices;
14646a818285SBarry Smith   }
14656a818285SBarry Smith   PetscFunctionReturn(0);
14666a818285SBarry Smith }
14676a818285SBarry Smith 
146807b52d57SBarry Smith /*@C
146907b52d57SBarry Smith     ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by ISLocalToGlobalMappingGetInfo()
147089d82c54SBarry Smith 
147107b52d57SBarry Smith     Collective on ISLocalToGlobalMapping
147207b52d57SBarry Smith 
147307b52d57SBarry Smith     Input Parameters:
147407b52d57SBarry Smith .   mapping - the mapping from local to global indexing
147507b52d57SBarry Smith 
147607b52d57SBarry Smith     Output Parameter:
147707b52d57SBarry Smith +   nproc - number of processors that are connected to this one
147807b52d57SBarry Smith .   proc - neighboring processors
147907b52d57SBarry Smith .   numproc - number of indices for each processor
148007b52d57SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
148107b52d57SBarry Smith 
148207b52d57SBarry Smith     Level: advanced
148307b52d57SBarry Smith 
148407b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
148507b52d57SBarry Smith           ISLocalToGlobalMappingGetInfo()
148607b52d57SBarry Smith @*/
14877087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
148807b52d57SBarry Smith {
14896849ba73SBarry Smith   PetscErrorCode ierr;
149007b52d57SBarry Smith 
149107b52d57SBarry Smith   PetscFunctionBegin;
14926a818285SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockInfo(mapping,nproc,procs,numprocs,indices);CHKERRQ(ierr);
149307b52d57SBarry Smith   PetscFunctionReturn(0);
149407b52d57SBarry Smith }
149586994e45SJed Brown 
149686994e45SJed Brown /*@C
1497107e9a97SBarry Smith    ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped
149886994e45SJed Brown 
149986994e45SJed Brown    Not Collective
150086994e45SJed Brown 
150186994e45SJed Brown    Input Arguments:
150286994e45SJed Brown . ltog - local to global mapping
150386994e45SJed Brown 
150486994e45SJed Brown    Output Arguments:
1505565245c5SBarry Smith . array - array of indices, the length of this array may be obtained with ISLocalToGlobalMappingGetSize()
150686994e45SJed Brown 
150786994e45SJed Brown    Level: advanced
150886994e45SJed Brown 
1509107e9a97SBarry Smith    Notes: ISLocalToGlobalMappingGetSize() returns the length the this array
1510107e9a97SBarry Smith 
1511107e9a97SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreIndices(), ISLocalToGlobalMappingGetBlockIndices(), ISLocalToGlobalMappingRestoreBlockIndices()
151286994e45SJed Brown @*/
15137087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
151486994e45SJed Brown {
151586994e45SJed Brown   PetscFunctionBegin;
151686994e45SJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
151786994e45SJed Brown   PetscValidPointer(array,2);
151845b6f7e9SBarry Smith   if (ltog->bs == 1) {
151986994e45SJed Brown     *array = ltog->indices;
152045b6f7e9SBarry Smith   } else {
152145b6f7e9SBarry Smith     PetscInt       *jj,k,i,j,n = ltog->n, bs = ltog->bs;
152245b6f7e9SBarry Smith     const PetscInt *ii;
152345b6f7e9SBarry Smith     PetscErrorCode ierr;
152445b6f7e9SBarry Smith 
152545b6f7e9SBarry Smith     ierr = PetscMalloc1(bs*n,&jj);CHKERRQ(ierr);
152645b6f7e9SBarry Smith     *array = jj;
152745b6f7e9SBarry Smith     k    = 0;
152845b6f7e9SBarry Smith     ii   = ltog->indices;
152945b6f7e9SBarry Smith     for (i=0; i<n; i++)
153045b6f7e9SBarry Smith       for (j=0; j<bs; j++)
153145b6f7e9SBarry Smith         jj[k++] = bs*ii[i] + j;
153245b6f7e9SBarry Smith   }
153386994e45SJed Brown   PetscFunctionReturn(0);
153486994e45SJed Brown }
153586994e45SJed Brown 
153686994e45SJed Brown /*@C
1537193a2b41SJulian Andrej    ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with ISLocalToGlobalMappingGetIndices()
153886994e45SJed Brown 
153986994e45SJed Brown    Not Collective
154086994e45SJed Brown 
154186994e45SJed Brown    Input Arguments:
154286994e45SJed Brown + ltog - local to global mapping
154386994e45SJed Brown - array - array of indices
154486994e45SJed Brown 
154586994e45SJed Brown    Level: advanced
154686994e45SJed Brown 
154786994e45SJed Brown .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices()
154886994e45SJed Brown @*/
15497087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
155086994e45SJed Brown {
155186994e45SJed Brown   PetscFunctionBegin;
155286994e45SJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
155386994e45SJed Brown   PetscValidPointer(array,2);
155445b6f7e9SBarry Smith   if (ltog->bs == 1 && *array != ltog->indices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
155545b6f7e9SBarry Smith 
155645b6f7e9SBarry Smith   if (ltog->bs > 1) {
155745b6f7e9SBarry Smith     PetscErrorCode ierr;
155845b6f7e9SBarry Smith     ierr = PetscFree(*(void**)array);CHKERRQ(ierr);
155945b6f7e9SBarry Smith   }
156045b6f7e9SBarry Smith   PetscFunctionReturn(0);
156145b6f7e9SBarry Smith }
156245b6f7e9SBarry Smith 
156345b6f7e9SBarry Smith /*@C
156445b6f7e9SBarry Smith    ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block
156545b6f7e9SBarry Smith 
156645b6f7e9SBarry Smith    Not Collective
156745b6f7e9SBarry Smith 
156845b6f7e9SBarry Smith    Input Arguments:
156945b6f7e9SBarry Smith . ltog - local to global mapping
157045b6f7e9SBarry Smith 
157145b6f7e9SBarry Smith    Output Arguments:
157245b6f7e9SBarry Smith . array - array of indices
157345b6f7e9SBarry Smith 
157445b6f7e9SBarry Smith    Level: advanced
157545b6f7e9SBarry Smith 
157645b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreBlockIndices()
157745b6f7e9SBarry Smith @*/
157845b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
157945b6f7e9SBarry Smith {
158045b6f7e9SBarry Smith   PetscFunctionBegin;
158145b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
158245b6f7e9SBarry Smith   PetscValidPointer(array,2);
158345b6f7e9SBarry Smith   *array = ltog->indices;
158445b6f7e9SBarry Smith   PetscFunctionReturn(0);
158545b6f7e9SBarry Smith }
158645b6f7e9SBarry Smith 
158745b6f7e9SBarry Smith /*@C
158845b6f7e9SBarry Smith    ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with ISLocalToGlobalMappingGetBlockIndices()
158945b6f7e9SBarry Smith 
159045b6f7e9SBarry Smith    Not Collective
159145b6f7e9SBarry Smith 
159245b6f7e9SBarry Smith    Input Arguments:
159345b6f7e9SBarry Smith + ltog - local to global mapping
159445b6f7e9SBarry Smith - array - array of indices
159545b6f7e9SBarry Smith 
159645b6f7e9SBarry Smith    Level: advanced
159745b6f7e9SBarry Smith 
159845b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices()
159945b6f7e9SBarry Smith @*/
160045b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
160145b6f7e9SBarry Smith {
160245b6f7e9SBarry Smith   PetscFunctionBegin;
160345b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
160445b6f7e9SBarry Smith   PetscValidPointer(array,2);
160586994e45SJed Brown   if (*array != ltog->indices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
16060298fd71SBarry Smith   *array = NULL;
160786994e45SJed Brown   PetscFunctionReturn(0);
160886994e45SJed Brown }
1609f7efa3c7SJed Brown 
1610f7efa3c7SJed Brown /*@C
1611f7efa3c7SJed Brown    ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings
1612f7efa3c7SJed Brown 
1613f7efa3c7SJed Brown    Not Collective
1614f7efa3c7SJed Brown 
1615f7efa3c7SJed Brown    Input Arguments:
1616f7efa3c7SJed Brown + comm - communicator for the new mapping, must contain the communicator of every mapping to concatenate
1617f7efa3c7SJed Brown . n - number of mappings to concatenate
1618f7efa3c7SJed Brown - ltogs - local to global mappings
1619f7efa3c7SJed Brown 
1620f7efa3c7SJed Brown    Output Arguments:
1621f7efa3c7SJed Brown . ltogcat - new mapping
1622f7efa3c7SJed Brown 
16239d90f715SBarry Smith    Note: this currently always returns a mapping with block size of 1
16249d90f715SBarry Smith 
16259d90f715SBarry Smith    Developer Note: If all the input mapping have the same block size we could easily handle that as a special case
16269d90f715SBarry Smith 
1627f7efa3c7SJed Brown    Level: advanced
1628f7efa3c7SJed Brown 
1629f7efa3c7SJed Brown .seealso: ISLocalToGlobalMappingCreate()
1630f7efa3c7SJed Brown @*/
1631f7efa3c7SJed Brown PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm,PetscInt n,const ISLocalToGlobalMapping ltogs[],ISLocalToGlobalMapping *ltogcat)
1632f7efa3c7SJed Brown {
1633f7efa3c7SJed Brown   PetscInt       i,cnt,m,*idx;
1634f7efa3c7SJed Brown   PetscErrorCode ierr;
1635f7efa3c7SJed Brown 
1636f7efa3c7SJed Brown   PetscFunctionBegin;
1637f7efa3c7SJed Brown   if (n < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have a non-negative number of mappings, given %D",n);
1638f7efa3c7SJed Brown   if (n > 0) PetscValidPointer(ltogs,3);
1639f7efa3c7SJed Brown   for (i=0; i<n; i++) PetscValidHeaderSpecific(ltogs[i],IS_LTOGM_CLASSID,3);
1640f7efa3c7SJed Brown   PetscValidPointer(ltogcat,4);
1641f7efa3c7SJed Brown   for (cnt=0,i=0; i<n; i++) {
1642f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
1643f7efa3c7SJed Brown     cnt += m;
1644f7efa3c7SJed Brown   }
1645785e854fSJed Brown   ierr = PetscMalloc1(cnt,&idx);CHKERRQ(ierr);
1646f7efa3c7SJed Brown   for (cnt=0,i=0; i<n; i++) {
1647f7efa3c7SJed Brown     const PetscInt *subidx;
1648f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
1649f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetIndices(ltogs[i],&subidx);CHKERRQ(ierr);
1650f7efa3c7SJed Brown     ierr = PetscMemcpy(&idx[cnt],subidx,m*sizeof(PetscInt));CHKERRQ(ierr);
1651f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingRestoreIndices(ltogs[i],&subidx);CHKERRQ(ierr);
1652f7efa3c7SJed Brown     cnt += m;
1653f7efa3c7SJed Brown   }
1654f0413b6fSBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,1,cnt,idx,PETSC_OWN_POINTER,ltogcat);CHKERRQ(ierr);
1655f7efa3c7SJed Brown   PetscFunctionReturn(0);
1656f7efa3c7SJed Brown }
165704a59952SBarry Smith 
165804a59952SBarry Smith 
1659