xref: /petsc/src/vec/is/utils/isltog.c (revision bbf7bc21ccc5b8343ec17e15d9927257b0fe4856)
12362add9SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/isimpl.h>    /*I "petscis.h"  I*/
3*bbf7bc21SLisandro Dalcin #include <petsc/private/hash.h>
40c312b8eSJed Brown #include <petscsf.h>
5665c2dedSJed Brown #include <petscviewer.h>
62362add9SBarry Smith 
77087cfbeSBarry Smith PetscClassId IS_LTOGM_CLASSID;
8268a049cSStefano Zampini static PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping,PetscInt*,PetscInt**,PetscInt**,PetscInt***);
98e58c17dSMatthew Knepley 
10413f72f0SBarry Smith typedef struct {
11413f72f0SBarry Smith   PetscInt    *globals;
12413f72f0SBarry Smith } ISLocalToGlobalMapping_Basic;
13413f72f0SBarry Smith 
14413f72f0SBarry Smith typedef struct {
15413f72f0SBarry Smith   PetscHashI globalht;
16413f72f0SBarry Smith } ISLocalToGlobalMapping_Hash;
17413f72f0SBarry Smith 
18413f72f0SBarry Smith 
19413f72f0SBarry Smith /* -----------------------------------------------------------------------------------------*/
20413f72f0SBarry Smith 
21413f72f0SBarry Smith /*
22413f72f0SBarry Smith     Creates the global mapping information in the ISLocalToGlobalMapping structure
23413f72f0SBarry Smith 
24413f72f0SBarry Smith     If the user has not selected how to handle the global to local mapping then use HASH for "large" problems
25413f72f0SBarry Smith */
26413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp(ISLocalToGlobalMapping mapping)
27413f72f0SBarry Smith {
28413f72f0SBarry Smith   PetscInt       i,*idx = mapping->indices,n = mapping->n,end,start;
29413f72f0SBarry Smith   PetscErrorCode ierr;
30413f72f0SBarry Smith 
31413f72f0SBarry Smith   PetscFunctionBegin;
32413f72f0SBarry Smith   if (mapping->data) PetscFunctionReturn(0);
33413f72f0SBarry Smith   end   = 0;
34413f72f0SBarry Smith   start = PETSC_MAX_INT;
35413f72f0SBarry Smith 
36413f72f0SBarry Smith   for (i=0; i<n; i++) {
37413f72f0SBarry Smith     if (idx[i] < 0) continue;
38413f72f0SBarry Smith     if (idx[i] < start) start = idx[i];
39413f72f0SBarry Smith     if (idx[i] > end)   end   = idx[i];
40413f72f0SBarry Smith   }
41413f72f0SBarry Smith   if (start > end) {start = 0; end = -1;}
42413f72f0SBarry Smith   mapping->globalstart = start;
43413f72f0SBarry Smith   mapping->globalend   = end;
44413f72f0SBarry Smith   if (!((PetscObject)mapping)->type_name) {
45413f72f0SBarry Smith     if ((end - start) > PetscMax(4*n,1000000)) {
46413f72f0SBarry Smith       ierr = ISLocalToGlobalMappingSetType(mapping,ISLOCALTOGLOBALMAPPINGHASH);
47413f72f0SBarry Smith     } else {
48413f72f0SBarry Smith       ierr = ISLocalToGlobalMappingSetType(mapping,ISLOCALTOGLOBALMAPPINGBASIC);
49413f72f0SBarry Smith     }
50413f72f0SBarry Smith   }
51413f72f0SBarry Smith   ierr = (*mapping->ops->globaltolocalmappingsetup)(mapping);CHKERRQ(ierr);
52413f72f0SBarry Smith   PetscFunctionReturn(0);
53413f72f0SBarry Smith }
54413f72f0SBarry Smith 
55413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp_Basic(ISLocalToGlobalMapping mapping)
56413f72f0SBarry Smith {
57413f72f0SBarry Smith   PetscErrorCode              ierr;
58413f72f0SBarry Smith   PetscInt                    i,*idx = mapping->indices,n = mapping->n,end,start,*globals;
59413f72f0SBarry Smith   ISLocalToGlobalMapping_Basic *map;
60413f72f0SBarry Smith 
61413f72f0SBarry Smith   PetscFunctionBegin;
62413f72f0SBarry Smith   start            = mapping->globalstart;
63413f72f0SBarry Smith   end              = mapping->globalend;
64413f72f0SBarry Smith   ierr             = PetscNew(&map);CHKERRQ(ierr);
65413f72f0SBarry Smith   ierr             = PetscMalloc1(end-start+2,&globals);CHKERRQ(ierr);
66413f72f0SBarry Smith   map->globals     = globals;
67413f72f0SBarry Smith   for (i=0; i<end-start+1; i++) globals[i] = -1;
68413f72f0SBarry Smith   for (i=0; i<n; i++) {
69413f72f0SBarry Smith     if (idx[i] < 0) continue;
70413f72f0SBarry Smith     globals[idx[i] - start] = i;
71413f72f0SBarry Smith   }
72413f72f0SBarry Smith   mapping->data = (void*)map;
73413f72f0SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)mapping,(end-start+1)*sizeof(PetscInt));CHKERRQ(ierr);
74413f72f0SBarry Smith   PetscFunctionReturn(0);
75413f72f0SBarry Smith }
76413f72f0SBarry Smith 
77413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp_Hash(ISLocalToGlobalMapping mapping)
78413f72f0SBarry Smith {
79413f72f0SBarry Smith   PetscErrorCode              ierr;
80413f72f0SBarry Smith   PetscInt                    i,*idx = mapping->indices,n = mapping->n;
81413f72f0SBarry Smith   ISLocalToGlobalMapping_Hash *map;
82413f72f0SBarry Smith 
83413f72f0SBarry Smith   PetscFunctionBegin;
84413f72f0SBarry Smith   ierr = PetscNew(&map);CHKERRQ(ierr);
85413f72f0SBarry Smith   PetscHashICreate(map->globalht);
86413f72f0SBarry Smith   for (i=0; i<n; i++ ) {
87413f72f0SBarry Smith     if (idx[i] < 0) continue;
88413f72f0SBarry Smith     PetscHashIAdd(map->globalht, idx[i], i);
89413f72f0SBarry Smith   }
90413f72f0SBarry Smith   mapping->data = (void*)map;
91413f72f0SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)mapping,2*n*sizeof(PetscInt));CHKERRQ(ierr);
92413f72f0SBarry Smith   PetscFunctionReturn(0);
93413f72f0SBarry Smith }
94413f72f0SBarry Smith 
95413f72f0SBarry Smith static PetscErrorCode ISLocalToGlobalMappingDestroy_Basic(ISLocalToGlobalMapping mapping)
96413f72f0SBarry Smith {
97413f72f0SBarry Smith   PetscErrorCode              ierr;
98413f72f0SBarry Smith   ISLocalToGlobalMapping_Basic *map  = (ISLocalToGlobalMapping_Basic *)mapping->data;
99413f72f0SBarry Smith 
100413f72f0SBarry Smith   PetscFunctionBegin;
101413f72f0SBarry Smith   if (!map) PetscFunctionReturn(0);
102413f72f0SBarry Smith   ierr = PetscFree(map->globals);CHKERRQ(ierr);
103413f72f0SBarry Smith   ierr = PetscFree(mapping->data);CHKERRQ(ierr);
104413f72f0SBarry Smith   PetscFunctionReturn(0);
105413f72f0SBarry Smith }
106413f72f0SBarry Smith 
107413f72f0SBarry Smith static PetscErrorCode ISLocalToGlobalMappingDestroy_Hash(ISLocalToGlobalMapping mapping)
108413f72f0SBarry Smith {
109413f72f0SBarry Smith   PetscErrorCode              ierr;
110413f72f0SBarry Smith   ISLocalToGlobalMapping_Hash *map  = (ISLocalToGlobalMapping_Hash*)mapping->data;
111413f72f0SBarry Smith 
112413f72f0SBarry Smith   PetscFunctionBegin;
113413f72f0SBarry Smith   if (!map) PetscFunctionReturn(0);
114413f72f0SBarry Smith   PetscHashIDestroy(map->globalht);
115413f72f0SBarry Smith   ierr = PetscFree(mapping->data);CHKERRQ(ierr);
116413f72f0SBarry Smith   PetscFunctionReturn(0);
117413f72f0SBarry Smith }
118413f72f0SBarry Smith 
119413f72f0SBarry Smith #define GTOLTYPE _Basic
120413f72f0SBarry Smith #define GTOLNAME _Basic
121413f72f0SBarry Smith #define GTOL(g, local) do {                  \
122413f72f0SBarry Smith     local = map->globals[g/bs - start];      \
123413f72f0SBarry Smith     local = bs*local + (g % bs);             \
124413f72f0SBarry Smith   } while (0)
125413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
126413f72f0SBarry Smith 
127413f72f0SBarry Smith #define GTOLTYPE _Basic
128413f72f0SBarry Smith #define GTOLNAME Block_Basic
129413f72f0SBarry Smith #define GTOL(g, local) do {                  \
130413f72f0SBarry Smith     local = map->globals[g - start];         \
131413f72f0SBarry Smith   } while (0)
132413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
133413f72f0SBarry Smith 
134413f72f0SBarry Smith #define GTOLTYPE _Hash
135413f72f0SBarry Smith #define GTOLNAME _Hash
136413f72f0SBarry Smith #define GTOL(g, local) do {                  \
137413f72f0SBarry Smith     PetscHashIMap(map->globalht,g/bs,local); \
138413f72f0SBarry Smith     local = bs*local + (g % bs);             \
139413f72f0SBarry Smith    } while (0)
140413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
141413f72f0SBarry Smith 
142413f72f0SBarry Smith #define GTOLTYPE _Hash
143413f72f0SBarry Smith #define GTOLNAME Block_Hash
144413f72f0SBarry Smith #define GTOL(g, local) do {                  \
145413f72f0SBarry Smith     PetscHashIMap(map->globalht,g,local);    \
146413f72f0SBarry Smith   } while (0)
147413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
148413f72f0SBarry Smith 
1496658fb44Sstefano_zampini /*@
1506658fb44Sstefano_zampini     ISLocalToGlobalMappingDuplicate - Duplicates the local to global mapping object
1516658fb44Sstefano_zampini 
1526658fb44Sstefano_zampini     Not Collective
1536658fb44Sstefano_zampini 
1546658fb44Sstefano_zampini     Input Parameter:
1556658fb44Sstefano_zampini .   ltog - local to global mapping
1566658fb44Sstefano_zampini 
1576658fb44Sstefano_zampini     Output Parameter:
1586658fb44Sstefano_zampini .   nltog - the duplicated local to global mapping
1596658fb44Sstefano_zampini 
1606658fb44Sstefano_zampini     Level: advanced
1616658fb44Sstefano_zampini 
1626658fb44Sstefano_zampini     Concepts: mapping^local to global
1636658fb44Sstefano_zampini 
1646658fb44Sstefano_zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
1656658fb44Sstefano_zampini @*/
1666658fb44Sstefano_zampini PetscErrorCode  ISLocalToGlobalMappingDuplicate(ISLocalToGlobalMapping ltog,ISLocalToGlobalMapping* nltog)
1676658fb44Sstefano_zampini {
1686658fb44Sstefano_zampini   PetscErrorCode ierr;
1696658fb44Sstefano_zampini 
1706658fb44Sstefano_zampini   PetscFunctionBegin;
1716658fb44Sstefano_zampini   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
1726658fb44Sstefano_zampini   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)ltog),ltog->bs,ltog->n,ltog->indices,PETSC_COPY_VALUES,nltog);CHKERRQ(ierr);
1736658fb44Sstefano_zampini   PetscFunctionReturn(0);
1746658fb44Sstefano_zampini }
1756658fb44Sstefano_zampini 
176565245c5SBarry Smith /*@
177107e9a97SBarry Smith     ISLocalToGlobalMappingGetSize - Gets the local size of a local to global mapping
1783b9aefa3SBarry Smith 
1793b9aefa3SBarry Smith     Not Collective
1803b9aefa3SBarry Smith 
1813b9aefa3SBarry Smith     Input Parameter:
1823b9aefa3SBarry Smith .   ltog - local to global mapping
1833b9aefa3SBarry Smith 
1843b9aefa3SBarry Smith     Output Parameter:
185107e9a97SBarry Smith .   n - the number of entries in the local mapping, ISLocalToGlobalMappingGetIndices() returns an array of this length
1863b9aefa3SBarry Smith 
1873b9aefa3SBarry Smith     Level: advanced
1883b9aefa3SBarry Smith 
189273d9f13SBarry Smith     Concepts: mapping^local to global
1903b9aefa3SBarry Smith 
1913b9aefa3SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
1923b9aefa3SBarry Smith @*/
1937087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping mapping,PetscInt *n)
1943b9aefa3SBarry Smith {
1953b9aefa3SBarry Smith   PetscFunctionBegin;
1960700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1974482741eSBarry Smith   PetscValidIntPointer(n,2);
198107e9a97SBarry Smith   *n = mapping->bs*mapping->n;
1993b9aefa3SBarry Smith   PetscFunctionReturn(0);
2003b9aefa3SBarry Smith }
2013b9aefa3SBarry Smith 
2025a5d4f66SBarry Smith /*@C
2035a5d4f66SBarry Smith     ISLocalToGlobalMappingView - View a local to global mapping
2045a5d4f66SBarry Smith 
205b9cd556bSLois Curfman McInnes     Not Collective
206b9cd556bSLois Curfman McInnes 
2075a5d4f66SBarry Smith     Input Parameters:
2083b9aefa3SBarry Smith +   ltog - local to global mapping
2093b9aefa3SBarry Smith -   viewer - viewer
2105a5d4f66SBarry Smith 
211a997ad1aSLois Curfman McInnes     Level: advanced
212a997ad1aSLois Curfman McInnes 
213273d9f13SBarry Smith     Concepts: mapping^local to global
2145a5d4f66SBarry Smith 
2155a5d4f66SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
2165a5d4f66SBarry Smith @*/
2177087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingView(ISLocalToGlobalMapping mapping,PetscViewer viewer)
2185a5d4f66SBarry Smith {
21932dcc486SBarry Smith   PetscInt       i;
22032dcc486SBarry Smith   PetscMPIInt    rank;
221ace3abfcSBarry Smith   PetscBool      iascii;
2226849ba73SBarry Smith   PetscErrorCode ierr;
2235a5d4f66SBarry Smith 
2245a5d4f66SBarry Smith   PetscFunctionBegin;
2250700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
2263050cee2SBarry Smith   if (!viewer) {
227ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mapping),&viewer);CHKERRQ(ierr);
2283050cee2SBarry Smith   }
2290700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
2305a5d4f66SBarry Smith 
231ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mapping),&rank);CHKERRQ(ierr);
232251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
23332077d6dSBarry Smith   if (iascii) {
23498c3331eSBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)mapping,viewer);CHKERRQ(ierr);
2351575c14dSBarry Smith     ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
2365a5d4f66SBarry Smith     for (i=0; i<mapping->n; i++) {
2377904a332SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D %D\n",rank,i,mapping->indices[i]);CHKERRQ(ierr);
2386831982aSBarry Smith     }
239b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
2401575c14dSBarry Smith     ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
2411575c14dSBarry Smith   }
2425a5d4f66SBarry Smith   PetscFunctionReturn(0);
2435a5d4f66SBarry Smith }
2445a5d4f66SBarry Smith 
2451f428162SBarry Smith /*@
2462bdab257SBarry Smith     ISLocalToGlobalMappingCreateIS - Creates a mapping between a local (0 to n)
2472bdab257SBarry Smith     ordering and a global parallel ordering.
2482bdab257SBarry Smith 
2490f5bd95cSBarry Smith     Not collective
250b9cd556bSLois Curfman McInnes 
251a997ad1aSLois Curfman McInnes     Input Parameter:
2528c03b21aSDmitry Karpeev .   is - index set containing the global numbers for each local number
2532bdab257SBarry Smith 
254a997ad1aSLois Curfman McInnes     Output Parameter:
2552bdab257SBarry Smith .   mapping - new mapping data structure
2562bdab257SBarry Smith 
257f0413b6fSBarry Smith     Notes: the block size of the IS determines the block size of the mapping
258a997ad1aSLois Curfman McInnes     Level: advanced
259a997ad1aSLois Curfman McInnes 
260273d9f13SBarry Smith     Concepts: mapping^local to global
2612bdab257SBarry Smith 
2627e99dc12SLawrence Mitchell .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetFromOptions()
2632bdab257SBarry Smith @*/
2647087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingCreateIS(IS is,ISLocalToGlobalMapping *mapping)
2652bdab257SBarry Smith {
2666849ba73SBarry Smith   PetscErrorCode ierr;
2673bbf0e92SBarry Smith   PetscInt       n,bs;
2685d0c19d7SBarry Smith   const PetscInt *indices;
2692bdab257SBarry Smith   MPI_Comm       comm;
2703bbf0e92SBarry Smith   PetscBool      isblock;
2713a40ed3dSBarry Smith 
2723a40ed3dSBarry Smith   PetscFunctionBegin;
2730700a824SBarry Smith   PetscValidHeaderSpecific(is,IS_CLASSID,1);
2744482741eSBarry Smith   PetscValidPointer(mapping,2);
2752bdab257SBarry Smith 
2762bdab257SBarry Smith   ierr = PetscObjectGetComm((PetscObject)is,&comm);CHKERRQ(ierr);
2773b9aefa3SBarry Smith   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
2783bbf0e92SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)is,ISBLOCK,&isblock);CHKERRQ(ierr);
2796006e8d2SBarry Smith   if (!isblock) {
280f0413b6fSBarry Smith     ierr = ISGetIndices(is,&indices);CHKERRQ(ierr);
281f0413b6fSBarry Smith     ierr = ISLocalToGlobalMappingCreate(comm,1,n,indices,PETSC_COPY_VALUES,mapping);CHKERRQ(ierr);
2822bdab257SBarry Smith     ierr = ISRestoreIndices(is,&indices);CHKERRQ(ierr);
2836006e8d2SBarry Smith   } else {
2846006e8d2SBarry Smith     ierr = ISGetBlockSize(is,&bs);CHKERRQ(ierr);
285f0413b6fSBarry Smith     ierr = ISBlockGetIndices(is,&indices);CHKERRQ(ierr);
28628bc9809SBarry Smith     ierr = ISLocalToGlobalMappingCreate(comm,bs,n/bs,indices,PETSC_COPY_VALUES,mapping);CHKERRQ(ierr);
287f0413b6fSBarry Smith     ierr = ISBlockRestoreIndices(is,&indices);CHKERRQ(ierr);
2886006e8d2SBarry Smith   }
2893a40ed3dSBarry Smith   PetscFunctionReturn(0);
2902bdab257SBarry Smith }
2915a5d4f66SBarry Smith 
292a4d96a55SJed Brown /*@C
293a4d96a55SJed Brown     ISLocalToGlobalMappingCreateSF - Creates a mapping between a local (0 to n)
294a4d96a55SJed Brown     ordering and a global parallel ordering.
295a4d96a55SJed Brown 
296a4d96a55SJed Brown     Collective
297a4d96a55SJed Brown 
298a4d96a55SJed Brown     Input Parameter:
299a4d96a55SJed Brown +   sf - star forest mapping contiguous local indices to (rank, offset)
300a4d96a55SJed Brown -   start - first global index on this process
301a4d96a55SJed Brown 
302a4d96a55SJed Brown     Output Parameter:
303a4d96a55SJed Brown .   mapping - new mapping data structure
304a4d96a55SJed Brown 
305a4d96a55SJed Brown     Level: advanced
306a4d96a55SJed Brown 
307a4d96a55SJed Brown     Concepts: mapping^local to global
308a4d96a55SJed Brown 
3097e99dc12SLawrence Mitchell .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingSetFromOptions()
310a4d96a55SJed Brown @*/
311a4d96a55SJed Brown PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF sf,PetscInt start,ISLocalToGlobalMapping *mapping)
312a4d96a55SJed Brown {
313a4d96a55SJed Brown   PetscErrorCode ierr;
314a4d96a55SJed Brown   PetscInt       i,maxlocal,nroots,nleaves,*globals,*ltog;
315a4d96a55SJed Brown   const PetscInt *ilocal;
316a4d96a55SJed Brown   MPI_Comm       comm;
317a4d96a55SJed Brown 
318a4d96a55SJed Brown   PetscFunctionBegin;
319a4d96a55SJed Brown   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
320a4d96a55SJed Brown   PetscValidPointer(mapping,3);
321a4d96a55SJed Brown 
322a4d96a55SJed Brown   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
3230298fd71SBarry Smith   ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr);
324f6e5521dSKarl Rupp   if (ilocal) {
325f6e5521dSKarl Rupp     for (i=0,maxlocal=0; i<nleaves; i++) maxlocal = PetscMax(maxlocal,ilocal[i]+1);
326f6e5521dSKarl Rupp   }
327a4d96a55SJed Brown   else maxlocal = nleaves;
328785e854fSJed Brown   ierr = PetscMalloc1(nroots,&globals);CHKERRQ(ierr);
329785e854fSJed Brown   ierr = PetscMalloc1(maxlocal,&ltog);CHKERRQ(ierr);
330a4d96a55SJed Brown   for (i=0; i<nroots; i++) globals[i] = start + i;
331a4d96a55SJed Brown   for (i=0; i<maxlocal; i++) ltog[i] = -1;
332a4d96a55SJed Brown   ierr = PetscSFBcastBegin(sf,MPIU_INT,globals,ltog);CHKERRQ(ierr);
333a4d96a55SJed Brown   ierr = PetscSFBcastEnd(sf,MPIU_INT,globals,ltog);CHKERRQ(ierr);
334f0413b6fSBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,1,maxlocal,ltog,PETSC_OWN_POINTER,mapping);CHKERRQ(ierr);
335a4d96a55SJed Brown   ierr = PetscFree(globals);CHKERRQ(ierr);
336a4d96a55SJed Brown   PetscFunctionReturn(0);
337a4d96a55SJed Brown }
338b46b645bSBarry Smith 
33963fa5c83Sstefano_zampini /*@
34063fa5c83Sstefano_zampini     ISLocalToGlobalMappingSetBlockSize - Sets the blocksize of the mapping
34163fa5c83Sstefano_zampini 
34263fa5c83Sstefano_zampini     Not collective
34363fa5c83Sstefano_zampini 
34463fa5c83Sstefano_zampini     Input Parameters:
34563fa5c83Sstefano_zampini .   mapping - mapping data structure
34663fa5c83Sstefano_zampini .   bs - the blocksize
34763fa5c83Sstefano_zampini 
34863fa5c83Sstefano_zampini     Level: advanced
34963fa5c83Sstefano_zampini 
35063fa5c83Sstefano_zampini     Concepts: mapping^local to global
35163fa5c83Sstefano_zampini 
35263fa5c83Sstefano_zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS()
35363fa5c83Sstefano_zampini @*/
35463fa5c83Sstefano_zampini PetscErrorCode  ISLocalToGlobalMappingSetBlockSize(ISLocalToGlobalMapping mapping,PetscInt bs)
35563fa5c83Sstefano_zampini {
356a59f3c4dSstefano_zampini   PetscInt       *nid;
357a59f3c4dSstefano_zampini   const PetscInt *oid;
358a59f3c4dSstefano_zampini   PetscInt       i,cn,on,obs,nn;
35963fa5c83Sstefano_zampini   PetscErrorCode ierr;
36063fa5c83Sstefano_zampini 
36163fa5c83Sstefano_zampini   PetscFunctionBegin;
36263fa5c83Sstefano_zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
36363fa5c83Sstefano_zampini   if (bs < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid block size %D",bs);
36463fa5c83Sstefano_zampini   if (bs == mapping->bs) PetscFunctionReturn(0);
36563fa5c83Sstefano_zampini   on  = mapping->n;
36663fa5c83Sstefano_zampini   obs = mapping->bs;
36763fa5c83Sstefano_zampini   oid = mapping->indices;
36863fa5c83Sstefano_zampini   nn  = (on*obs)/bs;
36963fa5c83Sstefano_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);
370a59f3c4dSstefano_zampini 
37163fa5c83Sstefano_zampini   ierr = PetscMalloc1(nn,&nid);CHKERRQ(ierr);
372a59f3c4dSstefano_zampini   ierr = ISLocalToGlobalMappingGetIndices(mapping,&oid);CHKERRQ(ierr);
373a59f3c4dSstefano_zampini   for (i=0;i<nn;i++) {
374a59f3c4dSstefano_zampini     PetscInt j;
375a59f3c4dSstefano_zampini     for (j=0,cn=0;j<bs-1;j++) {
376a59f3c4dSstefano_zampini       if (oid[i*bs+j] < 0) { cn++; continue; }
377a59f3c4dSstefano_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]);
378a59f3c4dSstefano_zampini     }
379a59f3c4dSstefano_zampini     if (oid[i*bs+j] < 0) cn++;
3808b7cb0e6Sstefano_zampini     if (cn) {
381a59f3c4dSstefano_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);
382a59f3c4dSstefano_zampini       nid[i] = -1;
3838b7cb0e6Sstefano_zampini     } else {
384a59f3c4dSstefano_zampini       nid[i] = oid[i*bs]/bs;
38563fa5c83Sstefano_zampini     }
38663fa5c83Sstefano_zampini   }
387a59f3c4dSstefano_zampini   ierr = ISLocalToGlobalMappingRestoreIndices(mapping,&oid);CHKERRQ(ierr);
388a59f3c4dSstefano_zampini 
38963fa5c83Sstefano_zampini   mapping->n           = nn;
39063fa5c83Sstefano_zampini   mapping->bs          = bs;
39163fa5c83Sstefano_zampini   ierr                 = PetscFree(mapping->indices);CHKERRQ(ierr);
39263fa5c83Sstefano_zampini   mapping->indices     = nid;
393c9345713Sstefano_zampini   mapping->globalstart = 0;
394c9345713Sstefano_zampini   mapping->globalend   = 0;
395413f72f0SBarry Smith   if (mapping->ops->destroy) {
396413f72f0SBarry Smith     ierr = (*mapping->ops->destroy)(mapping);CHKERRQ(ierr);
397413f72f0SBarry Smith   }
39863fa5c83Sstefano_zampini   PetscFunctionReturn(0);
39963fa5c83Sstefano_zampini }
40063fa5c83Sstefano_zampini 
40145b6f7e9SBarry Smith /*@
40245b6f7e9SBarry Smith     ISLocalToGlobalMappingGetBlockSize - Gets the blocksize of the mapping
40345b6f7e9SBarry Smith     ordering and a global parallel ordering.
40445b6f7e9SBarry Smith 
40545b6f7e9SBarry Smith     Not Collective
40645b6f7e9SBarry Smith 
40745b6f7e9SBarry Smith     Input Parameters:
40845b6f7e9SBarry Smith .   mapping - mapping data structure
40945b6f7e9SBarry Smith 
41045b6f7e9SBarry Smith     Output Parameter:
41145b6f7e9SBarry Smith .   bs - the blocksize
41245b6f7e9SBarry Smith 
41345b6f7e9SBarry Smith     Level: advanced
41445b6f7e9SBarry Smith 
41545b6f7e9SBarry Smith     Concepts: mapping^local to global
41645b6f7e9SBarry Smith 
41745b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS()
41845b6f7e9SBarry Smith @*/
41945b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping mapping,PetscInt *bs)
42045b6f7e9SBarry Smith {
42145b6f7e9SBarry Smith   PetscFunctionBegin;
422cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
42345b6f7e9SBarry Smith   *bs = mapping->bs;
42445b6f7e9SBarry Smith   PetscFunctionReturn(0);
42545b6f7e9SBarry Smith }
42645b6f7e9SBarry Smith 
427ba5bb76aSSatish Balay /*@
42890f02eecSBarry Smith     ISLocalToGlobalMappingCreate - Creates a mapping between a local (0 to n)
42990f02eecSBarry Smith     ordering and a global parallel ordering.
4302362add9SBarry Smith 
43189d82c54SBarry Smith     Not Collective, but communicator may have more than one process
432b9cd556bSLois Curfman McInnes 
4332362add9SBarry Smith     Input Parameters:
43489d82c54SBarry Smith +   comm - MPI communicator
435f0413b6fSBarry Smith .   bs - the block size
43628bc9809SBarry Smith .   n - the number of local elements divided by the block size, or equivalently the number of block indices
43728bc9809SBarry 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
438d5ad8652SBarry Smith -   mode - see PetscCopyMode
4392362add9SBarry Smith 
440a997ad1aSLois Curfman McInnes     Output Parameter:
44190f02eecSBarry Smith .   mapping - new mapping data structure
4422362add9SBarry Smith 
443f0413b6fSBarry 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
444413f72f0SBarry Smith 
445413f72f0SBarry Smith     For "small" problems when using ISGlobalToMappingApply() and ISGlobalToMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
446413f72f0SBarry Smith     this uses more memory but is faster; this approach is not scalable for extremely large mappings. For large problems ISLOCALTOGLOBALMAPPINGHASH is used, this is scalable.
447413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
448413f72f0SBarry Smith 
449a997ad1aSLois Curfman McInnes     Level: advanced
450a997ad1aSLois Curfman McInnes 
451273d9f13SBarry Smith     Concepts: mapping^local to global
4522362add9SBarry Smith 
453413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingSetFromOptions(), ISLOCALTOGLOBALMAPPINGBASIC, ISLOCALTOGLOBALMAPPINGHASH
454413f72f0SBarry Smith           ISLocalToGlobalMappingSetType(), ISLocalToGlobalMappingType
4552362add9SBarry Smith @*/
45660c7cefcSBarry Smith PetscErrorCode  ISLocalToGlobalMappingCreate(MPI_Comm comm,PetscInt bs,PetscInt n,const PetscInt indices[],PetscCopyMode mode,ISLocalToGlobalMapping *mapping)
4572362add9SBarry Smith {
4586849ba73SBarry Smith   PetscErrorCode ierr;
45932dcc486SBarry Smith   PetscInt       *in;
460b46b645bSBarry Smith 
461b46b645bSBarry Smith   PetscFunctionBegin;
46273911063SBarry Smith   if (n) PetscValidIntPointer(indices,3);
4634482741eSBarry Smith   PetscValidPointer(mapping,4);
464b46b645bSBarry Smith 
4650298fd71SBarry Smith   *mapping = NULL;
466607a6623SBarry Smith   ierr = ISInitializePackage();CHKERRQ(ierr);
4672362add9SBarry Smith 
46873107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(*mapping,IS_LTOGM_CLASSID,"ISLocalToGlobalMapping","Local to global mapping","IS",
46960c7cefcSBarry Smith                            comm,ISLocalToGlobalMappingDestroy,ISLocalToGlobalMappingView);CHKERRQ(ierr);
470d4bb536fSBarry Smith   (*mapping)->n             = n;
471f0413b6fSBarry Smith   (*mapping)->bs            = bs;
472268a049cSStefano Zampini   (*mapping)->info_cached   = PETSC_FALSE;
473268a049cSStefano Zampini   (*mapping)->info_free     = PETSC_FALSE;
474268a049cSStefano Zampini   (*mapping)->info_procs    = NULL;
475268a049cSStefano Zampini   (*mapping)->info_numprocs = NULL;
476268a049cSStefano Zampini   (*mapping)->info_indices  = NULL;
477413f72f0SBarry Smith 
478413f72f0SBarry Smith   (*mapping)->ops->globaltolocalmappingapply      = NULL;
479413f72f0SBarry Smith   (*mapping)->ops->globaltolocalmappingapplyblock = NULL;
480413f72f0SBarry Smith   (*mapping)->ops->destroy                        = NULL;
481d5ad8652SBarry Smith   if (mode == PETSC_COPY_VALUES) {
482785e854fSJed Brown     ierr = PetscMalloc1(n,&in);CHKERRQ(ierr);
483d5ad8652SBarry Smith     ierr = PetscMemcpy(in,indices,n*sizeof(PetscInt));CHKERRQ(ierr);
484d5ad8652SBarry Smith     (*mapping)->indices = in;
4856389a1a1SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt));CHKERRQ(ierr);
4866389a1a1SBarry Smith   } else if (mode == PETSC_OWN_POINTER) {
4876389a1a1SBarry Smith     (*mapping)->indices = (PetscInt*)indices;
4886389a1a1SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt));CHKERRQ(ierr);
4896389a1a1SBarry Smith   }
49060c7cefcSBarry Smith   else SETERRQ(comm,PETSC_ERR_SUP,"Cannot currently use PETSC_USE_POINTER");
4913a40ed3dSBarry Smith   PetscFunctionReturn(0);
4922362add9SBarry Smith }
4932362add9SBarry Smith 
494413f72f0SBarry Smith PetscFunctionList ISLocalToGlobalMappingList = NULL;
495413f72f0SBarry Smith 
496413f72f0SBarry Smith 
49790f02eecSBarry Smith /*@
4987e99dc12SLawrence Mitchell    ISLocalToGlobalMappingSetFromOptions - Set mapping options from the options database.
4997e99dc12SLawrence Mitchell 
5007e99dc12SLawrence Mitchell    Not collective
5017e99dc12SLawrence Mitchell 
5027e99dc12SLawrence Mitchell    Input Parameters:
5037e99dc12SLawrence Mitchell .  mapping - mapping data structure
5047e99dc12SLawrence Mitchell 
5057e99dc12SLawrence Mitchell    Level: advanced
5067e99dc12SLawrence Mitchell 
5077e99dc12SLawrence Mitchell @*/
5087e99dc12SLawrence Mitchell PetscErrorCode ISLocalToGlobalMappingSetFromOptions(ISLocalToGlobalMapping mapping)
5097e99dc12SLawrence Mitchell {
5107e99dc12SLawrence Mitchell   PetscErrorCode             ierr;
511413f72f0SBarry Smith   char                       type[256];
512413f72f0SBarry Smith   ISLocalToGlobalMappingType defaulttype = "Not set";
5137e99dc12SLawrence Mitchell   PetscBool                  flg;
5147e99dc12SLawrence Mitchell 
5157e99dc12SLawrence Mitchell   PetscFunctionBegin;
5167e99dc12SLawrence Mitchell   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
517413f72f0SBarry Smith   ierr = ISLocalToGlobalMappingRegisterAll();CHKERRQ(ierr);
5187e99dc12SLawrence Mitchell   ierr = PetscObjectOptionsBegin((PetscObject)mapping);CHKERRQ(ierr);
519413f72f0SBarry Smith   ierr = PetscOptionsFList("-islocaltoglobalmapping_type","ISLocalToGlobalMapping method","ISLocalToGlobalMappingSetType",ISLocalToGlobalMappingList,(char*)(((PetscObject)mapping)->type_name) ? ((PetscObject)mapping)->type_name : defaulttype,type,256,&flg);CHKERRQ(ierr);
520413f72f0SBarry Smith   if (flg) {
521413f72f0SBarry Smith     ierr = ISLocalToGlobalMappingSetType(mapping,type);CHKERRQ(ierr);
522413f72f0SBarry Smith   }
5237e99dc12SLawrence Mitchell   ierr = PetscOptionsEnd();CHKERRQ(ierr);
5247e99dc12SLawrence Mitchell   PetscFunctionReturn(0);
5257e99dc12SLawrence Mitchell }
5267e99dc12SLawrence Mitchell 
5277e99dc12SLawrence Mitchell /*@
52890f02eecSBarry Smith    ISLocalToGlobalMappingDestroy - Destroys a mapping between a local (0 to n)
52990f02eecSBarry Smith    ordering and a global parallel ordering.
53090f02eecSBarry Smith 
5310f5bd95cSBarry Smith    Note Collective
532b9cd556bSLois Curfman McInnes 
53390f02eecSBarry Smith    Input Parameters:
53490f02eecSBarry Smith .  mapping - mapping data structure
53590f02eecSBarry Smith 
536a997ad1aSLois Curfman McInnes    Level: advanced
537a997ad1aSLois Curfman McInnes 
5383acfe500SLois Curfman McInnes .seealso: ISLocalToGlobalMappingCreate()
53990f02eecSBarry Smith @*/
5406bf464f9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping *mapping)
54190f02eecSBarry Smith {
542dfbe8321SBarry Smith   PetscErrorCode ierr;
5435fd66863SKarl Rupp 
5443a40ed3dSBarry Smith   PetscFunctionBegin;
5456bf464f9SBarry Smith   if (!*mapping) PetscFunctionReturn(0);
5466bf464f9SBarry Smith   PetscValidHeaderSpecific((*mapping),IS_LTOGM_CLASSID,1);
547997056adSBarry Smith   if (--((PetscObject)(*mapping))->refct > 0) {*mapping = 0;PetscFunctionReturn(0);}
5486bf464f9SBarry Smith   ierr = PetscFree((*mapping)->indices);CHKERRQ(ierr);
549268a049cSStefano Zampini   ierr = PetscFree((*mapping)->info_procs);CHKERRQ(ierr);
550268a049cSStefano Zampini   ierr = PetscFree((*mapping)->info_numprocs);CHKERRQ(ierr);
551268a049cSStefano Zampini   if ((*mapping)->info_indices) {
552268a049cSStefano Zampini     PetscInt i;
553268a049cSStefano Zampini 
554268a049cSStefano Zampini     ierr = PetscFree(((*mapping)->info_indices)[0]);CHKERRQ(ierr);
555268a049cSStefano Zampini     for (i=1; i<(*mapping)->info_nproc; i++) {
556268a049cSStefano Zampini       ierr = PetscFree(((*mapping)->info_indices)[i]);CHKERRQ(ierr);
557268a049cSStefano Zampini     }
558268a049cSStefano Zampini     ierr = PetscFree((*mapping)->info_indices);CHKERRQ(ierr);
559268a049cSStefano Zampini   }
560413f72f0SBarry Smith   if ((*mapping)->ops->destroy) {
561413f72f0SBarry Smith     ierr = (*(*mapping)->ops->destroy)(*mapping);CHKERRQ(ierr);
562413f72f0SBarry Smith   }
563d38fa0fbSBarry Smith   ierr     = PetscHeaderDestroy(mapping);CHKERRQ(ierr);
564992144d0SBarry Smith   *mapping = 0;
5653a40ed3dSBarry Smith   PetscFunctionReturn(0);
56690f02eecSBarry Smith }
56790f02eecSBarry Smith 
56890f02eecSBarry Smith /*@
5693acfe500SLois Curfman McInnes     ISLocalToGlobalMappingApplyIS - Creates from an IS in the local numbering
5703acfe500SLois Curfman McInnes     a new index set using the global numbering defined in an ISLocalToGlobalMapping
5713acfe500SLois Curfman McInnes     context.
57290f02eecSBarry Smith 
573b9cd556bSLois Curfman McInnes     Not collective
574b9cd556bSLois Curfman McInnes 
57590f02eecSBarry Smith     Input Parameters:
576b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
577b9cd556bSLois Curfman McInnes -   is - index set in local numbering
57890f02eecSBarry Smith 
57990f02eecSBarry Smith     Output Parameters:
58090f02eecSBarry Smith .   newis - index set in global numbering
58190f02eecSBarry Smith 
582a997ad1aSLois Curfman McInnes     Level: advanced
583a997ad1aSLois Curfman McInnes 
584273d9f13SBarry Smith     Concepts: mapping^local to global
5853acfe500SLois Curfman McInnes 
58690f02eecSBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),
587d4bb536fSBarry Smith           ISLocalToGlobalMappingDestroy(), ISGlobalToLocalMappingApply()
58890f02eecSBarry Smith @*/
5897087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping,IS is,IS *newis)
59090f02eecSBarry Smith {
5916849ba73SBarry Smith   PetscErrorCode ierr;
592e24637baSBarry Smith   PetscInt       n,*idxout;
5935d0c19d7SBarry Smith   const PetscInt *idxin;
5943a40ed3dSBarry Smith 
5953a40ed3dSBarry Smith   PetscFunctionBegin;
5960700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
5970700a824SBarry Smith   PetscValidHeaderSpecific(is,IS_CLASSID,2);
5984482741eSBarry Smith   PetscValidPointer(newis,3);
59990f02eecSBarry Smith 
6003b9aefa3SBarry Smith   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
60190f02eecSBarry Smith   ierr = ISGetIndices(is,&idxin);CHKERRQ(ierr);
602785e854fSJed Brown   ierr = PetscMalloc1(n,&idxout);CHKERRQ(ierr);
603e24637baSBarry Smith   ierr = ISLocalToGlobalMappingApply(mapping,n,idxin,idxout);CHKERRQ(ierr);
6043b9aefa3SBarry Smith   ierr = ISRestoreIndices(is,&idxin);CHKERRQ(ierr);
605543f3098SMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),n,idxout,PETSC_OWN_POINTER,newis);CHKERRQ(ierr);
6063a40ed3dSBarry Smith   PetscFunctionReturn(0);
60790f02eecSBarry Smith }
60890f02eecSBarry Smith 
609b89cb25eSSatish Balay /*@
6103acfe500SLois Curfman McInnes    ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering
6113acfe500SLois Curfman McInnes    and converts them to the global numbering.
61290f02eecSBarry Smith 
613b9cd556bSLois Curfman McInnes    Not collective
614b9cd556bSLois Curfman McInnes 
615bb25748dSBarry Smith    Input Parameters:
616b9cd556bSLois Curfman McInnes +  mapping - the local to global mapping context
617bb25748dSBarry Smith .  N - number of integers
618b9cd556bSLois Curfman McInnes -  in - input indices in local numbering
619bb25748dSBarry Smith 
620bb25748dSBarry Smith    Output Parameter:
621bb25748dSBarry Smith .  out - indices in global numbering
622bb25748dSBarry Smith 
623b9cd556bSLois Curfman McInnes    Notes:
624b9cd556bSLois Curfman McInnes    The in and out array parameters may be identical.
625d4bb536fSBarry Smith 
626a997ad1aSLois Curfman McInnes    Level: advanced
627a997ad1aSLois Curfman McInnes 
62845b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
6290752156aSBarry Smith           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
630d4bb536fSBarry Smith           AOPetscToApplication(), ISGlobalToLocalMappingApply()
631bb25748dSBarry Smith 
632273d9f13SBarry Smith     Concepts: mapping^local to global
633afcb2eb5SJed Brown @*/
634afcb2eb5SJed Brown PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
635afcb2eb5SJed Brown {
636cbc1caf0SMatthew G. Knepley   PetscInt i,bs,Nmax;
63745b6f7e9SBarry Smith 
63845b6f7e9SBarry Smith   PetscFunctionBegin;
639cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
640cbc1caf0SMatthew G. Knepley   bs   = mapping->bs;
641cbc1caf0SMatthew G. Knepley   Nmax = bs*mapping->n;
64245b6f7e9SBarry Smith   if (bs == 1) {
643cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
64445b6f7e9SBarry Smith     for (i=0; i<N; i++) {
64545b6f7e9SBarry Smith       if (in[i] < 0) {
64645b6f7e9SBarry Smith         out[i] = in[i];
64745b6f7e9SBarry Smith         continue;
64845b6f7e9SBarry Smith       }
649e24637baSBarry 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);
65045b6f7e9SBarry Smith       out[i] = idx[in[i]];
65145b6f7e9SBarry Smith     }
65245b6f7e9SBarry Smith   } else {
653cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
65445b6f7e9SBarry Smith     for (i=0; i<N; i++) {
65545b6f7e9SBarry Smith       if (in[i] < 0) {
65645b6f7e9SBarry Smith         out[i] = in[i];
65745b6f7e9SBarry Smith         continue;
65845b6f7e9SBarry Smith       }
659e24637baSBarry 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);
66045b6f7e9SBarry Smith       out[i] = idx[in[i]/bs]*bs + (in[i] % bs);
66145b6f7e9SBarry Smith     }
66245b6f7e9SBarry Smith   }
66345b6f7e9SBarry Smith   PetscFunctionReturn(0);
66445b6f7e9SBarry Smith }
66545b6f7e9SBarry Smith 
66645b6f7e9SBarry Smith /*@
6676006e8d2SBarry Smith    ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local block numbering and converts them to the global block numbering
66845b6f7e9SBarry Smith 
66945b6f7e9SBarry Smith    Not collective
67045b6f7e9SBarry Smith 
67145b6f7e9SBarry Smith    Input Parameters:
67245b6f7e9SBarry Smith +  mapping - the local to global mapping context
67345b6f7e9SBarry Smith .  N - number of integers
6746006e8d2SBarry Smith -  in - input indices in local block numbering
67545b6f7e9SBarry Smith 
67645b6f7e9SBarry Smith    Output Parameter:
6776006e8d2SBarry Smith .  out - indices in global block numbering
67845b6f7e9SBarry Smith 
67945b6f7e9SBarry Smith    Notes:
68045b6f7e9SBarry Smith    The in and out array parameters may be identical.
68145b6f7e9SBarry Smith 
6826006e8d2SBarry Smith    Example:
6836006e8d2SBarry 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
6846006e8d2SBarry Smith      (the first block) would produce 0 and the mapping applied to 1 (the second block) would produce 3.
6856006e8d2SBarry Smith 
68645b6f7e9SBarry Smith    Level: advanced
68745b6f7e9SBarry Smith 
68845b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
68945b6f7e9SBarry Smith           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
69045b6f7e9SBarry Smith           AOPetscToApplication(), ISGlobalToLocalMappingApply()
69145b6f7e9SBarry Smith 
69245b6f7e9SBarry Smith     Concepts: mapping^local to global
69345b6f7e9SBarry Smith @*/
69445b6f7e9SBarry Smith PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
69545b6f7e9SBarry Smith {
696cbc1caf0SMatthew G. Knepley 
697cbc1caf0SMatthew G. Knepley   PetscFunctionBegin;
698cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
699cbc1caf0SMatthew G. Knepley   {
700afcb2eb5SJed Brown     PetscInt i,Nmax = mapping->n;
701afcb2eb5SJed Brown     const PetscInt *idx = mapping->indices;
702d4bb536fSBarry Smith 
703afcb2eb5SJed Brown     for (i=0; i<N; i++) {
704afcb2eb5SJed Brown       if (in[i] < 0) {
705afcb2eb5SJed Brown         out[i] = in[i];
706afcb2eb5SJed Brown         continue;
707afcb2eb5SJed Brown       }
708e24637baSBarry 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);
709afcb2eb5SJed Brown       out[i] = idx[in[i]];
710afcb2eb5SJed Brown     }
711cbc1caf0SMatthew G. Knepley   }
712afcb2eb5SJed Brown   PetscFunctionReturn(0);
713afcb2eb5SJed Brown }
714d4bb536fSBarry Smith 
7157e99dc12SLawrence Mitchell /*@
716a997ad1aSLois Curfman McInnes     ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers
717a997ad1aSLois Curfman McInnes     specified with a global numbering.
718d4bb536fSBarry Smith 
719b9cd556bSLois Curfman McInnes     Not collective
720b9cd556bSLois Curfman McInnes 
721d4bb536fSBarry Smith     Input Parameters:
722b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
723d4bb536fSBarry Smith .   type - IS_GTOLM_MASK - replaces global indices with no local value with -1
724d4bb536fSBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
725d4bb536fSBarry Smith .   n - number of global indices to map
726b9cd556bSLois Curfman McInnes -   idx - global indices to map
727d4bb536fSBarry Smith 
728d4bb536fSBarry Smith     Output Parameters:
729b9cd556bSLois Curfman McInnes +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
730b9cd556bSLois Curfman McInnes -   idxout - local index of each global index, one must pass in an array long enough
731e182c471SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApply() with
7320298fd71SBarry Smith              idxout == NULL to determine the required length (returned in nout)
733e182c471SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApply()
734e182c471SBarry Smith              a second time to set the values.
735d4bb536fSBarry Smith 
736b9cd556bSLois Curfman McInnes     Notes:
7370298fd71SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
738d4bb536fSBarry Smith 
739413f72f0SBarry Smith     For "small" problems when using ISGlobalToMappingApply() and ISGlobalToMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
740413f72f0SBarry Smith     this uses more memory but is faster; this approach is not scalable for extremely large mappings. For large problems ISLOCALTOGLOBALMAPPINGHASH is used, this is scalable.
741413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
7420f5bd95cSBarry Smith 
743a997ad1aSLois Curfman McInnes     Level: advanced
744a997ad1aSLois Curfman McInnes 
74532fd6b96SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
74632fd6b96SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
74732fd6b96SBarry Smith 
748273d9f13SBarry Smith     Concepts: mapping^global to local
749d4bb536fSBarry Smith 
7509d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),
751413f72f0SBarry Smith           ISLocalToGlobalMappingDestroy()
752d4bb536fSBarry Smith @*/
753413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
754d4bb536fSBarry Smith {
7559d90f715SBarry Smith   PetscErrorCode ierr;
7569d90f715SBarry Smith 
7579d90f715SBarry Smith   PetscFunctionBegin;
7589d90f715SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
759413f72f0SBarry Smith   if (!mapping->data) {
760413f72f0SBarry Smith     ierr = ISGlobalToLocalMappingSetUp(mapping);CHKERRQ(ierr);
7619d90f715SBarry Smith   }
762413f72f0SBarry Smith   ierr = (*mapping->ops->globaltolocalmappingapply)(mapping,type,n,idx,nout,idxout);CHKERRQ(ierr);
7639d90f715SBarry Smith   PetscFunctionReturn(0);
7649d90f715SBarry Smith }
7659d90f715SBarry Smith 
766d4fe737eSStefano Zampini /*@
767d4fe737eSStefano Zampini     ISGlobalToLocalMappingApplyIS - Creates from an IS in the global numbering
768d4fe737eSStefano Zampini     a new index set using the local numbering defined in an ISLocalToGlobalMapping
769d4fe737eSStefano Zampini     context.
770d4fe737eSStefano Zampini 
771d4fe737eSStefano Zampini     Not collective
772d4fe737eSStefano Zampini 
773d4fe737eSStefano Zampini     Input Parameters:
774d4fe737eSStefano Zampini +   mapping - mapping between local and global numbering
775d4fe737eSStefano Zampini -   is - index set in global numbering
776d4fe737eSStefano Zampini 
777d4fe737eSStefano Zampini     Output Parameters:
778d4fe737eSStefano Zampini .   newis - index set in local numbering
779d4fe737eSStefano Zampini 
780d4fe737eSStefano Zampini     Level: advanced
781d4fe737eSStefano Zampini 
782d4fe737eSStefano Zampini     Concepts: mapping^local to global
783d4fe737eSStefano Zampini 
784d4fe737eSStefano Zampini .seealso: ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(),
785d4fe737eSStefano Zampini           ISLocalToGlobalMappingDestroy()
786d4fe737eSStefano Zampini @*/
787413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type, IS is,IS *newis)
788d4fe737eSStefano Zampini {
789d4fe737eSStefano Zampini   PetscErrorCode ierr;
790d4fe737eSStefano Zampini   PetscInt       n,nout,*idxout;
791d4fe737eSStefano Zampini   const PetscInt *idxin;
792d4fe737eSStefano Zampini 
793d4fe737eSStefano Zampini   PetscFunctionBegin;
794d4fe737eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
795d4fe737eSStefano Zampini   PetscValidHeaderSpecific(is,IS_CLASSID,3);
796d4fe737eSStefano Zampini   PetscValidPointer(newis,4);
797d4fe737eSStefano Zampini 
798d4fe737eSStefano Zampini   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
799d4fe737eSStefano Zampini   ierr = ISGetIndices(is,&idxin);CHKERRQ(ierr);
800d4fe737eSStefano Zampini   if (type == IS_GTOLM_MASK) {
801d4fe737eSStefano Zampini     ierr = PetscMalloc1(n,&idxout);CHKERRQ(ierr);
802d4fe737eSStefano Zampini   } else {
803d4fe737eSStefano Zampini     ierr = ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,NULL);CHKERRQ(ierr);
804d4fe737eSStefano Zampini     ierr = PetscMalloc1(nout,&idxout);CHKERRQ(ierr);
805d4fe737eSStefano Zampini   }
806d4fe737eSStefano Zampini   ierr = ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,idxout);CHKERRQ(ierr);
807d4fe737eSStefano Zampini   ierr = ISRestoreIndices(is,&idxin);CHKERRQ(ierr);
808d4fe737eSStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),nout,idxout,PETSC_OWN_POINTER,newis);CHKERRQ(ierr);
809d4fe737eSStefano Zampini   PetscFunctionReturn(0);
810d4fe737eSStefano Zampini }
811d4fe737eSStefano Zampini 
8129d90f715SBarry Smith /*@
8139d90f715SBarry Smith     ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers
8149d90f715SBarry Smith     specified with a block global numbering.
8159d90f715SBarry Smith 
8169d90f715SBarry Smith     Not collective
8179d90f715SBarry Smith 
8189d90f715SBarry Smith     Input Parameters:
8199d90f715SBarry Smith +   mapping - mapping between local and global numbering
8209d90f715SBarry Smith .   type - IS_GTOLM_MASK - replaces global indices with no local value with -1
8219d90f715SBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
8229d90f715SBarry Smith .   n - number of global indices to map
8239d90f715SBarry Smith -   idx - global indices to map
8249d90f715SBarry Smith 
8259d90f715SBarry Smith     Output Parameters:
8269d90f715SBarry Smith +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
8279d90f715SBarry Smith -   idxout - local index of each global index, one must pass in an array long enough
8289d90f715SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApplyBlock() with
8299d90f715SBarry Smith              idxout == NULL to determine the required length (returned in nout)
8309d90f715SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApplyBlock()
8319d90f715SBarry Smith              a second time to set the values.
8329d90f715SBarry Smith 
8339d90f715SBarry Smith     Notes:
8349d90f715SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
8359d90f715SBarry Smith 
836413f72f0SBarry Smith     For "small" problems when using ISGlobalToMappingApply() and ISGlobalToMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
837413f72f0SBarry Smith     this uses more memory but is faster; this approach is not scalable for extremely large mappings. For large problems ISLOCALTOGLOBALMAPPINGHASH is used, this is scalable.
838413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
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(),
848413f72f0SBarry Smith           ISLocalToGlobalMappingDestroy()
8499d90f715SBarry Smith @*/
850413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,
8519d90f715SBarry Smith                                                  PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
8529d90f715SBarry Smith {
8536849ba73SBarry Smith   PetscErrorCode ierr;
854d4bb536fSBarry Smith 
8553a40ed3dSBarry Smith   PetscFunctionBegin;
8560700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
857413f72f0SBarry Smith   if (!mapping->data) {
858413f72f0SBarry Smith     ierr = ISGlobalToLocalMappingSetUp(mapping);CHKERRQ(ierr);
859d4bb536fSBarry Smith   }
860413f72f0SBarry Smith   ierr = (*mapping->ops->globaltolocalmappingapplyblock)(mapping,type,n,idx,nout,idxout);CHKERRQ(ierr);
8613a40ed3dSBarry Smith   PetscFunctionReturn(0);
862d4bb536fSBarry Smith }
86390f02eecSBarry Smith 
864413f72f0SBarry Smith 
86589d82c54SBarry Smith /*@C
8666a818285SBarry Smith     ISLocalToGlobalMappingGetBlockInfo - Gets the neighbor information for each processor and
86789d82c54SBarry Smith      each index shared by more than one processor
86889d82c54SBarry Smith 
86989d82c54SBarry Smith     Collective on ISLocalToGlobalMapping
87089d82c54SBarry Smith 
87189d82c54SBarry Smith     Input Parameters:
87289d82c54SBarry Smith .   mapping - the mapping from local to global indexing
87389d82c54SBarry Smith 
87489d82c54SBarry Smith     Output Parameter:
87589d82c54SBarry Smith +   nproc - number of processors that are connected to this one
87689d82c54SBarry Smith .   proc - neighboring processors
87707b52d57SBarry Smith .   numproc - number of indices for each subdomain (processor)
8783463a7baSJed Brown -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
87989d82c54SBarry Smith 
88089d82c54SBarry Smith     Level: advanced
88189d82c54SBarry Smith 
882273d9f13SBarry Smith     Concepts: mapping^local to global
88389d82c54SBarry Smith 
8842cfcea29SBarry Smith     Fortran Usage:
8852cfcea29SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
8862cfcea29SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
8872cfcea29SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
8882cfcea29SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
8892cfcea29SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
8902cfcea29SBarry Smith 
8912cfcea29SBarry Smith 
89207b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
89307b52d57SBarry Smith           ISLocalToGlobalMappingRestoreInfo()
89489d82c54SBarry Smith @*/
8956a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
89689d82c54SBarry Smith {
8976849ba73SBarry Smith   PetscErrorCode ierr;
898268a049cSStefano Zampini 
899268a049cSStefano Zampini   PetscFunctionBegin;
900268a049cSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
901268a049cSStefano Zampini   if (mapping->info_cached) {
902268a049cSStefano Zampini     *nproc    = mapping->info_nproc;
903268a049cSStefano Zampini     *procs    = mapping->info_procs;
904268a049cSStefano Zampini     *numprocs = mapping->info_numprocs;
905268a049cSStefano Zampini     *indices  = mapping->info_indices;
906268a049cSStefano Zampini   } else {
907268a049cSStefano Zampini     ierr = ISLocalToGlobalMappingGetBlockInfo_Private(mapping,nproc,procs,numprocs,indices);CHKERRQ(ierr);
908268a049cSStefano Zampini   }
909268a049cSStefano Zampini   PetscFunctionReturn(0);
910268a049cSStefano Zampini }
911268a049cSStefano Zampini 
912268a049cSStefano Zampini static PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
913268a049cSStefano Zampini {
914268a049cSStefano Zampini   PetscErrorCode ierr;
91597f1f81fSBarry Smith   PetscMPIInt    size,rank,tag1,tag2,tag3,*len,*source,imdex;
91632dcc486SBarry Smith   PetscInt       i,n = mapping->n,Ng,ng,max = 0,*lindices = mapping->indices;
91732dcc486SBarry Smith   PetscInt       *nprocs,*owner,nsends,*sends,j,*starts,nmax,nrecvs,*recvs,proc;
91897f1f81fSBarry Smith   PetscInt       cnt,scale,*ownedsenders,*nownedsenders,rstart,nowned;
91932dcc486SBarry Smith   PetscInt       node,nownedm,nt,*sends2,nsends2,*starts2,*lens2,*dest,nrecvs2,*starts3,*recvs2,k,*bprocs,*tmp;
92032dcc486SBarry Smith   PetscInt       first_procs,first_numprocs,*first_indices;
92189d82c54SBarry Smith   MPI_Request    *recv_waits,*send_waits;
92230dcb7c9SBarry Smith   MPI_Status     recv_status,*send_status,*recv_statuses;
923ce94432eSBarry Smith   MPI_Comm       comm;
924ace3abfcSBarry Smith   PetscBool      debug = PETSC_FALSE;
92589d82c54SBarry Smith 
92689d82c54SBarry Smith   PetscFunctionBegin;
927ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)mapping,&comm);CHKERRQ(ierr);
92824cf384cSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
92924cf384cSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
93024cf384cSBarry Smith   if (size == 1) {
93124cf384cSBarry Smith     *nproc         = 0;
9320298fd71SBarry Smith     *procs         = NULL;
93395dccacaSBarry Smith     ierr           = PetscNew(numprocs);CHKERRQ(ierr);
9341e2105dcSBarry Smith     (*numprocs)[0] = 0;
93595dccacaSBarry Smith     ierr           = PetscNew(indices);CHKERRQ(ierr);
9360298fd71SBarry Smith     (*indices)[0]  = NULL;
937268a049cSStefano Zampini     /* save info for reuse */
938268a049cSStefano Zampini     mapping->info_nproc = *nproc;
939268a049cSStefano Zampini     mapping->info_procs = *procs;
940268a049cSStefano Zampini     mapping->info_numprocs = *numprocs;
941268a049cSStefano Zampini     mapping->info_indices = *indices;
942268a049cSStefano Zampini     mapping->info_cached = PETSC_TRUE;
94324cf384cSBarry Smith     PetscFunctionReturn(0);
94424cf384cSBarry Smith   }
94524cf384cSBarry Smith 
946c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)mapping)->options,NULL,"-islocaltoglobalmappinggetinfo_debug",&debug,NULL);CHKERRQ(ierr);
94707b52d57SBarry Smith 
9483677ff5aSBarry Smith   /*
9496a818285SBarry Smith     Notes on ISLocalToGlobalMappingGetBlockInfo
9503677ff5aSBarry Smith 
9513677ff5aSBarry Smith     globally owned node - the nodes that have been assigned to this processor in global
9523677ff5aSBarry Smith            numbering, just for this routine.
9533677ff5aSBarry Smith 
9543677ff5aSBarry Smith     nontrivial globally owned node - node assigned to this processor that is on a subdomain
9553677ff5aSBarry Smith            boundary (i.e. is has more than one local owner)
9563677ff5aSBarry Smith 
9573677ff5aSBarry Smith     locally owned node - node that exists on this processors subdomain
9583677ff5aSBarry Smith 
9593677ff5aSBarry Smith     nontrivial locally owned node - node that is not in the interior (i.e. has more than one
9603677ff5aSBarry Smith            local subdomain
9613677ff5aSBarry Smith   */
96224cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag1);CHKERRQ(ierr);
96324cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag2);CHKERRQ(ierr);
96424cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag3);CHKERRQ(ierr);
96589d82c54SBarry Smith 
96689d82c54SBarry Smith   for (i=0; i<n; i++) {
96789d82c54SBarry Smith     if (lindices[i] > max) max = lindices[i];
96889d82c54SBarry Smith   }
969b2566f29SBarry Smith   ierr   = MPIU_Allreduce(&max,&Ng,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr);
97078058e43SBarry Smith   Ng++;
97189d82c54SBarry Smith   ierr   = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
97289d82c54SBarry Smith   ierr   = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
973bc8ff85bSBarry Smith   scale  = Ng/size + 1;
974a2e34c3dSBarry Smith   ng     = scale; if (rank == size-1) ng = Ng - scale*(size-1); ng = PetscMax(1,ng);
975caba0dd0SBarry Smith   rstart = scale*rank;
97689d82c54SBarry Smith 
97789d82c54SBarry Smith   /* determine ownership ranges of global indices */
978785e854fSJed Brown   ierr = PetscMalloc1(2*size,&nprocs);CHKERRQ(ierr);
97932dcc486SBarry Smith   ierr = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr);
98089d82c54SBarry Smith 
98189d82c54SBarry Smith   /* determine owners of each local node  */
982785e854fSJed Brown   ierr = PetscMalloc1(n,&owner);CHKERRQ(ierr);
98389d82c54SBarry Smith   for (i=0; i<n; i++) {
9843677ff5aSBarry Smith     proc             = lindices[i]/scale; /* processor that globally owns this index */
98527c402fcSBarry Smith     nprocs[2*proc+1] = 1;                 /* processor globally owns at least one of ours */
9863677ff5aSBarry Smith     owner[i]         = proc;
98727c402fcSBarry Smith     nprocs[2*proc]++;                     /* count of how many that processor globally owns of ours */
98889d82c54SBarry Smith   }
98927c402fcSBarry Smith   nsends = 0; for (i=0; i<size; i++) nsends += nprocs[2*i+1];
9907904a332SBarry Smith   ierr = PetscInfo1(mapping,"Number of global owners for my local data %D\n",nsends);CHKERRQ(ierr);
99189d82c54SBarry Smith 
99289d82c54SBarry Smith   /* inform other processors of number of messages and max length*/
99327c402fcSBarry Smith   ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr);
9947904a332SBarry Smith   ierr = PetscInfo1(mapping,"Number of local owners for my global data %D\n",nrecvs);CHKERRQ(ierr);
99589d82c54SBarry Smith 
99689d82c54SBarry Smith   /* post receives for owned rows */
997785e854fSJed Brown   ierr = PetscMalloc1((2*nrecvs+1)*(nmax+1),&recvs);CHKERRQ(ierr);
998854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&recv_waits);CHKERRQ(ierr);
99989d82c54SBarry Smith   for (i=0; i<nrecvs; i++) {
100032dcc486SBarry Smith     ierr = MPI_Irecv(recvs+2*nmax*i,2*nmax,MPIU_INT,MPI_ANY_SOURCE,tag1,comm,recv_waits+i);CHKERRQ(ierr);
100189d82c54SBarry Smith   }
100289d82c54SBarry Smith 
100389d82c54SBarry Smith   /* pack messages containing lists of local nodes to owners */
1004854ce69bSBarry Smith   ierr      = PetscMalloc1(2*n+1,&sends);CHKERRQ(ierr);
1005854ce69bSBarry Smith   ierr      = PetscMalloc1(size+1,&starts);CHKERRQ(ierr);
100689d82c54SBarry Smith   starts[0] = 0;
1007f6e5521dSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
100889d82c54SBarry Smith   for (i=0; i<n; i++) {
100989d82c54SBarry Smith     sends[starts[owner[i]]++] = lindices[i];
101030dcb7c9SBarry Smith     sends[starts[owner[i]]++] = i;
101189d82c54SBarry Smith   }
101289d82c54SBarry Smith   ierr = PetscFree(owner);CHKERRQ(ierr);
101389d82c54SBarry Smith   starts[0] = 0;
1014f6e5521dSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
101589d82c54SBarry Smith 
101689d82c54SBarry Smith   /* send the messages */
1017854ce69bSBarry Smith   ierr = PetscMalloc1(nsends+1,&send_waits);CHKERRQ(ierr);
1018854ce69bSBarry Smith   ierr = PetscMalloc1(nsends+1,&dest);CHKERRQ(ierr);
101989d82c54SBarry Smith   cnt = 0;
102089d82c54SBarry Smith   for (i=0; i<size; i++) {
102127c402fcSBarry Smith     if (nprocs[2*i]) {
102232dcc486SBarry Smith       ierr      = MPI_Isend(sends+starts[i],2*nprocs[2*i],MPIU_INT,i,tag1,comm,send_waits+cnt);CHKERRQ(ierr);
102330dcb7c9SBarry Smith       dest[cnt] = i;
102489d82c54SBarry Smith       cnt++;
102589d82c54SBarry Smith     }
102689d82c54SBarry Smith   }
102789d82c54SBarry Smith   ierr = PetscFree(starts);CHKERRQ(ierr);
102889d82c54SBarry Smith 
102989d82c54SBarry Smith   /* wait on receives */
1030854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&source);CHKERRQ(ierr);
1031854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&len);CHKERRQ(ierr);
103289d82c54SBarry Smith   cnt  = nrecvs;
1033854ce69bSBarry Smith   ierr = PetscMalloc1(ng+1,&nownedsenders);CHKERRQ(ierr);
103432dcc486SBarry Smith   ierr = PetscMemzero(nownedsenders,ng*sizeof(PetscInt));CHKERRQ(ierr);
103589d82c54SBarry Smith   while (cnt) {
103689d82c54SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
103789d82c54SBarry Smith     /* unpack receives into our local space */
103832dcc486SBarry Smith     ierr          = MPI_Get_count(&recv_status,MPIU_INT,&len[imdex]);CHKERRQ(ierr);
103989d82c54SBarry Smith     source[imdex] = recv_status.MPI_SOURCE;
104030dcb7c9SBarry Smith     len[imdex]    = len[imdex]/2;
1041caba0dd0SBarry Smith     /* count how many local owners for each of my global owned indices */
104230dcb7c9SBarry Smith     for (i=0; i<len[imdex]; i++) nownedsenders[recvs[2*imdex*nmax+2*i]-rstart]++;
104389d82c54SBarry Smith     cnt--;
104489d82c54SBarry Smith   }
104589d82c54SBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
104689d82c54SBarry Smith 
104730dcb7c9SBarry Smith   /* count how many globally owned indices are on an edge multiplied by how many processors own them. */
1048bc8ff85bSBarry Smith   nowned  = 0;
1049bc8ff85bSBarry Smith   nownedm = 0;
1050bc8ff85bSBarry Smith   for (i=0; i<ng; i++) {
1051bc8ff85bSBarry Smith     if (nownedsenders[i] > 1) {nownedm += nownedsenders[i]; nowned++;}
1052bc8ff85bSBarry Smith   }
1053bc8ff85bSBarry Smith 
1054bc8ff85bSBarry Smith   /* create single array to contain rank of all local owners of each globally owned index */
1055854ce69bSBarry Smith   ierr      = PetscMalloc1(nownedm+1,&ownedsenders);CHKERRQ(ierr);
1056854ce69bSBarry Smith   ierr      = PetscMalloc1(ng+1,&starts);CHKERRQ(ierr);
1057bc8ff85bSBarry Smith   starts[0] = 0;
1058bc8ff85bSBarry Smith   for (i=1; i<ng; i++) {
1059bc8ff85bSBarry Smith     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
1060bc8ff85bSBarry Smith     else starts[i] = starts[i-1];
1061bc8ff85bSBarry Smith   }
1062bc8ff85bSBarry Smith 
106330dcb7c9SBarry Smith   /* for each nontrival globally owned node list all arriving processors */
1064bc8ff85bSBarry Smith   for (i=0; i<nrecvs; i++) {
1065bc8ff85bSBarry Smith     for (j=0; j<len[i]; j++) {
106630dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
1067f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) ownedsenders[starts[node]++] = source[i];
1068bc8ff85bSBarry Smith     }
1069bc8ff85bSBarry Smith   }
1070bc8ff85bSBarry Smith 
107107b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
107230dcb7c9SBarry Smith     starts[0] = 0;
107330dcb7c9SBarry Smith     for (i=1; i<ng; i++) {
107430dcb7c9SBarry Smith       if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
107530dcb7c9SBarry Smith       else starts[i] = starts[i-1];
107630dcb7c9SBarry Smith     }
107730dcb7c9SBarry Smith     for (i=0; i<ng; i++) {
107830dcb7c9SBarry Smith       if (nownedsenders[i] > 1) {
10797904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"[%d] global node %D local owner processors: ",rank,i+rstart);CHKERRQ(ierr);
108030dcb7c9SBarry Smith         for (j=0; j<nownedsenders[i]; j++) {
10817904a332SBarry Smith           ierr = PetscSynchronizedPrintf(comm,"%D ",ownedsenders[starts[i]+j]);CHKERRQ(ierr);
108230dcb7c9SBarry Smith         }
108330dcb7c9SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
108430dcb7c9SBarry Smith       }
108530dcb7c9SBarry Smith     }
10860ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
108707b52d57SBarry Smith   } /* -----------------------------------  */
108830dcb7c9SBarry Smith 
10893677ff5aSBarry Smith   /* wait on original sends */
10903a96401aSBarry Smith   if (nsends) {
1091785e854fSJed Brown     ierr = PetscMalloc1(nsends,&send_status);CHKERRQ(ierr);
10923a96401aSBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
10933a96401aSBarry Smith     ierr = PetscFree(send_status);CHKERRQ(ierr);
10943a96401aSBarry Smith   }
109589d82c54SBarry Smith   ierr = PetscFree(send_waits);CHKERRQ(ierr);
10963a96401aSBarry Smith   ierr = PetscFree(sends);CHKERRQ(ierr);
10973677ff5aSBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
10983677ff5aSBarry Smith 
10993677ff5aSBarry Smith   /* pack messages to send back to local owners */
110030dcb7c9SBarry Smith   starts[0] = 0;
110130dcb7c9SBarry Smith   for (i=1; i<ng; i++) {
110230dcb7c9SBarry Smith     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
110330dcb7c9SBarry Smith     else starts[i] = starts[i-1];
110430dcb7c9SBarry Smith   }
110530dcb7c9SBarry Smith   nsends2 = nrecvs;
1106854ce69bSBarry Smith   ierr    = PetscMalloc1(nsends2+1,&nprocs);CHKERRQ(ierr); /* length of each message */
110730dcb7c9SBarry Smith   for (i=0; i<nrecvs; i++) {
110830dcb7c9SBarry Smith     nprocs[i] = 1;
110930dcb7c9SBarry Smith     for (j=0; j<len[i]; j++) {
111030dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
1111f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) nprocs[i] += 2 + nownedsenders[node];
111230dcb7c9SBarry Smith     }
111330dcb7c9SBarry Smith   }
1114f6e5521dSKarl Rupp   nt = 0;
1115f6e5521dSKarl Rupp   for (i=0; i<nsends2; i++) nt += nprocs[i];
1116f6e5521dSKarl Rupp 
1117854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,&sends2);CHKERRQ(ierr);
1118854ce69bSBarry Smith   ierr = PetscMalloc1(nsends2+1,&starts2);CHKERRQ(ierr);
1119f6e5521dSKarl Rupp 
1120f6e5521dSKarl Rupp   starts2[0] = 0;
1121f6e5521dSKarl Rupp   for (i=1; i<nsends2; i++) starts2[i] = starts2[i-1] + nprocs[i-1];
112230dcb7c9SBarry Smith   /*
112330dcb7c9SBarry Smith      Each message is 1 + nprocs[i] long, and consists of
112430dcb7c9SBarry Smith        (0) the number of nodes being sent back
112530dcb7c9SBarry Smith        (1) the local node number,
112630dcb7c9SBarry Smith        (2) the number of processors sharing it,
112730dcb7c9SBarry Smith        (3) the processors sharing it
112830dcb7c9SBarry Smith   */
112930dcb7c9SBarry Smith   for (i=0; i<nsends2; i++) {
113030dcb7c9SBarry Smith     cnt = 1;
113130dcb7c9SBarry Smith     sends2[starts2[i]] = 0;
113230dcb7c9SBarry Smith     for (j=0; j<len[i]; j++) {
113330dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
113430dcb7c9SBarry Smith       if (nownedsenders[node] > 1) {
113530dcb7c9SBarry Smith         sends2[starts2[i]]++;
113630dcb7c9SBarry Smith         sends2[starts2[i]+cnt++] = recvs[2*i*nmax+2*j+1];
113730dcb7c9SBarry Smith         sends2[starts2[i]+cnt++] = nownedsenders[node];
113832dcc486SBarry Smith         ierr = PetscMemcpy(&sends2[starts2[i]+cnt],&ownedsenders[starts[node]],nownedsenders[node]*sizeof(PetscInt));CHKERRQ(ierr);
113930dcb7c9SBarry Smith         cnt += nownedsenders[node];
114030dcb7c9SBarry Smith       }
114130dcb7c9SBarry Smith     }
114230dcb7c9SBarry Smith   }
114330dcb7c9SBarry Smith 
114430dcb7c9SBarry Smith   /* receive the message lengths */
114530dcb7c9SBarry Smith   nrecvs2 = nsends;
1146854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&lens2);CHKERRQ(ierr);
1147854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&starts3);CHKERRQ(ierr);
1148854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&recv_waits);CHKERRQ(ierr);
114930dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
1150d44834fbSBarry Smith     ierr = MPI_Irecv(&lens2[i],1,MPIU_INT,dest[i],tag2,comm,recv_waits+i);CHKERRQ(ierr);
115130dcb7c9SBarry Smith   }
1152d44834fbSBarry Smith 
11538a8e0b3aSBarry Smith   /* send the message lengths */
11548a8e0b3aSBarry Smith   for (i=0; i<nsends2; i++) {
11558a8e0b3aSBarry Smith     ierr = MPI_Send(&nprocs[i],1,MPIU_INT,source[i],tag2,comm);CHKERRQ(ierr);
11568a8e0b3aSBarry Smith   }
11578a8e0b3aSBarry Smith 
1158d44834fbSBarry Smith   /* wait on receives of lens */
11590c468ba9SBarry Smith   if (nrecvs2) {
1160785e854fSJed Brown     ierr = PetscMalloc1(nrecvs2,&recv_statuses);CHKERRQ(ierr);
1161d44834fbSBarry Smith     ierr = MPI_Waitall(nrecvs2,recv_waits,recv_statuses);CHKERRQ(ierr);
1162d44834fbSBarry Smith     ierr = PetscFree(recv_statuses);CHKERRQ(ierr);
11630c468ba9SBarry Smith   }
1164a2ea699eSBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
1165d44834fbSBarry Smith 
116630dcb7c9SBarry Smith   starts3[0] = 0;
1167d44834fbSBarry Smith   nt         = 0;
116830dcb7c9SBarry Smith   for (i=0; i<nrecvs2-1; i++) {
116930dcb7c9SBarry Smith     starts3[i+1] = starts3[i] + lens2[i];
1170d44834fbSBarry Smith     nt          += lens2[i];
117130dcb7c9SBarry Smith   }
117276466f69SStefano Zampini   if (nrecvs2) nt += lens2[nrecvs2-1];
1173d44834fbSBarry Smith 
1174854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,&recvs2);CHKERRQ(ierr);
1175854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs2+1,&recv_waits);CHKERRQ(ierr);
117652b72c4aSBarry Smith   for (i=0; i<nrecvs2; i++) {
117732dcc486SBarry Smith     ierr = MPI_Irecv(recvs2+starts3[i],lens2[i],MPIU_INT,dest[i],tag3,comm,recv_waits+i);CHKERRQ(ierr);
117830dcb7c9SBarry Smith   }
117930dcb7c9SBarry Smith 
118030dcb7c9SBarry Smith   /* send the messages */
1181854ce69bSBarry Smith   ierr = PetscMalloc1(nsends2+1,&send_waits);CHKERRQ(ierr);
118230dcb7c9SBarry Smith   for (i=0; i<nsends2; i++) {
118332dcc486SBarry Smith     ierr = MPI_Isend(sends2+starts2[i],nprocs[i],MPIU_INT,source[i],tag3,comm,send_waits+i);CHKERRQ(ierr);
118430dcb7c9SBarry Smith   }
118530dcb7c9SBarry Smith 
118630dcb7c9SBarry Smith   /* wait on receives */
11870c468ba9SBarry Smith   if (nrecvs2) {
1188785e854fSJed Brown     ierr = PetscMalloc1(nrecvs2,&recv_statuses);CHKERRQ(ierr);
118930dcb7c9SBarry Smith     ierr = MPI_Waitall(nrecvs2,recv_waits,recv_statuses);CHKERRQ(ierr);
119030dcb7c9SBarry Smith     ierr = PetscFree(recv_statuses);CHKERRQ(ierr);
11910c468ba9SBarry Smith   }
119230dcb7c9SBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
119330dcb7c9SBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
119430dcb7c9SBarry Smith 
119507b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
119630dcb7c9SBarry Smith     cnt = 0;
119730dcb7c9SBarry Smith     for (i=0; i<nrecvs2; i++) {
119830dcb7c9SBarry Smith       nt = recvs2[cnt++];
119930dcb7c9SBarry Smith       for (j=0; j<nt; j++) {
12007904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"[%d] local node %D number of subdomains %D: ",rank,recvs2[cnt],recvs2[cnt+1]);CHKERRQ(ierr);
120130dcb7c9SBarry Smith         for (k=0; k<recvs2[cnt+1]; k++) {
12027904a332SBarry Smith           ierr = PetscSynchronizedPrintf(comm,"%D ",recvs2[cnt+2+k]);CHKERRQ(ierr);
120330dcb7c9SBarry Smith         }
120430dcb7c9SBarry Smith         cnt += 2 + recvs2[cnt+1];
120530dcb7c9SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
120630dcb7c9SBarry Smith       }
120730dcb7c9SBarry Smith     }
12080ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
120907b52d57SBarry Smith   } /* -----------------------------------  */
121030dcb7c9SBarry Smith 
121130dcb7c9SBarry Smith   /* count number subdomains for each local node */
1212785e854fSJed Brown   ierr = PetscMalloc1(size,&nprocs);CHKERRQ(ierr);
121332dcc486SBarry Smith   ierr = PetscMemzero(nprocs,size*sizeof(PetscInt));CHKERRQ(ierr);
121430dcb7c9SBarry Smith   cnt  = 0;
121530dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
121630dcb7c9SBarry Smith     nt = recvs2[cnt++];
121730dcb7c9SBarry Smith     for (j=0; j<nt; j++) {
1218f6e5521dSKarl Rupp       for (k=0; k<recvs2[cnt+1]; k++) nprocs[recvs2[cnt+2+k]]++;
121930dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt+1];
122030dcb7c9SBarry Smith     }
122130dcb7c9SBarry Smith   }
122230dcb7c9SBarry Smith   nt = 0; for (i=0; i<size; i++) nt += (nprocs[i] > 0);
122330dcb7c9SBarry Smith   *nproc    = nt;
1224854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,procs);CHKERRQ(ierr);
1225854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,numprocs);CHKERRQ(ierr);
1226854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,indices);CHKERRQ(ierr);
12270298fd71SBarry Smith   for (i=0;i<nt+1;i++) (*indices)[i]=NULL;
1228785e854fSJed Brown   ierr = PetscMalloc1(size,&bprocs);CHKERRQ(ierr);
122930dcb7c9SBarry Smith   cnt  = 0;
123030dcb7c9SBarry Smith   for (i=0; i<size; i++) {
123130dcb7c9SBarry Smith     if (nprocs[i] > 0) {
123230dcb7c9SBarry Smith       bprocs[i]        = cnt;
123330dcb7c9SBarry Smith       (*procs)[cnt]    = i;
123430dcb7c9SBarry Smith       (*numprocs)[cnt] = nprocs[i];
1235785e854fSJed Brown       ierr             = PetscMalloc1(nprocs[i],&(*indices)[cnt]);CHKERRQ(ierr);
123630dcb7c9SBarry Smith       cnt++;
123730dcb7c9SBarry Smith     }
123830dcb7c9SBarry Smith   }
123930dcb7c9SBarry Smith 
124030dcb7c9SBarry Smith   /* make the list of subdomains for each nontrivial local node */
124132dcc486SBarry Smith   ierr = PetscMemzero(*numprocs,nt*sizeof(PetscInt));CHKERRQ(ierr);
124230dcb7c9SBarry Smith   cnt  = 0;
124330dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
124430dcb7c9SBarry Smith     nt = recvs2[cnt++];
124530dcb7c9SBarry Smith     for (j=0; j<nt; j++) {
1246f6e5521dSKarl Rupp       for (k=0; k<recvs2[cnt+1]; k++) (*indices)[bprocs[recvs2[cnt+2+k]]][(*numprocs)[bprocs[recvs2[cnt+2+k]]]++] = recvs2[cnt];
124730dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt+1];
124830dcb7c9SBarry Smith     }
124930dcb7c9SBarry Smith   }
125030dcb7c9SBarry Smith   ierr = PetscFree(bprocs);CHKERRQ(ierr);
125107b52d57SBarry Smith   ierr = PetscFree(recvs2);CHKERRQ(ierr);
125230dcb7c9SBarry Smith 
125307b52d57SBarry Smith   /* sort the node indexing by their global numbers */
125407b52d57SBarry Smith   nt = *nproc;
125507b52d57SBarry Smith   for (i=0; i<nt; i++) {
1256854ce69bSBarry Smith     ierr = PetscMalloc1((*numprocs)[i],&tmp);CHKERRQ(ierr);
1257f6e5521dSKarl Rupp     for (j=0; j<(*numprocs)[i]; j++) tmp[j] = lindices[(*indices)[i][j]];
125807b52d57SBarry Smith     ierr = PetscSortIntWithArray((*numprocs)[i],tmp,(*indices)[i]);CHKERRQ(ierr);
125907b52d57SBarry Smith     ierr = PetscFree(tmp);CHKERRQ(ierr);
126007b52d57SBarry Smith   }
126107b52d57SBarry Smith 
126207b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
126330dcb7c9SBarry Smith     nt = *nproc;
126430dcb7c9SBarry Smith     for (i=0; i<nt; i++) {
12657904a332SBarry Smith       ierr = PetscSynchronizedPrintf(comm,"[%d] subdomain %D number of indices %D: ",rank,(*procs)[i],(*numprocs)[i]);CHKERRQ(ierr);
126630dcb7c9SBarry Smith       for (j=0; j<(*numprocs)[i]; j++) {
12677904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"%D ",(*indices)[i][j]);CHKERRQ(ierr);
126830dcb7c9SBarry Smith       }
126930dcb7c9SBarry Smith       ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
127030dcb7c9SBarry Smith     }
12710ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
127207b52d57SBarry Smith   } /* -----------------------------------  */
127330dcb7c9SBarry Smith 
127430dcb7c9SBarry Smith   /* wait on sends */
127530dcb7c9SBarry Smith   if (nsends2) {
1276785e854fSJed Brown     ierr = PetscMalloc1(nsends2,&send_status);CHKERRQ(ierr);
127730dcb7c9SBarry Smith     ierr = MPI_Waitall(nsends2,send_waits,send_status);CHKERRQ(ierr);
127830dcb7c9SBarry Smith     ierr = PetscFree(send_status);CHKERRQ(ierr);
127930dcb7c9SBarry Smith   }
128030dcb7c9SBarry Smith 
128130dcb7c9SBarry Smith   ierr = PetscFree(starts3);CHKERRQ(ierr);
128230dcb7c9SBarry Smith   ierr = PetscFree(dest);CHKERRQ(ierr);
128330dcb7c9SBarry Smith   ierr = PetscFree(send_waits);CHKERRQ(ierr);
12843677ff5aSBarry Smith 
1285bc8ff85bSBarry Smith   ierr = PetscFree(nownedsenders);CHKERRQ(ierr);
1286bc8ff85bSBarry Smith   ierr = PetscFree(ownedsenders);CHKERRQ(ierr);
1287bc8ff85bSBarry Smith   ierr = PetscFree(starts);CHKERRQ(ierr);
128830dcb7c9SBarry Smith   ierr = PetscFree(starts2);CHKERRQ(ierr);
128930dcb7c9SBarry Smith   ierr = PetscFree(lens2);CHKERRQ(ierr);
129089d82c54SBarry Smith 
129189d82c54SBarry Smith   ierr = PetscFree(source);CHKERRQ(ierr);
129297f1f81fSBarry Smith   ierr = PetscFree(len);CHKERRQ(ierr);
129389d82c54SBarry Smith   ierr = PetscFree(recvs);CHKERRQ(ierr);
12943a96401aSBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
129530dcb7c9SBarry Smith   ierr = PetscFree(sends2);CHKERRQ(ierr);
129624cf384cSBarry Smith 
129724cf384cSBarry Smith   /* put the information about myself as the first entry in the list */
129824cf384cSBarry Smith   first_procs    = (*procs)[0];
129924cf384cSBarry Smith   first_numprocs = (*numprocs)[0];
130024cf384cSBarry Smith   first_indices  = (*indices)[0];
130124cf384cSBarry Smith   for (i=0; i<*nproc; i++) {
130224cf384cSBarry Smith     if ((*procs)[i] == rank) {
130324cf384cSBarry Smith       (*procs)[0]    = (*procs)[i];
130424cf384cSBarry Smith       (*numprocs)[0] = (*numprocs)[i];
130524cf384cSBarry Smith       (*indices)[0]  = (*indices)[i];
130624cf384cSBarry Smith       (*procs)[i]    = first_procs;
130724cf384cSBarry Smith       (*numprocs)[i] = first_numprocs;
130824cf384cSBarry Smith       (*indices)[i]  = first_indices;
130924cf384cSBarry Smith       break;
131024cf384cSBarry Smith     }
131124cf384cSBarry Smith   }
1312268a049cSStefano Zampini 
1313268a049cSStefano Zampini   /* save info for reuse */
1314268a049cSStefano Zampini   mapping->info_nproc = *nproc;
1315268a049cSStefano Zampini   mapping->info_procs = *procs;
1316268a049cSStefano Zampini   mapping->info_numprocs = *numprocs;
1317268a049cSStefano Zampini   mapping->info_indices = *indices;
1318268a049cSStefano Zampini   mapping->info_cached = PETSC_TRUE;
131989d82c54SBarry Smith   PetscFunctionReturn(0);
132089d82c54SBarry Smith }
132189d82c54SBarry Smith 
13226a818285SBarry Smith /*@C
13236a818285SBarry Smith     ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by ISLocalToGlobalMappingGetBlockInfo()
13246a818285SBarry Smith 
13256a818285SBarry Smith     Collective on ISLocalToGlobalMapping
13266a818285SBarry Smith 
13276a818285SBarry Smith     Input Parameters:
13286a818285SBarry Smith .   mapping - the mapping from local to global indexing
13296a818285SBarry Smith 
13306a818285SBarry Smith     Output Parameter:
13316a818285SBarry Smith +   nproc - number of processors that are connected to this one
13326a818285SBarry Smith .   proc - neighboring processors
13336a818285SBarry Smith .   numproc - number of indices for each processor
13346a818285SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
13356a818285SBarry Smith 
13366a818285SBarry Smith     Level: advanced
13376a818285SBarry Smith 
13386a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
13396a818285SBarry Smith           ISLocalToGlobalMappingGetInfo()
13406a818285SBarry Smith @*/
13416a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
13426a818285SBarry Smith {
13436a818285SBarry Smith   PetscErrorCode ierr;
13446a818285SBarry Smith 
13456a818285SBarry Smith   PetscFunctionBegin;
1346cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1347268a049cSStefano Zampini   if (mapping->info_free) {
13486a818285SBarry Smith     ierr = PetscFree(*numprocs);CHKERRQ(ierr);
13496a818285SBarry Smith     if (*indices) {
1350268a049cSStefano Zampini       PetscInt i;
1351268a049cSStefano Zampini 
13526a818285SBarry Smith       ierr = PetscFree((*indices)[0]);CHKERRQ(ierr);
13536a818285SBarry Smith       for (i=1; i<*nproc; i++) {
13546a818285SBarry Smith         ierr = PetscFree((*indices)[i]);CHKERRQ(ierr);
13556a818285SBarry Smith       }
13566a818285SBarry Smith       ierr = PetscFree(*indices);CHKERRQ(ierr);
13576a818285SBarry Smith     }
1358268a049cSStefano Zampini   }
1359268a049cSStefano Zampini   *nproc    = 0;
1360268a049cSStefano Zampini   *procs    = NULL;
1361268a049cSStefano Zampini   *numprocs = NULL;
1362268a049cSStefano Zampini   *indices  = NULL;
13636a818285SBarry Smith   PetscFunctionReturn(0);
13646a818285SBarry Smith }
13656a818285SBarry Smith 
13666a818285SBarry Smith /*@C
13676a818285SBarry Smith     ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each processor and
13686a818285SBarry Smith      each index shared by more than one processor
13696a818285SBarry Smith 
13706a818285SBarry Smith     Collective on ISLocalToGlobalMapping
13716a818285SBarry Smith 
13726a818285SBarry Smith     Input Parameters:
13736a818285SBarry Smith .   mapping - the mapping from local to global indexing
13746a818285SBarry Smith 
13756a818285SBarry Smith     Output Parameter:
13766a818285SBarry Smith +   nproc - number of processors that are connected to this one
13776a818285SBarry Smith .   proc - neighboring processors
13786a818285SBarry Smith .   numproc - number of indices for each subdomain (processor)
13796a818285SBarry Smith -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
13806a818285SBarry Smith 
13816a818285SBarry Smith     Level: advanced
13826a818285SBarry Smith 
13836a818285SBarry Smith     Concepts: mapping^local to global
13846a818285SBarry Smith 
13856a818285SBarry Smith     Fortran Usage:
13866a818285SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
13876a818285SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
13886a818285SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
13896a818285SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
13906a818285SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
13916a818285SBarry Smith 
13926a818285SBarry Smith 
13936a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
13946a818285SBarry Smith           ISLocalToGlobalMappingRestoreInfo()
13956a818285SBarry Smith @*/
13966a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
13976a818285SBarry Smith {
13986a818285SBarry Smith   PetscErrorCode ierr;
1399268a049cSStefano Zampini   PetscInt       **bindices = NULL,*bnumprocs = NULL,bs = mapping->bs,i,j,k;
14006a818285SBarry Smith 
14016a818285SBarry Smith   PetscFunctionBegin;
14026a818285SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1403268a049cSStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockInfo(mapping,nproc,procs,&bnumprocs,&bindices);CHKERRQ(ierr);
1404268a049cSStefano Zampini   if (bs > 1) { /* we need to expand the cached info */
1405732f65e3SBarry Smith     ierr = PetscCalloc1(*nproc,&*indices);CHKERRQ(ierr);
1406268a049cSStefano Zampini     ierr = PetscCalloc1(*nproc,&*numprocs);CHKERRQ(ierr);
14076a818285SBarry Smith     for (i=0; i<*nproc; i++) {
1408268a049cSStefano Zampini       ierr = PetscMalloc1(bs*bnumprocs[i],&(*indices)[i]);CHKERRQ(ierr);
1409268a049cSStefano Zampini       for (j=0; j<bnumprocs[i]; j++) {
14106a818285SBarry Smith         for (k=0; k<bs; k++) {
14116a818285SBarry Smith           (*indices)[i][j*bs+k] = bs*bindices[i][j] + k;
14126a818285SBarry Smith         }
14136a818285SBarry Smith       }
1414268a049cSStefano Zampini       (*numprocs)[i] = bnumprocs[i]*bs;
14156a818285SBarry Smith     }
1416268a049cSStefano Zampini     mapping->info_free = PETSC_TRUE;
1417268a049cSStefano Zampini   } else {
1418268a049cSStefano Zampini     *numprocs = bnumprocs;
1419268a049cSStefano Zampini     *indices  = bindices;
14206a818285SBarry Smith   }
14216a818285SBarry Smith   PetscFunctionReturn(0);
14226a818285SBarry Smith }
14236a818285SBarry Smith 
142407b52d57SBarry Smith /*@C
142507b52d57SBarry Smith     ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by ISLocalToGlobalMappingGetInfo()
142689d82c54SBarry Smith 
142707b52d57SBarry Smith     Collective on ISLocalToGlobalMapping
142807b52d57SBarry Smith 
142907b52d57SBarry Smith     Input Parameters:
143007b52d57SBarry Smith .   mapping - the mapping from local to global indexing
143107b52d57SBarry Smith 
143207b52d57SBarry Smith     Output Parameter:
143307b52d57SBarry Smith +   nproc - number of processors that are connected to this one
143407b52d57SBarry Smith .   proc - neighboring processors
143507b52d57SBarry Smith .   numproc - number of indices for each processor
143607b52d57SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
143707b52d57SBarry Smith 
143807b52d57SBarry Smith     Level: advanced
143907b52d57SBarry Smith 
144007b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
144107b52d57SBarry Smith           ISLocalToGlobalMappingGetInfo()
144207b52d57SBarry Smith @*/
14437087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
144407b52d57SBarry Smith {
14456849ba73SBarry Smith   PetscErrorCode ierr;
144607b52d57SBarry Smith 
144707b52d57SBarry Smith   PetscFunctionBegin;
14486a818285SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockInfo(mapping,nproc,procs,numprocs,indices);CHKERRQ(ierr);
144907b52d57SBarry Smith   PetscFunctionReturn(0);
145007b52d57SBarry Smith }
145186994e45SJed Brown 
145286994e45SJed Brown /*@C
1453107e9a97SBarry Smith    ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped
145486994e45SJed Brown 
145586994e45SJed Brown    Not Collective
145686994e45SJed Brown 
145786994e45SJed Brown    Input Arguments:
145886994e45SJed Brown . ltog - local to global mapping
145986994e45SJed Brown 
146086994e45SJed Brown    Output Arguments:
1461565245c5SBarry Smith . array - array of indices, the length of this array may be obtained with ISLocalToGlobalMappingGetSize()
146286994e45SJed Brown 
146386994e45SJed Brown    Level: advanced
146486994e45SJed Brown 
1465107e9a97SBarry Smith    Notes: ISLocalToGlobalMappingGetSize() returns the length the this array
1466107e9a97SBarry Smith 
1467107e9a97SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreIndices(), ISLocalToGlobalMappingGetBlockIndices(), ISLocalToGlobalMappingRestoreBlockIndices()
146886994e45SJed Brown @*/
14697087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
147086994e45SJed Brown {
147186994e45SJed Brown   PetscFunctionBegin;
147286994e45SJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
147386994e45SJed Brown   PetscValidPointer(array,2);
147445b6f7e9SBarry Smith   if (ltog->bs == 1) {
147586994e45SJed Brown     *array = ltog->indices;
147645b6f7e9SBarry Smith   } else {
147745b6f7e9SBarry Smith     PetscInt       *jj,k,i,j,n = ltog->n, bs = ltog->bs;
147845b6f7e9SBarry Smith     const PetscInt *ii;
147945b6f7e9SBarry Smith     PetscErrorCode ierr;
148045b6f7e9SBarry Smith 
148145b6f7e9SBarry Smith     ierr = PetscMalloc1(bs*n,&jj);CHKERRQ(ierr);
148245b6f7e9SBarry Smith     *array = jj;
148345b6f7e9SBarry Smith     k    = 0;
148445b6f7e9SBarry Smith     ii   = ltog->indices;
148545b6f7e9SBarry Smith     for (i=0; i<n; i++)
148645b6f7e9SBarry Smith       for (j=0; j<bs; j++)
148745b6f7e9SBarry Smith         jj[k++] = bs*ii[i] + j;
148845b6f7e9SBarry Smith   }
148986994e45SJed Brown   PetscFunctionReturn(0);
149086994e45SJed Brown }
149186994e45SJed Brown 
149286994e45SJed Brown /*@C
1493193a2b41SJulian Andrej    ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with ISLocalToGlobalMappingGetIndices()
149486994e45SJed Brown 
149586994e45SJed Brown    Not Collective
149686994e45SJed Brown 
149786994e45SJed Brown    Input Arguments:
149886994e45SJed Brown + ltog - local to global mapping
149986994e45SJed Brown - array - array of indices
150086994e45SJed Brown 
150186994e45SJed Brown    Level: advanced
150286994e45SJed Brown 
150386994e45SJed Brown .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices()
150486994e45SJed Brown @*/
15057087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
150686994e45SJed Brown {
150786994e45SJed Brown   PetscFunctionBegin;
150886994e45SJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
150986994e45SJed Brown   PetscValidPointer(array,2);
151045b6f7e9SBarry Smith   if (ltog->bs == 1 && *array != ltog->indices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
151145b6f7e9SBarry Smith 
151245b6f7e9SBarry Smith   if (ltog->bs > 1) {
151345b6f7e9SBarry Smith     PetscErrorCode ierr;
151445b6f7e9SBarry Smith     ierr = PetscFree(*(void**)array);CHKERRQ(ierr);
151545b6f7e9SBarry Smith   }
151645b6f7e9SBarry Smith   PetscFunctionReturn(0);
151745b6f7e9SBarry Smith }
151845b6f7e9SBarry Smith 
151945b6f7e9SBarry Smith /*@C
152045b6f7e9SBarry Smith    ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block
152145b6f7e9SBarry Smith 
152245b6f7e9SBarry Smith    Not Collective
152345b6f7e9SBarry Smith 
152445b6f7e9SBarry Smith    Input Arguments:
152545b6f7e9SBarry Smith . ltog - local to global mapping
152645b6f7e9SBarry Smith 
152745b6f7e9SBarry Smith    Output Arguments:
152845b6f7e9SBarry Smith . array - array of indices
152945b6f7e9SBarry Smith 
153045b6f7e9SBarry Smith    Level: advanced
153145b6f7e9SBarry Smith 
153245b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreBlockIndices()
153345b6f7e9SBarry Smith @*/
153445b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
153545b6f7e9SBarry Smith {
153645b6f7e9SBarry Smith   PetscFunctionBegin;
153745b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
153845b6f7e9SBarry Smith   PetscValidPointer(array,2);
153945b6f7e9SBarry Smith   *array = ltog->indices;
154045b6f7e9SBarry Smith   PetscFunctionReturn(0);
154145b6f7e9SBarry Smith }
154245b6f7e9SBarry Smith 
154345b6f7e9SBarry Smith /*@C
154445b6f7e9SBarry Smith    ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with ISLocalToGlobalMappingGetBlockIndices()
154545b6f7e9SBarry Smith 
154645b6f7e9SBarry Smith    Not Collective
154745b6f7e9SBarry Smith 
154845b6f7e9SBarry Smith    Input Arguments:
154945b6f7e9SBarry Smith + ltog - local to global mapping
155045b6f7e9SBarry Smith - array - array of indices
155145b6f7e9SBarry Smith 
155245b6f7e9SBarry Smith    Level: advanced
155345b6f7e9SBarry Smith 
155445b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices()
155545b6f7e9SBarry Smith @*/
155645b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
155745b6f7e9SBarry Smith {
155845b6f7e9SBarry Smith   PetscFunctionBegin;
155945b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
156045b6f7e9SBarry Smith   PetscValidPointer(array,2);
156186994e45SJed Brown   if (*array != ltog->indices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
15620298fd71SBarry Smith   *array = NULL;
156386994e45SJed Brown   PetscFunctionReturn(0);
156486994e45SJed Brown }
1565f7efa3c7SJed Brown 
1566f7efa3c7SJed Brown /*@C
1567f7efa3c7SJed Brown    ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings
1568f7efa3c7SJed Brown 
1569f7efa3c7SJed Brown    Not Collective
1570f7efa3c7SJed Brown 
1571f7efa3c7SJed Brown    Input Arguments:
1572f7efa3c7SJed Brown + comm - communicator for the new mapping, must contain the communicator of every mapping to concatenate
1573f7efa3c7SJed Brown . n - number of mappings to concatenate
1574f7efa3c7SJed Brown - ltogs - local to global mappings
1575f7efa3c7SJed Brown 
1576f7efa3c7SJed Brown    Output Arguments:
1577f7efa3c7SJed Brown . ltogcat - new mapping
1578f7efa3c7SJed Brown 
15799d90f715SBarry Smith    Note: this currently always returns a mapping with block size of 1
15809d90f715SBarry Smith 
15819d90f715SBarry Smith    Developer Note: If all the input mapping have the same block size we could easily handle that as a special case
15829d90f715SBarry Smith 
1583f7efa3c7SJed Brown    Level: advanced
1584f7efa3c7SJed Brown 
1585f7efa3c7SJed Brown .seealso: ISLocalToGlobalMappingCreate()
1586f7efa3c7SJed Brown @*/
1587f7efa3c7SJed Brown PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm,PetscInt n,const ISLocalToGlobalMapping ltogs[],ISLocalToGlobalMapping *ltogcat)
1588f7efa3c7SJed Brown {
1589f7efa3c7SJed Brown   PetscInt       i,cnt,m,*idx;
1590f7efa3c7SJed Brown   PetscErrorCode ierr;
1591f7efa3c7SJed Brown 
1592f7efa3c7SJed Brown   PetscFunctionBegin;
1593f7efa3c7SJed Brown   if (n < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have a non-negative number of mappings, given %D",n);
1594f7efa3c7SJed Brown   if (n > 0) PetscValidPointer(ltogs,3);
1595f7efa3c7SJed Brown   for (i=0; i<n; i++) PetscValidHeaderSpecific(ltogs[i],IS_LTOGM_CLASSID,3);
1596f7efa3c7SJed Brown   PetscValidPointer(ltogcat,4);
1597f7efa3c7SJed Brown   for (cnt=0,i=0; i<n; i++) {
1598f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
1599f7efa3c7SJed Brown     cnt += m;
1600f7efa3c7SJed Brown   }
1601785e854fSJed Brown   ierr = PetscMalloc1(cnt,&idx);CHKERRQ(ierr);
1602f7efa3c7SJed Brown   for (cnt=0,i=0; i<n; i++) {
1603f7efa3c7SJed Brown     const PetscInt *subidx;
1604f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
1605f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetIndices(ltogs[i],&subidx);CHKERRQ(ierr);
1606f7efa3c7SJed Brown     ierr = PetscMemcpy(&idx[cnt],subidx,m*sizeof(PetscInt));CHKERRQ(ierr);
1607f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingRestoreIndices(ltogs[i],&subidx);CHKERRQ(ierr);
1608f7efa3c7SJed Brown     cnt += m;
1609f7efa3c7SJed Brown   }
1610f0413b6fSBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,1,cnt,idx,PETSC_OWN_POINTER,ltogcat);CHKERRQ(ierr);
1611f7efa3c7SJed Brown   PetscFunctionReturn(0);
1612f7efa3c7SJed Brown }
161304a59952SBarry Smith 
1614413f72f0SBarry Smith /*MC
1615413f72f0SBarry Smith       ISLOCALTOGLOBALMAPPINGBASIC - basic implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is
1616413f72f0SBarry Smith                                     used this is good for only small and moderate size problems.
1617413f72f0SBarry Smith 
1618413f72f0SBarry Smith    Options Database Keys:
1619413f72f0SBarry Smith +   -islocaltoglobalmapping_type basic - select this method
1620413f72f0SBarry Smith 
1621413f72f0SBarry Smith    Level: beginner
1622413f72f0SBarry Smith 
1623413f72f0SBarry Smith .seealso:  ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType(), ISLOCALTOGLOBALMAPPINGHASH
1624413f72f0SBarry Smith M*/
1625413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Basic(ISLocalToGlobalMapping ltog)
1626413f72f0SBarry Smith {
1627413f72f0SBarry Smith   PetscFunctionBegin;
1628413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Basic;
1629413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Basic;
1630413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Basic;
1631413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Basic;
1632413f72f0SBarry Smith   PetscFunctionReturn(0);
1633413f72f0SBarry Smith }
1634413f72f0SBarry Smith 
1635413f72f0SBarry Smith /*MC
1636413f72f0SBarry Smith       ISLOCALTOGLOBALMAPPINGHASH - hash implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is
1637ed56e8eaSBarry Smith                                     used this is good for large memory problems.
1638413f72f0SBarry Smith 
1639413f72f0SBarry Smith    Options Database Keys:
1640413f72f0SBarry Smith +   -islocaltoglobalmapping_type hash - select this method
1641413f72f0SBarry Smith 
1642ed56e8eaSBarry Smith 
1643ed56e8eaSBarry Smith    Notes: This is selected automatically for large problems if the user does not set the type.
1644ed56e8eaSBarry Smith 
1645413f72f0SBarry Smith    Level: beginner
1646413f72f0SBarry Smith 
1647413f72f0SBarry Smith .seealso:  ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType(), ISLOCALTOGLOBALMAPPINGHASH
1648413f72f0SBarry Smith M*/
1649413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Hash(ISLocalToGlobalMapping ltog)
1650413f72f0SBarry Smith {
1651413f72f0SBarry Smith   PetscFunctionBegin;
1652413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Hash;
1653413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Hash;
1654413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Hash;
1655413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Hash;
1656413f72f0SBarry Smith   PetscFunctionReturn(0);
1657413f72f0SBarry Smith }
1658413f72f0SBarry Smith 
1659413f72f0SBarry Smith 
1660413f72f0SBarry Smith /*@C
1661413f72f0SBarry Smith     ISLocalToGlobalMappingRegister -  Adds a method for applying a global to local mapping with an ISLocalToGlobalMapping
1662413f72f0SBarry Smith 
1663413f72f0SBarry Smith    Not Collective
1664413f72f0SBarry Smith 
1665413f72f0SBarry Smith    Input Parameters:
1666413f72f0SBarry Smith +  sname - name of a new method
1667413f72f0SBarry Smith -  routine_create - routine to create method context
1668413f72f0SBarry Smith 
1669413f72f0SBarry Smith    Notes:
1670ed56e8eaSBarry Smith    ISLocalToGlobalMappingRegister() may be called multiple times to add several user-defined mappings.
1671413f72f0SBarry Smith 
1672413f72f0SBarry Smith    Sample usage:
1673413f72f0SBarry Smith .vb
1674413f72f0SBarry Smith    ISLocalToGlobalMappingRegister("my_mapper",MyCreate);
1675413f72f0SBarry Smith .ve
1676413f72f0SBarry Smith 
1677ed56e8eaSBarry Smith    Then, your mapping can be chosen with the procedural interface via
1678413f72f0SBarry Smith $     ISLocalToGlobalMappingSetType(ltog,"my_mapper")
1679413f72f0SBarry Smith    or at runtime via the option
1680ed56e8eaSBarry Smith $     -islocaltoglobalmapping_type my_mapper
1681413f72f0SBarry Smith 
1682413f72f0SBarry Smith    Level: advanced
1683413f72f0SBarry Smith 
1684413f72f0SBarry Smith .keywords: ISLocalToGlobalMappingType, register
1685413f72f0SBarry Smith 
1686413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingRegisterAll(), ISLocalToGlobalMappingRegisterDestroy(), ISLOCALTOGLOBALMAPPINGBASIC, ISLOCALTOGLOBALMAPPINGHASH
1687413f72f0SBarry Smith 
1688413f72f0SBarry Smith @*/
1689413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRegister(const char sname[],PetscErrorCode (*function)(ISLocalToGlobalMapping))
1690413f72f0SBarry Smith {
1691413f72f0SBarry Smith   PetscErrorCode ierr;
1692413f72f0SBarry Smith 
1693413f72f0SBarry Smith   PetscFunctionBegin;
1694413f72f0SBarry Smith   ierr = PetscFunctionListAdd(&ISLocalToGlobalMappingList,sname,function);CHKERRQ(ierr);
1695413f72f0SBarry Smith   PetscFunctionReturn(0);
1696413f72f0SBarry Smith }
1697413f72f0SBarry Smith 
1698413f72f0SBarry Smith /*@C
1699ed56e8eaSBarry Smith    ISLocalToGlobalMappingSetType - Builds ISLocalToGlobalMapping for a particular global to local mapping approach.
1700413f72f0SBarry Smith 
1701413f72f0SBarry Smith    Logically Collective on ISLocalToGlobalMapping
1702413f72f0SBarry Smith 
1703413f72f0SBarry Smith    Input Parameters:
1704413f72f0SBarry Smith +  ltog - the ISLocalToGlobalMapping object
1705413f72f0SBarry Smith -  type - a known method
1706413f72f0SBarry Smith 
1707413f72f0SBarry Smith    Options Database Key:
1708ed56e8eaSBarry Smith .  -islocaltoglobalmapping_type  <method> - Sets the method; use -help for a list
1709413f72f0SBarry Smith     of available methods (for instance, basic or hash)
1710413f72f0SBarry Smith 
1711413f72f0SBarry Smith    Notes:
1712413f72f0SBarry Smith    See "petsc/include/petscis.h" for available methods
1713413f72f0SBarry Smith 
1714413f72f0SBarry Smith   Normally, it is best to use the ISLocalToGlobalMappingSetFromOptions() command and
1715413f72f0SBarry Smith   then set the ISLocalToGlobalMapping type from the options database rather than by using
1716413f72f0SBarry Smith   this routine.
1717413f72f0SBarry Smith 
1718413f72f0SBarry Smith   Level: intermediate
1719413f72f0SBarry Smith 
1720413f72f0SBarry Smith   Developer Note: ISLocalToGlobalMappingRegister() is used to add new types to ISLocalToGlobalMappingList from which they
1721413f72f0SBarry Smith   are accessed by ISLocalToGlobalMappingSetType().
1722413f72f0SBarry Smith 
1723413f72f0SBarry Smith .keywords: ISLocalToGlobalMapping, set, method
1724413f72f0SBarry Smith 
1725413f72f0SBarry Smith .seealso: PCSetType(), ISLocalToGlobalMappingType, ISLocalToGlobalMappingRegister(), ISLocalToGlobalMappingCreate()
1726413f72f0SBarry Smith 
1727413f72f0SBarry Smith @*/
1728413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType type)
1729413f72f0SBarry Smith {
1730413f72f0SBarry Smith   PetscErrorCode ierr,(*r)(ISLocalToGlobalMapping);
1731413f72f0SBarry Smith   PetscBool      match;
1732413f72f0SBarry Smith 
1733413f72f0SBarry Smith   PetscFunctionBegin;
1734413f72f0SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
1735413f72f0SBarry Smith   PetscValidCharPointer(type,2);
1736413f72f0SBarry Smith 
1737413f72f0SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)ltog,type,&match);CHKERRQ(ierr);
1738413f72f0SBarry Smith   if (match) PetscFunctionReturn(0);
1739413f72f0SBarry Smith 
1740413f72f0SBarry Smith   ierr =  PetscFunctionListFind(ISLocalToGlobalMappingList,type,&r);CHKERRQ(ierr);
1741413f72f0SBarry Smith   if (!r) SETERRQ1(PetscObjectComm((PetscObject)ltog),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested ISLocalToGlobalMapping type %s",type);
1742413f72f0SBarry Smith   /* Destroy the previous private LTOG context */
1743413f72f0SBarry Smith   if (ltog->ops->destroy) {
1744413f72f0SBarry Smith     ierr              = (*ltog->ops->destroy)(ltog);CHKERRQ(ierr);
1745413f72f0SBarry Smith     ltog->ops->destroy = NULL;
1746413f72f0SBarry Smith   }
1747413f72f0SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)ltog,type);CHKERRQ(ierr);
1748413f72f0SBarry Smith   ierr = (*r)(ltog);CHKERRQ(ierr);
1749413f72f0SBarry Smith   PetscFunctionReturn(0);
1750413f72f0SBarry Smith }
1751413f72f0SBarry Smith 
1752413f72f0SBarry Smith PetscBool ISLocalToGlobalMappingRegisterAllCalled = PETSC_FALSE;
1753413f72f0SBarry Smith 
1754413f72f0SBarry Smith /*@C
1755413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAll - Registers all of the local to global mapping components in the IS package.
1756413f72f0SBarry Smith 
1757413f72f0SBarry Smith   Not Collective
1758413f72f0SBarry Smith 
1759413f72f0SBarry Smith   Level: advanced
1760413f72f0SBarry Smith 
1761413f72f0SBarry Smith .keywords: IS, register, all
1762413f72f0SBarry Smith .seealso:  ISRegister(),  ISLocalToGlobalRegister()
1763413f72f0SBarry Smith @*/
1764413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRegisterAll(void)
1765413f72f0SBarry Smith {
1766413f72f0SBarry Smith   PetscErrorCode ierr;
1767413f72f0SBarry Smith 
1768413f72f0SBarry Smith   PetscFunctionBegin;
1769413f72f0SBarry Smith   if (ISLocalToGlobalMappingRegisterAllCalled) PetscFunctionReturn(0);
1770413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAllCalled = PETSC_TRUE;
1771413f72f0SBarry Smith 
1772413f72f0SBarry Smith   ierr = ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGBASIC, ISLocalToGlobalMappingCreate_Basic);CHKERRQ(ierr);
1773413f72f0SBarry Smith   ierr = ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGHASH, ISLocalToGlobalMappingCreate_Hash);CHKERRQ(ierr);
1774413f72f0SBarry Smith   PetscFunctionReturn(0);
1775413f72f0SBarry Smith }
177604a59952SBarry Smith 
1777