xref: /petsc/src/vec/is/utils/isltog.c (revision ed56e8eadf73fdc0691184ae99f1f70c6e63eeb7)
12362add9SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/isimpl.h>    /*I "petscis.h"  I*/
30c312b8eSJed Brown #include <petscsf.h>
4665c2dedSJed Brown #include <petscviewer.h>
52362add9SBarry Smith 
67087cfbeSBarry Smith PetscClassId IS_LTOGM_CLASSID;
7268a049cSStefano Zampini static PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping,PetscInt*,PetscInt**,PetscInt**,PetscInt***);
88e58c17dSMatthew Knepley 
9413f72f0SBarry Smith typedef struct {
10413f72f0SBarry Smith   PetscInt    *globals;
11413f72f0SBarry Smith } ISLocalToGlobalMapping_Basic;
12413f72f0SBarry Smith 
13413f72f0SBarry Smith typedef struct {
14413f72f0SBarry Smith   PetscHashI globalht;
15413f72f0SBarry Smith } ISLocalToGlobalMapping_Hash;
16413f72f0SBarry Smith 
17413f72f0SBarry Smith 
18413f72f0SBarry Smith /* -----------------------------------------------------------------------------------------*/
19413f72f0SBarry Smith 
20413f72f0SBarry Smith /*
21413f72f0SBarry Smith     Creates the global mapping information in the ISLocalToGlobalMapping structure
22413f72f0SBarry Smith 
23413f72f0SBarry Smith     If the user has not selected how to handle the global to local mapping then use HASH for "large" problems
24413f72f0SBarry Smith */
25413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp(ISLocalToGlobalMapping mapping)
26413f72f0SBarry Smith {
27413f72f0SBarry Smith   PetscInt       i,*idx = mapping->indices,n = mapping->n,end,start;
28413f72f0SBarry Smith   PetscErrorCode ierr;
29413f72f0SBarry Smith 
30413f72f0SBarry Smith   PetscFunctionBegin;
31413f72f0SBarry Smith   if (mapping->data) PetscFunctionReturn(0);
32413f72f0SBarry Smith   end   = 0;
33413f72f0SBarry Smith   start = PETSC_MAX_INT;
34413f72f0SBarry Smith 
35413f72f0SBarry Smith   for (i=0; i<n; i++) {
36413f72f0SBarry Smith     if (idx[i] < 0) continue;
37413f72f0SBarry Smith     if (idx[i] < start) start = idx[i];
38413f72f0SBarry Smith     if (idx[i] > end)   end   = idx[i];
39413f72f0SBarry Smith   }
40413f72f0SBarry Smith   if (start > end) {start = 0; end = -1;}
41413f72f0SBarry Smith   mapping->globalstart = start;
42413f72f0SBarry Smith   mapping->globalend   = end;
43413f72f0SBarry Smith   if (!((PetscObject)mapping)->type_name) {
44413f72f0SBarry Smith     if ((end - start) > PetscMax(4*n,1000000)) {
45413f72f0SBarry Smith       ierr = ISLocalToGlobalMappingSetType(mapping,ISLOCALTOGLOBALMAPPINGHASH);
46413f72f0SBarry Smith     } else {
47413f72f0SBarry Smith       ierr = ISLocalToGlobalMappingSetType(mapping,ISLOCALTOGLOBALMAPPINGBASIC);
48413f72f0SBarry Smith     }
49413f72f0SBarry Smith   }
50413f72f0SBarry Smith   ierr = (*mapping->ops->globaltolocalmappingsetup)(mapping);CHKERRQ(ierr);
51413f72f0SBarry Smith   PetscFunctionReturn(0);
52413f72f0SBarry Smith }
53413f72f0SBarry Smith 
54413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp_Basic(ISLocalToGlobalMapping mapping)
55413f72f0SBarry Smith {
56413f72f0SBarry Smith   PetscErrorCode              ierr;
57413f72f0SBarry Smith   PetscInt                    i,*idx = mapping->indices,n = mapping->n,end,start,*globals;
58413f72f0SBarry Smith   ISLocalToGlobalMapping_Basic *map;
59413f72f0SBarry Smith 
60413f72f0SBarry Smith   PetscFunctionBegin;
61413f72f0SBarry Smith   start            = mapping->globalstart;
62413f72f0SBarry Smith   end              = mapping->globalend;
63413f72f0SBarry Smith   ierr             = PetscNew(&map);CHKERRQ(ierr);
64413f72f0SBarry Smith   ierr             = PetscMalloc1(end-start+2,&globals);CHKERRQ(ierr);
65413f72f0SBarry Smith   map->globals     = globals;
66413f72f0SBarry Smith   for (i=0; i<end-start+1; i++) globals[i] = -1;
67413f72f0SBarry Smith   for (i=0; i<n; i++) {
68413f72f0SBarry Smith     if (idx[i] < 0) continue;
69413f72f0SBarry Smith     globals[idx[i] - start] = i;
70413f72f0SBarry Smith   }
71413f72f0SBarry Smith   mapping->data = (void*)map;
72413f72f0SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)mapping,(end-start+1)*sizeof(PetscInt));CHKERRQ(ierr);
73413f72f0SBarry Smith   PetscFunctionReturn(0);
74413f72f0SBarry Smith }
75413f72f0SBarry Smith 
76413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp_Hash(ISLocalToGlobalMapping mapping)
77413f72f0SBarry Smith {
78413f72f0SBarry Smith   PetscErrorCode              ierr;
79413f72f0SBarry Smith   PetscInt                    i,*idx = mapping->indices,n = mapping->n;
80413f72f0SBarry Smith   ISLocalToGlobalMapping_Hash *map;
81413f72f0SBarry Smith 
82413f72f0SBarry Smith   PetscFunctionBegin;
83413f72f0SBarry Smith   ierr = PetscNew(&map);CHKERRQ(ierr);
84413f72f0SBarry Smith   PetscHashICreate(map->globalht);
85413f72f0SBarry Smith   for (i=0; i<n; i++ ) {
86413f72f0SBarry Smith     if (idx[i] < 0) continue;
87413f72f0SBarry Smith     PetscHashIAdd(map->globalht, idx[i], i);
88413f72f0SBarry Smith   }
89413f72f0SBarry Smith   mapping->data = (void*)map;
90413f72f0SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)mapping,2*n*sizeof(PetscInt));CHKERRQ(ierr);
91413f72f0SBarry Smith   PetscFunctionReturn(0);
92413f72f0SBarry Smith }
93413f72f0SBarry Smith 
94413f72f0SBarry Smith static PetscErrorCode ISLocalToGlobalMappingDestroy_Basic(ISLocalToGlobalMapping mapping)
95413f72f0SBarry Smith {
96413f72f0SBarry Smith   PetscErrorCode              ierr;
97413f72f0SBarry Smith   ISLocalToGlobalMapping_Basic *map  = (ISLocalToGlobalMapping_Basic *)mapping->data;
98413f72f0SBarry Smith 
99413f72f0SBarry Smith   PetscFunctionBegin;
100413f72f0SBarry Smith   if (!map) PetscFunctionReturn(0);
101413f72f0SBarry Smith   ierr = PetscFree(map->globals);CHKERRQ(ierr);
102413f72f0SBarry Smith   ierr = PetscFree(mapping->data);CHKERRQ(ierr);
103413f72f0SBarry Smith   PetscFunctionReturn(0);
104413f72f0SBarry Smith }
105413f72f0SBarry Smith 
106413f72f0SBarry Smith static PetscErrorCode ISLocalToGlobalMappingDestroy_Hash(ISLocalToGlobalMapping mapping)
107413f72f0SBarry Smith {
108413f72f0SBarry Smith   PetscErrorCode              ierr;
109413f72f0SBarry Smith   ISLocalToGlobalMapping_Hash *map  = (ISLocalToGlobalMapping_Hash*)mapping->data;
110413f72f0SBarry Smith 
111413f72f0SBarry Smith   PetscFunctionBegin;
112413f72f0SBarry Smith   if (!map) PetscFunctionReturn(0);
113413f72f0SBarry Smith   PetscHashIDestroy(map->globalht);
114413f72f0SBarry Smith   ierr = PetscFree(mapping->data);CHKERRQ(ierr);
115413f72f0SBarry Smith   PetscFunctionReturn(0);
116413f72f0SBarry Smith }
117413f72f0SBarry Smith 
118413f72f0SBarry Smith #define GTOLTYPE _Basic
119413f72f0SBarry Smith #define GTOLNAME _Basic
120413f72f0SBarry Smith #define GTOL(g, local) do {                  \
121413f72f0SBarry Smith     local = map->globals[g/bs - start];      \
122413f72f0SBarry Smith     local = bs*local + (g % bs);             \
123413f72f0SBarry Smith   } while (0)
124413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
125413f72f0SBarry Smith 
126413f72f0SBarry Smith #define GTOLTYPE _Basic
127413f72f0SBarry Smith #define GTOLNAME Block_Basic
128413f72f0SBarry Smith #define GTOL(g, local) do {                  \
129413f72f0SBarry Smith     local = map->globals[g - start];         \
130413f72f0SBarry Smith   } while (0)
131413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
132413f72f0SBarry Smith 
133413f72f0SBarry Smith #define GTOLTYPE _Hash
134413f72f0SBarry Smith #define GTOLNAME _Hash
135413f72f0SBarry Smith #define GTOL(g, local) do {                  \
136413f72f0SBarry Smith     PetscHashIMap(map->globalht,g/bs,local); \
137413f72f0SBarry Smith     local = bs*local + (g % bs);             \
138413f72f0SBarry Smith    } while (0)
139413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
140413f72f0SBarry Smith 
141413f72f0SBarry Smith #define GTOLTYPE _Hash
142413f72f0SBarry Smith #define GTOLNAME Block_Hash
143413f72f0SBarry Smith #define GTOL(g, local) do {                  \
144413f72f0SBarry Smith     PetscHashIMap(map->globalht,g,local);    \
145413f72f0SBarry Smith   } while (0)
146413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
147413f72f0SBarry Smith 
1486658fb44Sstefano_zampini /*@
1496658fb44Sstefano_zampini     ISLocalToGlobalMappingDuplicate - Duplicates the local to global mapping object
1506658fb44Sstefano_zampini 
1516658fb44Sstefano_zampini     Not Collective
1526658fb44Sstefano_zampini 
1536658fb44Sstefano_zampini     Input Parameter:
1546658fb44Sstefano_zampini .   ltog - local to global mapping
1556658fb44Sstefano_zampini 
1566658fb44Sstefano_zampini     Output Parameter:
1576658fb44Sstefano_zampini .   nltog - the duplicated local to global mapping
1586658fb44Sstefano_zampini 
1596658fb44Sstefano_zampini     Level: advanced
1606658fb44Sstefano_zampini 
1616658fb44Sstefano_zampini     Concepts: mapping^local to global
1626658fb44Sstefano_zampini 
1636658fb44Sstefano_zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
1646658fb44Sstefano_zampini @*/
1656658fb44Sstefano_zampini PetscErrorCode  ISLocalToGlobalMappingDuplicate(ISLocalToGlobalMapping ltog,ISLocalToGlobalMapping* nltog)
1666658fb44Sstefano_zampini {
1676658fb44Sstefano_zampini   PetscErrorCode ierr;
1686658fb44Sstefano_zampini 
1696658fb44Sstefano_zampini   PetscFunctionBegin;
1706658fb44Sstefano_zampini   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
1716658fb44Sstefano_zampini   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)ltog),ltog->bs,ltog->n,ltog->indices,PETSC_COPY_VALUES,nltog);CHKERRQ(ierr);
1726658fb44Sstefano_zampini   PetscFunctionReturn(0);
1736658fb44Sstefano_zampini }
1746658fb44Sstefano_zampini 
175565245c5SBarry Smith /*@
176107e9a97SBarry Smith     ISLocalToGlobalMappingGetSize - Gets the local size of a local to global mapping
1773b9aefa3SBarry Smith 
1783b9aefa3SBarry Smith     Not Collective
1793b9aefa3SBarry Smith 
1803b9aefa3SBarry Smith     Input Parameter:
1813b9aefa3SBarry Smith .   ltog - local to global mapping
1823b9aefa3SBarry Smith 
1833b9aefa3SBarry Smith     Output Parameter:
184107e9a97SBarry Smith .   n - the number of entries in the local mapping, ISLocalToGlobalMappingGetIndices() returns an array of this length
1853b9aefa3SBarry Smith 
1863b9aefa3SBarry Smith     Level: advanced
1873b9aefa3SBarry Smith 
188273d9f13SBarry Smith     Concepts: mapping^local to global
1893b9aefa3SBarry Smith 
1903b9aefa3SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
1913b9aefa3SBarry Smith @*/
1927087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping mapping,PetscInt *n)
1933b9aefa3SBarry Smith {
1943b9aefa3SBarry Smith   PetscFunctionBegin;
1950700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1964482741eSBarry Smith   PetscValidIntPointer(n,2);
197107e9a97SBarry Smith   *n = mapping->bs*mapping->n;
1983b9aefa3SBarry Smith   PetscFunctionReturn(0);
1993b9aefa3SBarry Smith }
2003b9aefa3SBarry Smith 
2015a5d4f66SBarry Smith /*@C
2025a5d4f66SBarry Smith     ISLocalToGlobalMappingView - View a local to global mapping
2035a5d4f66SBarry Smith 
204b9cd556bSLois Curfman McInnes     Not Collective
205b9cd556bSLois Curfman McInnes 
2065a5d4f66SBarry Smith     Input Parameters:
2073b9aefa3SBarry Smith +   ltog - local to global mapping
2083b9aefa3SBarry Smith -   viewer - viewer
2095a5d4f66SBarry Smith 
210a997ad1aSLois Curfman McInnes     Level: advanced
211a997ad1aSLois Curfman McInnes 
212273d9f13SBarry Smith     Concepts: mapping^local to global
2135a5d4f66SBarry Smith 
2145a5d4f66SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
2155a5d4f66SBarry Smith @*/
2167087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingView(ISLocalToGlobalMapping mapping,PetscViewer viewer)
2175a5d4f66SBarry Smith {
21832dcc486SBarry Smith   PetscInt       i;
21932dcc486SBarry Smith   PetscMPIInt    rank;
220ace3abfcSBarry Smith   PetscBool      iascii;
2216849ba73SBarry Smith   PetscErrorCode ierr;
2225a5d4f66SBarry Smith 
2235a5d4f66SBarry Smith   PetscFunctionBegin;
2240700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
2253050cee2SBarry Smith   if (!viewer) {
226ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mapping),&viewer);CHKERRQ(ierr);
2273050cee2SBarry Smith   }
2280700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
2295a5d4f66SBarry Smith 
230ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mapping),&rank);CHKERRQ(ierr);
231251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
23232077d6dSBarry Smith   if (iascii) {
23398c3331eSBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)mapping,viewer);CHKERRQ(ierr);
2341575c14dSBarry Smith     ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
2355a5d4f66SBarry Smith     for (i=0; i<mapping->n; i++) {
2367904a332SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D %D\n",rank,i,mapping->indices[i]);CHKERRQ(ierr);
2376831982aSBarry Smith     }
238b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
2391575c14dSBarry Smith     ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
2401575c14dSBarry Smith   }
2415a5d4f66SBarry Smith   PetscFunctionReturn(0);
2425a5d4f66SBarry Smith }
2435a5d4f66SBarry Smith 
2441f428162SBarry Smith /*@
2452bdab257SBarry Smith     ISLocalToGlobalMappingCreateIS - Creates a mapping between a local (0 to n)
2462bdab257SBarry Smith     ordering and a global parallel ordering.
2472bdab257SBarry Smith 
2480f5bd95cSBarry Smith     Not collective
249b9cd556bSLois Curfman McInnes 
250a997ad1aSLois Curfman McInnes     Input Parameter:
2518c03b21aSDmitry Karpeev .   is - index set containing the global numbers for each local number
2522bdab257SBarry Smith 
253a997ad1aSLois Curfman McInnes     Output Parameter:
2542bdab257SBarry Smith .   mapping - new mapping data structure
2552bdab257SBarry Smith 
256f0413b6fSBarry Smith     Notes: the block size of the IS determines the block size of the mapping
257a997ad1aSLois Curfman McInnes     Level: advanced
258a997ad1aSLois Curfman McInnes 
259273d9f13SBarry Smith     Concepts: mapping^local to global
2602bdab257SBarry Smith 
2617e99dc12SLawrence Mitchell .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetFromOptions()
2622bdab257SBarry Smith @*/
2637087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingCreateIS(IS is,ISLocalToGlobalMapping *mapping)
2642bdab257SBarry Smith {
2656849ba73SBarry Smith   PetscErrorCode ierr;
2663bbf0e92SBarry Smith   PetscInt       n,bs;
2675d0c19d7SBarry Smith   const PetscInt *indices;
2682bdab257SBarry Smith   MPI_Comm       comm;
2693bbf0e92SBarry Smith   PetscBool      isblock;
2703a40ed3dSBarry Smith 
2713a40ed3dSBarry Smith   PetscFunctionBegin;
2720700a824SBarry Smith   PetscValidHeaderSpecific(is,IS_CLASSID,1);
2734482741eSBarry Smith   PetscValidPointer(mapping,2);
2742bdab257SBarry Smith 
2752bdab257SBarry Smith   ierr = PetscObjectGetComm((PetscObject)is,&comm);CHKERRQ(ierr);
2763b9aefa3SBarry Smith   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
2773bbf0e92SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)is,ISBLOCK,&isblock);CHKERRQ(ierr);
2786006e8d2SBarry Smith   if (!isblock) {
279f0413b6fSBarry Smith     ierr = ISGetIndices(is,&indices);CHKERRQ(ierr);
280f0413b6fSBarry Smith     ierr = ISLocalToGlobalMappingCreate(comm,1,n,indices,PETSC_COPY_VALUES,mapping);CHKERRQ(ierr);
2812bdab257SBarry Smith     ierr = ISRestoreIndices(is,&indices);CHKERRQ(ierr);
2826006e8d2SBarry Smith   } else {
2836006e8d2SBarry Smith     ierr = ISGetBlockSize(is,&bs);CHKERRQ(ierr);
284f0413b6fSBarry Smith     ierr = ISBlockGetIndices(is,&indices);CHKERRQ(ierr);
28528bc9809SBarry Smith     ierr = ISLocalToGlobalMappingCreate(comm,bs,n/bs,indices,PETSC_COPY_VALUES,mapping);CHKERRQ(ierr);
286f0413b6fSBarry Smith     ierr = ISBlockRestoreIndices(is,&indices);CHKERRQ(ierr);
2876006e8d2SBarry Smith   }
2883a40ed3dSBarry Smith   PetscFunctionReturn(0);
2892bdab257SBarry Smith }
2905a5d4f66SBarry Smith 
291a4d96a55SJed Brown /*@C
292a4d96a55SJed Brown     ISLocalToGlobalMappingCreateSF - Creates a mapping between a local (0 to n)
293a4d96a55SJed Brown     ordering and a global parallel ordering.
294a4d96a55SJed Brown 
295a4d96a55SJed Brown     Collective
296a4d96a55SJed Brown 
297a4d96a55SJed Brown     Input Parameter:
298a4d96a55SJed Brown +   sf - star forest mapping contiguous local indices to (rank, offset)
299a4d96a55SJed Brown -   start - first global index on this process
300a4d96a55SJed Brown 
301a4d96a55SJed Brown     Output Parameter:
302a4d96a55SJed Brown .   mapping - new mapping data structure
303a4d96a55SJed Brown 
304a4d96a55SJed Brown     Level: advanced
305a4d96a55SJed Brown 
306a4d96a55SJed Brown     Concepts: mapping^local to global
307a4d96a55SJed Brown 
3087e99dc12SLawrence Mitchell .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingSetFromOptions()
309a4d96a55SJed Brown @*/
310a4d96a55SJed Brown PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF sf,PetscInt start,ISLocalToGlobalMapping *mapping)
311a4d96a55SJed Brown {
312a4d96a55SJed Brown   PetscErrorCode ierr;
313a4d96a55SJed Brown   PetscInt       i,maxlocal,nroots,nleaves,*globals,*ltog;
314a4d96a55SJed Brown   const PetscInt *ilocal;
315a4d96a55SJed Brown   MPI_Comm       comm;
316a4d96a55SJed Brown 
317a4d96a55SJed Brown   PetscFunctionBegin;
318a4d96a55SJed Brown   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
319a4d96a55SJed Brown   PetscValidPointer(mapping,3);
320a4d96a55SJed Brown 
321a4d96a55SJed Brown   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
3220298fd71SBarry Smith   ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr);
323f6e5521dSKarl Rupp   if (ilocal) {
324f6e5521dSKarl Rupp     for (i=0,maxlocal=0; i<nleaves; i++) maxlocal = PetscMax(maxlocal,ilocal[i]+1);
325f6e5521dSKarl Rupp   }
326a4d96a55SJed Brown   else maxlocal = nleaves;
327785e854fSJed Brown   ierr = PetscMalloc1(nroots,&globals);CHKERRQ(ierr);
328785e854fSJed Brown   ierr = PetscMalloc1(maxlocal,&ltog);CHKERRQ(ierr);
329a4d96a55SJed Brown   for (i=0; i<nroots; i++) globals[i] = start + i;
330a4d96a55SJed Brown   for (i=0; i<maxlocal; i++) ltog[i] = -1;
331a4d96a55SJed Brown   ierr = PetscSFBcastBegin(sf,MPIU_INT,globals,ltog);CHKERRQ(ierr);
332a4d96a55SJed Brown   ierr = PetscSFBcastEnd(sf,MPIU_INT,globals,ltog);CHKERRQ(ierr);
333f0413b6fSBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,1,maxlocal,ltog,PETSC_OWN_POINTER,mapping);CHKERRQ(ierr);
334a4d96a55SJed Brown   ierr = PetscFree(globals);CHKERRQ(ierr);
335a4d96a55SJed Brown   PetscFunctionReturn(0);
336a4d96a55SJed Brown }
337b46b645bSBarry Smith 
33863fa5c83Sstefano_zampini /*@
33963fa5c83Sstefano_zampini     ISLocalToGlobalMappingSetBlockSize - Sets the blocksize of the mapping
34063fa5c83Sstefano_zampini 
34163fa5c83Sstefano_zampini     Not collective
34263fa5c83Sstefano_zampini 
34363fa5c83Sstefano_zampini     Input Parameters:
34463fa5c83Sstefano_zampini .   mapping - mapping data structure
34563fa5c83Sstefano_zampini .   bs - the blocksize
34663fa5c83Sstefano_zampini 
34763fa5c83Sstefano_zampini     Level: advanced
34863fa5c83Sstefano_zampini 
34963fa5c83Sstefano_zampini     Concepts: mapping^local to global
35063fa5c83Sstefano_zampini 
35163fa5c83Sstefano_zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS()
35263fa5c83Sstefano_zampini @*/
35363fa5c83Sstefano_zampini PetscErrorCode  ISLocalToGlobalMappingSetBlockSize(ISLocalToGlobalMapping mapping,PetscInt bs)
35463fa5c83Sstefano_zampini {
355a59f3c4dSstefano_zampini   PetscInt       *nid;
356a59f3c4dSstefano_zampini   const PetscInt *oid;
357a59f3c4dSstefano_zampini   PetscInt       i,cn,on,obs,nn;
35863fa5c83Sstefano_zampini   PetscErrorCode ierr;
35963fa5c83Sstefano_zampini 
36063fa5c83Sstefano_zampini   PetscFunctionBegin;
36163fa5c83Sstefano_zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
36263fa5c83Sstefano_zampini   if (bs < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid block size %D",bs);
36363fa5c83Sstefano_zampini   if (bs == mapping->bs) PetscFunctionReturn(0);
36463fa5c83Sstefano_zampini   on  = mapping->n;
36563fa5c83Sstefano_zampini   obs = mapping->bs;
36663fa5c83Sstefano_zampini   oid = mapping->indices;
36763fa5c83Sstefano_zampini   nn  = (on*obs)/bs;
36863fa5c83Sstefano_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);
369a59f3c4dSstefano_zampini 
37063fa5c83Sstefano_zampini   ierr = PetscMalloc1(nn,&nid);CHKERRQ(ierr);
371a59f3c4dSstefano_zampini   ierr = ISLocalToGlobalMappingGetIndices(mapping,&oid);CHKERRQ(ierr);
372a59f3c4dSstefano_zampini   for (i=0;i<nn;i++) {
373a59f3c4dSstefano_zampini     PetscInt j;
374a59f3c4dSstefano_zampini     for (j=0,cn=0;j<bs-1;j++) {
375a59f3c4dSstefano_zampini       if (oid[i*bs+j] < 0) { cn++; continue; }
376a59f3c4dSstefano_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]);
377a59f3c4dSstefano_zampini     }
378a59f3c4dSstefano_zampini     if (oid[i*bs+j] < 0) cn++;
3798b7cb0e6Sstefano_zampini     if (cn) {
380a59f3c4dSstefano_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);
381a59f3c4dSstefano_zampini       nid[i] = -1;
3828b7cb0e6Sstefano_zampini     } else {
383a59f3c4dSstefano_zampini       nid[i] = oid[i*bs]/bs;
38463fa5c83Sstefano_zampini     }
38563fa5c83Sstefano_zampini   }
386a59f3c4dSstefano_zampini   ierr = ISLocalToGlobalMappingRestoreIndices(mapping,&oid);CHKERRQ(ierr);
387a59f3c4dSstefano_zampini 
38863fa5c83Sstefano_zampini   mapping->n           = nn;
38963fa5c83Sstefano_zampini   mapping->bs          = bs;
39063fa5c83Sstefano_zampini   ierr                 = PetscFree(mapping->indices);CHKERRQ(ierr);
39163fa5c83Sstefano_zampini   mapping->indices     = nid;
392c9345713Sstefano_zampini   mapping->globalstart = 0;
393c9345713Sstefano_zampini   mapping->globalend   = 0;
394413f72f0SBarry Smith   if (mapping->ops->destroy) {
395413f72f0SBarry Smith     ierr = (*mapping->ops->destroy)(mapping);CHKERRQ(ierr);
396413f72f0SBarry Smith   }
39763fa5c83Sstefano_zampini   PetscFunctionReturn(0);
39863fa5c83Sstefano_zampini }
39963fa5c83Sstefano_zampini 
40045b6f7e9SBarry Smith /*@
40145b6f7e9SBarry Smith     ISLocalToGlobalMappingGetBlockSize - Gets the blocksize of the mapping
40245b6f7e9SBarry Smith     ordering and a global parallel ordering.
40345b6f7e9SBarry Smith 
40445b6f7e9SBarry Smith     Not Collective
40545b6f7e9SBarry Smith 
40645b6f7e9SBarry Smith     Input Parameters:
40745b6f7e9SBarry Smith .   mapping - mapping data structure
40845b6f7e9SBarry Smith 
40945b6f7e9SBarry Smith     Output Parameter:
41045b6f7e9SBarry Smith .   bs - the blocksize
41145b6f7e9SBarry Smith 
41245b6f7e9SBarry Smith     Level: advanced
41345b6f7e9SBarry Smith 
41445b6f7e9SBarry Smith     Concepts: mapping^local to global
41545b6f7e9SBarry Smith 
41645b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS()
41745b6f7e9SBarry Smith @*/
41845b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping mapping,PetscInt *bs)
41945b6f7e9SBarry Smith {
42045b6f7e9SBarry Smith   PetscFunctionBegin;
421cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
42245b6f7e9SBarry Smith   *bs = mapping->bs;
42345b6f7e9SBarry Smith   PetscFunctionReturn(0);
42445b6f7e9SBarry Smith }
42545b6f7e9SBarry Smith 
426ba5bb76aSSatish Balay /*@
42790f02eecSBarry Smith     ISLocalToGlobalMappingCreate - Creates a mapping between a local (0 to n)
42890f02eecSBarry Smith     ordering and a global parallel ordering.
4292362add9SBarry Smith 
43089d82c54SBarry Smith     Not Collective, but communicator may have more than one process
431b9cd556bSLois Curfman McInnes 
4322362add9SBarry Smith     Input Parameters:
43389d82c54SBarry Smith +   comm - MPI communicator
434f0413b6fSBarry Smith .   bs - the block size
43528bc9809SBarry Smith .   n - the number of local elements divided by the block size, or equivalently the number of block indices
43628bc9809SBarry 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
437d5ad8652SBarry Smith -   mode - see PetscCopyMode
4382362add9SBarry Smith 
439a997ad1aSLois Curfman McInnes     Output Parameter:
44090f02eecSBarry Smith .   mapping - new mapping data structure
4412362add9SBarry Smith 
442f0413b6fSBarry 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
443413f72f0SBarry Smith 
444413f72f0SBarry Smith     For "small" problems when using ISGlobalToMappingApply() and ISGlobalToMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
445413f72f0SBarry 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.
446413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
447413f72f0SBarry Smith 
448a997ad1aSLois Curfman McInnes     Level: advanced
449a997ad1aSLois Curfman McInnes 
450273d9f13SBarry Smith     Concepts: mapping^local to global
4512362add9SBarry Smith 
452413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingSetFromOptions(), ISLOCALTOGLOBALMAPPINGBASIC, ISLOCALTOGLOBALMAPPINGHASH
453413f72f0SBarry Smith           ISLocalToGlobalMappingSetType(), ISLocalToGlobalMappingType
4542362add9SBarry Smith @*/
45560c7cefcSBarry Smith PetscErrorCode  ISLocalToGlobalMappingCreate(MPI_Comm comm,PetscInt bs,PetscInt n,const PetscInt indices[],PetscCopyMode mode,ISLocalToGlobalMapping *mapping)
4562362add9SBarry Smith {
4576849ba73SBarry Smith   PetscErrorCode ierr;
45832dcc486SBarry Smith   PetscInt       *in;
459b46b645bSBarry Smith 
460b46b645bSBarry Smith   PetscFunctionBegin;
46173911063SBarry Smith   if (n) PetscValidIntPointer(indices,3);
4624482741eSBarry Smith   PetscValidPointer(mapping,4);
463b46b645bSBarry Smith 
4640298fd71SBarry Smith   *mapping = NULL;
465607a6623SBarry Smith   ierr = ISInitializePackage();CHKERRQ(ierr);
4662362add9SBarry Smith 
46773107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(*mapping,IS_LTOGM_CLASSID,"ISLocalToGlobalMapping","Local to global mapping","IS",
46860c7cefcSBarry Smith                            comm,ISLocalToGlobalMappingDestroy,ISLocalToGlobalMappingView);CHKERRQ(ierr);
469d4bb536fSBarry Smith   (*mapping)->n             = n;
470f0413b6fSBarry Smith   (*mapping)->bs            = bs;
471268a049cSStefano Zampini   (*mapping)->info_cached   = PETSC_FALSE;
472268a049cSStefano Zampini   (*mapping)->info_free     = PETSC_FALSE;
473268a049cSStefano Zampini   (*mapping)->info_procs    = NULL;
474268a049cSStefano Zampini   (*mapping)->info_numprocs = NULL;
475268a049cSStefano Zampini   (*mapping)->info_indices  = NULL;
476413f72f0SBarry Smith 
477413f72f0SBarry Smith   (*mapping)->ops->globaltolocalmappingapply      = NULL;
478413f72f0SBarry Smith   (*mapping)->ops->globaltolocalmappingapplyblock = NULL;
479413f72f0SBarry Smith   (*mapping)->ops->destroy                        = NULL;
480d5ad8652SBarry Smith   if (mode == PETSC_COPY_VALUES) {
481785e854fSJed Brown     ierr = PetscMalloc1(n,&in);CHKERRQ(ierr);
482d5ad8652SBarry Smith     ierr = PetscMemcpy(in,indices,n*sizeof(PetscInt));CHKERRQ(ierr);
483d5ad8652SBarry Smith     (*mapping)->indices = in;
4846389a1a1SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt));CHKERRQ(ierr);
4856389a1a1SBarry Smith   } else if (mode == PETSC_OWN_POINTER) {
4866389a1a1SBarry Smith     (*mapping)->indices = (PetscInt*)indices;
4876389a1a1SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt));CHKERRQ(ierr);
4886389a1a1SBarry Smith   }
48960c7cefcSBarry Smith   else SETERRQ(comm,PETSC_ERR_SUP,"Cannot currently use PETSC_USE_POINTER");
4903a40ed3dSBarry Smith   PetscFunctionReturn(0);
4912362add9SBarry Smith }
4922362add9SBarry Smith 
493413f72f0SBarry Smith PetscFunctionList ISLocalToGlobalMappingList = NULL;
494413f72f0SBarry Smith 
495413f72f0SBarry Smith 
49690f02eecSBarry Smith /*@
4977e99dc12SLawrence Mitchell    ISLocalToGlobalMappingSetFromOptions - Set mapping options from the options database.
4987e99dc12SLawrence Mitchell 
4997e99dc12SLawrence Mitchell    Not collective
5007e99dc12SLawrence Mitchell 
5017e99dc12SLawrence Mitchell    Input Parameters:
5027e99dc12SLawrence Mitchell .  mapping - mapping data structure
5037e99dc12SLawrence Mitchell 
5047e99dc12SLawrence Mitchell    Level: advanced
5057e99dc12SLawrence Mitchell 
5067e99dc12SLawrence Mitchell @*/
5077e99dc12SLawrence Mitchell PetscErrorCode ISLocalToGlobalMappingSetFromOptions(ISLocalToGlobalMapping mapping)
5087e99dc12SLawrence Mitchell {
5097e99dc12SLawrence Mitchell   PetscErrorCode             ierr;
510413f72f0SBarry Smith   char                       type[256];
511413f72f0SBarry Smith   ISLocalToGlobalMappingType defaulttype = "Not set";
5127e99dc12SLawrence Mitchell   PetscBool                  flg;
5137e99dc12SLawrence Mitchell 
5147e99dc12SLawrence Mitchell   PetscFunctionBegin;
5157e99dc12SLawrence Mitchell   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
516413f72f0SBarry Smith   ierr = ISLocalToGlobalMappingRegisterAll();CHKERRQ(ierr);
5177e99dc12SLawrence Mitchell   ierr = PetscObjectOptionsBegin((PetscObject)mapping);CHKERRQ(ierr);
518413f72f0SBarry Smith   ierr = PetscOptionsFList("-islocaltoglobalmapping_type","ISLocalToGlobalMapping method","ISLocalToGlobalMappingSetType",ISLocalToGlobalMappingList,(char*)(((PetscObject)mapping)->type_name) ? ((PetscObject)mapping)->type_name : defaulttype,type,256,&flg);CHKERRQ(ierr);
519413f72f0SBarry Smith   if (flg) {
520413f72f0SBarry Smith     ierr = ISLocalToGlobalMappingSetType(mapping,type);CHKERRQ(ierr);
521413f72f0SBarry Smith   }
5227e99dc12SLawrence Mitchell   ierr = PetscOptionsEnd();CHKERRQ(ierr);
5237e99dc12SLawrence Mitchell   PetscFunctionReturn(0);
5247e99dc12SLawrence Mitchell }
5257e99dc12SLawrence Mitchell 
5267e99dc12SLawrence Mitchell /*@
52790f02eecSBarry Smith    ISLocalToGlobalMappingDestroy - Destroys a mapping between a local (0 to n)
52890f02eecSBarry Smith    ordering and a global parallel ordering.
52990f02eecSBarry Smith 
5300f5bd95cSBarry Smith    Note Collective
531b9cd556bSLois Curfman McInnes 
53290f02eecSBarry Smith    Input Parameters:
53390f02eecSBarry Smith .  mapping - mapping data structure
53490f02eecSBarry Smith 
535a997ad1aSLois Curfman McInnes    Level: advanced
536a997ad1aSLois Curfman McInnes 
5373acfe500SLois Curfman McInnes .seealso: ISLocalToGlobalMappingCreate()
53890f02eecSBarry Smith @*/
5396bf464f9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping *mapping)
54090f02eecSBarry Smith {
541dfbe8321SBarry Smith   PetscErrorCode ierr;
5425fd66863SKarl Rupp 
5433a40ed3dSBarry Smith   PetscFunctionBegin;
5446bf464f9SBarry Smith   if (!*mapping) PetscFunctionReturn(0);
5456bf464f9SBarry Smith   PetscValidHeaderSpecific((*mapping),IS_LTOGM_CLASSID,1);
546997056adSBarry Smith   if (--((PetscObject)(*mapping))->refct > 0) {*mapping = 0;PetscFunctionReturn(0);}
5476bf464f9SBarry Smith   ierr = PetscFree((*mapping)->indices);CHKERRQ(ierr);
548268a049cSStefano Zampini   ierr = PetscFree((*mapping)->info_procs);CHKERRQ(ierr);
549268a049cSStefano Zampini   ierr = PetscFree((*mapping)->info_numprocs);CHKERRQ(ierr);
550268a049cSStefano Zampini   if ((*mapping)->info_indices) {
551268a049cSStefano Zampini     PetscInt i;
552268a049cSStefano Zampini 
553268a049cSStefano Zampini     ierr = PetscFree(((*mapping)->info_indices)[0]);CHKERRQ(ierr);
554268a049cSStefano Zampini     for (i=1; i<(*mapping)->info_nproc; i++) {
555268a049cSStefano Zampini       ierr = PetscFree(((*mapping)->info_indices)[i]);CHKERRQ(ierr);
556268a049cSStefano Zampini     }
557268a049cSStefano Zampini     ierr = PetscFree((*mapping)->info_indices);CHKERRQ(ierr);
558268a049cSStefano Zampini   }
559413f72f0SBarry Smith   if ((*mapping)->ops->destroy) {
560413f72f0SBarry Smith     ierr = (*(*mapping)->ops->destroy)(*mapping);CHKERRQ(ierr);
561413f72f0SBarry Smith   }
562d38fa0fbSBarry Smith   ierr     = PetscHeaderDestroy(mapping);CHKERRQ(ierr);
563992144d0SBarry Smith   *mapping = 0;
5643a40ed3dSBarry Smith   PetscFunctionReturn(0);
56590f02eecSBarry Smith }
56690f02eecSBarry Smith 
56790f02eecSBarry Smith /*@
5683acfe500SLois Curfman McInnes     ISLocalToGlobalMappingApplyIS - Creates from an IS in the local numbering
5693acfe500SLois Curfman McInnes     a new index set using the global numbering defined in an ISLocalToGlobalMapping
5703acfe500SLois Curfman McInnes     context.
57190f02eecSBarry Smith 
572b9cd556bSLois Curfman McInnes     Not collective
573b9cd556bSLois Curfman McInnes 
57490f02eecSBarry Smith     Input Parameters:
575b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
576b9cd556bSLois Curfman McInnes -   is - index set in local numbering
57790f02eecSBarry Smith 
57890f02eecSBarry Smith     Output Parameters:
57990f02eecSBarry Smith .   newis - index set in global numbering
58090f02eecSBarry Smith 
581a997ad1aSLois Curfman McInnes     Level: advanced
582a997ad1aSLois Curfman McInnes 
583273d9f13SBarry Smith     Concepts: mapping^local to global
5843acfe500SLois Curfman McInnes 
58590f02eecSBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),
586d4bb536fSBarry Smith           ISLocalToGlobalMappingDestroy(), ISGlobalToLocalMappingApply()
58790f02eecSBarry Smith @*/
5887087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping,IS is,IS *newis)
58990f02eecSBarry Smith {
5906849ba73SBarry Smith   PetscErrorCode ierr;
591e24637baSBarry Smith   PetscInt       n,*idxout;
5925d0c19d7SBarry Smith   const PetscInt *idxin;
5933a40ed3dSBarry Smith 
5943a40ed3dSBarry Smith   PetscFunctionBegin;
5950700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
5960700a824SBarry Smith   PetscValidHeaderSpecific(is,IS_CLASSID,2);
5974482741eSBarry Smith   PetscValidPointer(newis,3);
59890f02eecSBarry Smith 
5993b9aefa3SBarry Smith   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
60090f02eecSBarry Smith   ierr = ISGetIndices(is,&idxin);CHKERRQ(ierr);
601785e854fSJed Brown   ierr = PetscMalloc1(n,&idxout);CHKERRQ(ierr);
602e24637baSBarry Smith   ierr = ISLocalToGlobalMappingApply(mapping,n,idxin,idxout);CHKERRQ(ierr);
6033b9aefa3SBarry Smith   ierr = ISRestoreIndices(is,&idxin);CHKERRQ(ierr);
604543f3098SMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),n,idxout,PETSC_OWN_POINTER,newis);CHKERRQ(ierr);
6053a40ed3dSBarry Smith   PetscFunctionReturn(0);
60690f02eecSBarry Smith }
60790f02eecSBarry Smith 
608b89cb25eSSatish Balay /*@
6093acfe500SLois Curfman McInnes    ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering
6103acfe500SLois Curfman McInnes    and converts them to the global numbering.
61190f02eecSBarry Smith 
612b9cd556bSLois Curfman McInnes    Not collective
613b9cd556bSLois Curfman McInnes 
614bb25748dSBarry Smith    Input Parameters:
615b9cd556bSLois Curfman McInnes +  mapping - the local to global mapping context
616bb25748dSBarry Smith .  N - number of integers
617b9cd556bSLois Curfman McInnes -  in - input indices in local numbering
618bb25748dSBarry Smith 
619bb25748dSBarry Smith    Output Parameter:
620bb25748dSBarry Smith .  out - indices in global numbering
621bb25748dSBarry Smith 
622b9cd556bSLois Curfman McInnes    Notes:
623b9cd556bSLois Curfman McInnes    The in and out array parameters may be identical.
624d4bb536fSBarry Smith 
625a997ad1aSLois Curfman McInnes    Level: advanced
626a997ad1aSLois Curfman McInnes 
62745b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
6280752156aSBarry Smith           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
629d4bb536fSBarry Smith           AOPetscToApplication(), ISGlobalToLocalMappingApply()
630bb25748dSBarry Smith 
631273d9f13SBarry Smith     Concepts: mapping^local to global
632afcb2eb5SJed Brown @*/
633afcb2eb5SJed Brown PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
634afcb2eb5SJed Brown {
635cbc1caf0SMatthew G. Knepley   PetscInt i,bs,Nmax;
63645b6f7e9SBarry Smith 
63745b6f7e9SBarry Smith   PetscFunctionBegin;
638cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
639cbc1caf0SMatthew G. Knepley   bs   = mapping->bs;
640cbc1caf0SMatthew G. Knepley   Nmax = bs*mapping->n;
64145b6f7e9SBarry Smith   if (bs == 1) {
642cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
64345b6f7e9SBarry Smith     for (i=0; i<N; i++) {
64445b6f7e9SBarry Smith       if (in[i] < 0) {
64545b6f7e9SBarry Smith         out[i] = in[i];
64645b6f7e9SBarry Smith         continue;
64745b6f7e9SBarry Smith       }
648e24637baSBarry 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);
64945b6f7e9SBarry Smith       out[i] = idx[in[i]];
65045b6f7e9SBarry Smith     }
65145b6f7e9SBarry Smith   } else {
652cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
65345b6f7e9SBarry Smith     for (i=0; i<N; i++) {
65445b6f7e9SBarry Smith       if (in[i] < 0) {
65545b6f7e9SBarry Smith         out[i] = in[i];
65645b6f7e9SBarry Smith         continue;
65745b6f7e9SBarry Smith       }
658e24637baSBarry 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);
65945b6f7e9SBarry Smith       out[i] = idx[in[i]/bs]*bs + (in[i] % bs);
66045b6f7e9SBarry Smith     }
66145b6f7e9SBarry Smith   }
66245b6f7e9SBarry Smith   PetscFunctionReturn(0);
66345b6f7e9SBarry Smith }
66445b6f7e9SBarry Smith 
66545b6f7e9SBarry Smith /*@
6666006e8d2SBarry Smith    ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local block numbering and converts them to the global block numbering
66745b6f7e9SBarry Smith 
66845b6f7e9SBarry Smith    Not collective
66945b6f7e9SBarry Smith 
67045b6f7e9SBarry Smith    Input Parameters:
67145b6f7e9SBarry Smith +  mapping - the local to global mapping context
67245b6f7e9SBarry Smith .  N - number of integers
6736006e8d2SBarry Smith -  in - input indices in local block numbering
67445b6f7e9SBarry Smith 
67545b6f7e9SBarry Smith    Output Parameter:
6766006e8d2SBarry Smith .  out - indices in global block numbering
67745b6f7e9SBarry Smith 
67845b6f7e9SBarry Smith    Notes:
67945b6f7e9SBarry Smith    The in and out array parameters may be identical.
68045b6f7e9SBarry Smith 
6816006e8d2SBarry Smith    Example:
6826006e8d2SBarry 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
6836006e8d2SBarry Smith      (the first block) would produce 0 and the mapping applied to 1 (the second block) would produce 3.
6846006e8d2SBarry Smith 
68545b6f7e9SBarry Smith    Level: advanced
68645b6f7e9SBarry Smith 
68745b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
68845b6f7e9SBarry Smith           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
68945b6f7e9SBarry Smith           AOPetscToApplication(), ISGlobalToLocalMappingApply()
69045b6f7e9SBarry Smith 
69145b6f7e9SBarry Smith     Concepts: mapping^local to global
69245b6f7e9SBarry Smith @*/
69345b6f7e9SBarry Smith PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
69445b6f7e9SBarry Smith {
695cbc1caf0SMatthew G. Knepley 
696cbc1caf0SMatthew G. Knepley   PetscFunctionBegin;
697cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
698cbc1caf0SMatthew G. Knepley   {
699afcb2eb5SJed Brown     PetscInt i,Nmax = mapping->n;
700afcb2eb5SJed Brown     const PetscInt *idx = mapping->indices;
701d4bb536fSBarry Smith 
702afcb2eb5SJed Brown     for (i=0; i<N; i++) {
703afcb2eb5SJed Brown       if (in[i] < 0) {
704afcb2eb5SJed Brown         out[i] = in[i];
705afcb2eb5SJed Brown         continue;
706afcb2eb5SJed Brown       }
707e24637baSBarry 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);
708afcb2eb5SJed Brown       out[i] = idx[in[i]];
709afcb2eb5SJed Brown     }
710cbc1caf0SMatthew G. Knepley   }
711afcb2eb5SJed Brown   PetscFunctionReturn(0);
712afcb2eb5SJed Brown }
713d4bb536fSBarry Smith 
7147e99dc12SLawrence Mitchell /*@
715a997ad1aSLois Curfman McInnes     ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers
716a997ad1aSLois Curfman McInnes     specified with a global numbering.
717d4bb536fSBarry Smith 
718b9cd556bSLois Curfman McInnes     Not collective
719b9cd556bSLois Curfman McInnes 
720d4bb536fSBarry Smith     Input Parameters:
721b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
722d4bb536fSBarry Smith .   type - IS_GTOLM_MASK - replaces global indices with no local value with -1
723d4bb536fSBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
724d4bb536fSBarry Smith .   n - number of global indices to map
725b9cd556bSLois Curfman McInnes -   idx - global indices to map
726d4bb536fSBarry Smith 
727d4bb536fSBarry Smith     Output Parameters:
728b9cd556bSLois Curfman McInnes +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
729b9cd556bSLois Curfman McInnes -   idxout - local index of each global index, one must pass in an array long enough
730e182c471SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApply() with
7310298fd71SBarry Smith              idxout == NULL to determine the required length (returned in nout)
732e182c471SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApply()
733e182c471SBarry Smith              a second time to set the values.
734d4bb536fSBarry Smith 
735b9cd556bSLois Curfman McInnes     Notes:
7360298fd71SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
737d4bb536fSBarry Smith 
738413f72f0SBarry Smith     For "small" problems when using ISGlobalToMappingApply() and ISGlobalToMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
739413f72f0SBarry 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.
740413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
7410f5bd95cSBarry Smith 
742a997ad1aSLois Curfman McInnes     Level: advanced
743a997ad1aSLois Curfman McInnes 
74432fd6b96SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
74532fd6b96SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
74632fd6b96SBarry Smith 
747273d9f13SBarry Smith     Concepts: mapping^global to local
748d4bb536fSBarry Smith 
7499d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),
750413f72f0SBarry Smith           ISLocalToGlobalMappingDestroy()
751d4bb536fSBarry Smith @*/
752413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
753d4bb536fSBarry Smith {
7549d90f715SBarry Smith   PetscErrorCode ierr;
7559d90f715SBarry Smith 
7569d90f715SBarry Smith   PetscFunctionBegin;
7579d90f715SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
758413f72f0SBarry Smith   if (!mapping->data) {
759413f72f0SBarry Smith     ierr = ISGlobalToLocalMappingSetUp(mapping);CHKERRQ(ierr);
7609d90f715SBarry Smith   }
761413f72f0SBarry Smith   ierr = (*mapping->ops->globaltolocalmappingapply)(mapping,type,n,idx,nout,idxout);CHKERRQ(ierr);
7629d90f715SBarry Smith   PetscFunctionReturn(0);
7639d90f715SBarry Smith }
7649d90f715SBarry Smith 
765d4fe737eSStefano Zampini /*@
766d4fe737eSStefano Zampini     ISGlobalToLocalMappingApplyIS - Creates from an IS in the global numbering
767d4fe737eSStefano Zampini     a new index set using the local numbering defined in an ISLocalToGlobalMapping
768d4fe737eSStefano Zampini     context.
769d4fe737eSStefano Zampini 
770d4fe737eSStefano Zampini     Not collective
771d4fe737eSStefano Zampini 
772d4fe737eSStefano Zampini     Input Parameters:
773d4fe737eSStefano Zampini +   mapping - mapping between local and global numbering
774d4fe737eSStefano Zampini -   is - index set in global numbering
775d4fe737eSStefano Zampini 
776d4fe737eSStefano Zampini     Output Parameters:
777d4fe737eSStefano Zampini .   newis - index set in local numbering
778d4fe737eSStefano Zampini 
779d4fe737eSStefano Zampini     Level: advanced
780d4fe737eSStefano Zampini 
781d4fe737eSStefano Zampini     Concepts: mapping^local to global
782d4fe737eSStefano Zampini 
783d4fe737eSStefano Zampini .seealso: ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(),
784d4fe737eSStefano Zampini           ISLocalToGlobalMappingDestroy()
785d4fe737eSStefano Zampini @*/
786413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type, IS is,IS *newis)
787d4fe737eSStefano Zampini {
788d4fe737eSStefano Zampini   PetscErrorCode ierr;
789d4fe737eSStefano Zampini   PetscInt       n,nout,*idxout;
790d4fe737eSStefano Zampini   const PetscInt *idxin;
791d4fe737eSStefano Zampini 
792d4fe737eSStefano Zampini   PetscFunctionBegin;
793d4fe737eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
794d4fe737eSStefano Zampini   PetscValidHeaderSpecific(is,IS_CLASSID,3);
795d4fe737eSStefano Zampini   PetscValidPointer(newis,4);
796d4fe737eSStefano Zampini 
797d4fe737eSStefano Zampini   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
798d4fe737eSStefano Zampini   ierr = ISGetIndices(is,&idxin);CHKERRQ(ierr);
799d4fe737eSStefano Zampini   if (type == IS_GTOLM_MASK) {
800d4fe737eSStefano Zampini     ierr = PetscMalloc1(n,&idxout);CHKERRQ(ierr);
801d4fe737eSStefano Zampini   } else {
802d4fe737eSStefano Zampini     ierr = ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,NULL);CHKERRQ(ierr);
803d4fe737eSStefano Zampini     ierr = PetscMalloc1(nout,&idxout);CHKERRQ(ierr);
804d4fe737eSStefano Zampini   }
805d4fe737eSStefano Zampini   ierr = ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,idxout);CHKERRQ(ierr);
806d4fe737eSStefano Zampini   ierr = ISRestoreIndices(is,&idxin);CHKERRQ(ierr);
807d4fe737eSStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),nout,idxout,PETSC_OWN_POINTER,newis);CHKERRQ(ierr);
808d4fe737eSStefano Zampini   PetscFunctionReturn(0);
809d4fe737eSStefano Zampini }
810d4fe737eSStefano Zampini 
8119d90f715SBarry Smith /*@
8129d90f715SBarry Smith     ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers
8139d90f715SBarry Smith     specified with a block global numbering.
8149d90f715SBarry Smith 
8159d90f715SBarry Smith     Not collective
8169d90f715SBarry Smith 
8179d90f715SBarry Smith     Input Parameters:
8189d90f715SBarry Smith +   mapping - mapping between local and global numbering
8199d90f715SBarry Smith .   type - IS_GTOLM_MASK - replaces global indices with no local value with -1
8209d90f715SBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
8219d90f715SBarry Smith .   n - number of global indices to map
8229d90f715SBarry Smith -   idx - global indices to map
8239d90f715SBarry Smith 
8249d90f715SBarry Smith     Output Parameters:
8259d90f715SBarry Smith +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
8269d90f715SBarry Smith -   idxout - local index of each global index, one must pass in an array long enough
8279d90f715SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApplyBlock() with
8289d90f715SBarry Smith              idxout == NULL to determine the required length (returned in nout)
8299d90f715SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApplyBlock()
8309d90f715SBarry Smith              a second time to set the values.
8319d90f715SBarry Smith 
8329d90f715SBarry Smith     Notes:
8339d90f715SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
8349d90f715SBarry Smith 
835413f72f0SBarry Smith     For "small" problems when using ISGlobalToMappingApply() and ISGlobalToMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
836413f72f0SBarry 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.
837413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
8389d90f715SBarry Smith 
8399d90f715SBarry Smith     Level: advanced
8409d90f715SBarry Smith 
8419d90f715SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
8429d90f715SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
8439d90f715SBarry Smith 
8449d90f715SBarry Smith     Concepts: mapping^global to local
8459d90f715SBarry Smith 
8469d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(),
847413f72f0SBarry Smith           ISLocalToGlobalMappingDestroy()
8489d90f715SBarry Smith @*/
849413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,
8509d90f715SBarry Smith                                                  PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
8519d90f715SBarry Smith {
8526849ba73SBarry Smith   PetscErrorCode ierr;
853d4bb536fSBarry Smith 
8543a40ed3dSBarry Smith   PetscFunctionBegin;
8550700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
856413f72f0SBarry Smith   if (!mapping->data) {
857413f72f0SBarry Smith     ierr = ISGlobalToLocalMappingSetUp(mapping);CHKERRQ(ierr);
858d4bb536fSBarry Smith   }
859413f72f0SBarry Smith   ierr = (*mapping->ops->globaltolocalmappingapplyblock)(mapping,type,n,idx,nout,idxout);CHKERRQ(ierr);
8603a40ed3dSBarry Smith   PetscFunctionReturn(0);
861d4bb536fSBarry Smith }
86290f02eecSBarry Smith 
863413f72f0SBarry Smith 
86489d82c54SBarry Smith /*@C
8656a818285SBarry Smith     ISLocalToGlobalMappingGetBlockInfo - Gets the neighbor information for each processor and
86689d82c54SBarry Smith      each index shared by more than one processor
86789d82c54SBarry Smith 
86889d82c54SBarry Smith     Collective on ISLocalToGlobalMapping
86989d82c54SBarry Smith 
87089d82c54SBarry Smith     Input Parameters:
87189d82c54SBarry Smith .   mapping - the mapping from local to global indexing
87289d82c54SBarry Smith 
87389d82c54SBarry Smith     Output Parameter:
87489d82c54SBarry Smith +   nproc - number of processors that are connected to this one
87589d82c54SBarry Smith .   proc - neighboring processors
87607b52d57SBarry Smith .   numproc - number of indices for each subdomain (processor)
8773463a7baSJed Brown -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
87889d82c54SBarry Smith 
87989d82c54SBarry Smith     Level: advanced
88089d82c54SBarry Smith 
881273d9f13SBarry Smith     Concepts: mapping^local to global
88289d82c54SBarry Smith 
8832cfcea29SBarry Smith     Fortran Usage:
8842cfcea29SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
8852cfcea29SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
8862cfcea29SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
8872cfcea29SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
8882cfcea29SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
8892cfcea29SBarry Smith 
8902cfcea29SBarry Smith 
89107b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
89207b52d57SBarry Smith           ISLocalToGlobalMappingRestoreInfo()
89389d82c54SBarry Smith @*/
8946a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
89589d82c54SBarry Smith {
8966849ba73SBarry Smith   PetscErrorCode ierr;
897268a049cSStefano Zampini 
898268a049cSStefano Zampini   PetscFunctionBegin;
899268a049cSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
900268a049cSStefano Zampini   if (mapping->info_cached) {
901268a049cSStefano Zampini     *nproc    = mapping->info_nproc;
902268a049cSStefano Zampini     *procs    = mapping->info_procs;
903268a049cSStefano Zampini     *numprocs = mapping->info_numprocs;
904268a049cSStefano Zampini     *indices  = mapping->info_indices;
905268a049cSStefano Zampini   } else {
906268a049cSStefano Zampini     ierr = ISLocalToGlobalMappingGetBlockInfo_Private(mapping,nproc,procs,numprocs,indices);CHKERRQ(ierr);
907268a049cSStefano Zampini   }
908268a049cSStefano Zampini   PetscFunctionReturn(0);
909268a049cSStefano Zampini }
910268a049cSStefano Zampini 
911268a049cSStefano Zampini static PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
912268a049cSStefano Zampini {
913268a049cSStefano Zampini   PetscErrorCode ierr;
91497f1f81fSBarry Smith   PetscMPIInt    size,rank,tag1,tag2,tag3,*len,*source,imdex;
91532dcc486SBarry Smith   PetscInt       i,n = mapping->n,Ng,ng,max = 0,*lindices = mapping->indices;
91632dcc486SBarry Smith   PetscInt       *nprocs,*owner,nsends,*sends,j,*starts,nmax,nrecvs,*recvs,proc;
91797f1f81fSBarry Smith   PetscInt       cnt,scale,*ownedsenders,*nownedsenders,rstart,nowned;
91832dcc486SBarry Smith   PetscInt       node,nownedm,nt,*sends2,nsends2,*starts2,*lens2,*dest,nrecvs2,*starts3,*recvs2,k,*bprocs,*tmp;
91932dcc486SBarry Smith   PetscInt       first_procs,first_numprocs,*first_indices;
92089d82c54SBarry Smith   MPI_Request    *recv_waits,*send_waits;
92130dcb7c9SBarry Smith   MPI_Status     recv_status,*send_status,*recv_statuses;
922ce94432eSBarry Smith   MPI_Comm       comm;
923ace3abfcSBarry Smith   PetscBool      debug = PETSC_FALSE;
92489d82c54SBarry Smith 
92589d82c54SBarry Smith   PetscFunctionBegin;
926ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)mapping,&comm);CHKERRQ(ierr);
92724cf384cSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
92824cf384cSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
92924cf384cSBarry Smith   if (size == 1) {
93024cf384cSBarry Smith     *nproc         = 0;
9310298fd71SBarry Smith     *procs         = NULL;
93295dccacaSBarry Smith     ierr           = PetscNew(numprocs);CHKERRQ(ierr);
9331e2105dcSBarry Smith     (*numprocs)[0] = 0;
93495dccacaSBarry Smith     ierr           = PetscNew(indices);CHKERRQ(ierr);
9350298fd71SBarry Smith     (*indices)[0]  = NULL;
936268a049cSStefano Zampini     /* save info for reuse */
937268a049cSStefano Zampini     mapping->info_nproc = *nproc;
938268a049cSStefano Zampini     mapping->info_procs = *procs;
939268a049cSStefano Zampini     mapping->info_numprocs = *numprocs;
940268a049cSStefano Zampini     mapping->info_indices = *indices;
941268a049cSStefano Zampini     mapping->info_cached = PETSC_TRUE;
94224cf384cSBarry Smith     PetscFunctionReturn(0);
94324cf384cSBarry Smith   }
94424cf384cSBarry Smith 
945c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)mapping)->options,NULL,"-islocaltoglobalmappinggetinfo_debug",&debug,NULL);CHKERRQ(ierr);
94607b52d57SBarry Smith 
9473677ff5aSBarry Smith   /*
9486a818285SBarry Smith     Notes on ISLocalToGlobalMappingGetBlockInfo
9493677ff5aSBarry Smith 
9503677ff5aSBarry Smith     globally owned node - the nodes that have been assigned to this processor in global
9513677ff5aSBarry Smith            numbering, just for this routine.
9523677ff5aSBarry Smith 
9533677ff5aSBarry Smith     nontrivial globally owned node - node assigned to this processor that is on a subdomain
9543677ff5aSBarry Smith            boundary (i.e. is has more than one local owner)
9553677ff5aSBarry Smith 
9563677ff5aSBarry Smith     locally owned node - node that exists on this processors subdomain
9573677ff5aSBarry Smith 
9583677ff5aSBarry Smith     nontrivial locally owned node - node that is not in the interior (i.e. has more than one
9593677ff5aSBarry Smith            local subdomain
9603677ff5aSBarry Smith   */
96124cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag1);CHKERRQ(ierr);
96224cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag2);CHKERRQ(ierr);
96324cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag3);CHKERRQ(ierr);
96489d82c54SBarry Smith 
96589d82c54SBarry Smith   for (i=0; i<n; i++) {
96689d82c54SBarry Smith     if (lindices[i] > max) max = lindices[i];
96789d82c54SBarry Smith   }
968b2566f29SBarry Smith   ierr   = MPIU_Allreduce(&max,&Ng,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr);
96978058e43SBarry Smith   Ng++;
97089d82c54SBarry Smith   ierr   = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
97189d82c54SBarry Smith   ierr   = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
972bc8ff85bSBarry Smith   scale  = Ng/size + 1;
973a2e34c3dSBarry Smith   ng     = scale; if (rank == size-1) ng = Ng - scale*(size-1); ng = PetscMax(1,ng);
974caba0dd0SBarry Smith   rstart = scale*rank;
97589d82c54SBarry Smith 
97689d82c54SBarry Smith   /* determine ownership ranges of global indices */
977785e854fSJed Brown   ierr = PetscMalloc1(2*size,&nprocs);CHKERRQ(ierr);
97832dcc486SBarry Smith   ierr = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr);
97989d82c54SBarry Smith 
98089d82c54SBarry Smith   /* determine owners of each local node  */
981785e854fSJed Brown   ierr = PetscMalloc1(n,&owner);CHKERRQ(ierr);
98289d82c54SBarry Smith   for (i=0; i<n; i++) {
9833677ff5aSBarry Smith     proc             = lindices[i]/scale; /* processor that globally owns this index */
98427c402fcSBarry Smith     nprocs[2*proc+1] = 1;                 /* processor globally owns at least one of ours */
9853677ff5aSBarry Smith     owner[i]         = proc;
98627c402fcSBarry Smith     nprocs[2*proc]++;                     /* count of how many that processor globally owns of ours */
98789d82c54SBarry Smith   }
98827c402fcSBarry Smith   nsends = 0; for (i=0; i<size; i++) nsends += nprocs[2*i+1];
9897904a332SBarry Smith   ierr = PetscInfo1(mapping,"Number of global owners for my local data %D\n",nsends);CHKERRQ(ierr);
99089d82c54SBarry Smith 
99189d82c54SBarry Smith   /* inform other processors of number of messages and max length*/
99227c402fcSBarry Smith   ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr);
9937904a332SBarry Smith   ierr = PetscInfo1(mapping,"Number of local owners for my global data %D\n",nrecvs);CHKERRQ(ierr);
99489d82c54SBarry Smith 
99589d82c54SBarry Smith   /* post receives for owned rows */
996785e854fSJed Brown   ierr = PetscMalloc1((2*nrecvs+1)*(nmax+1),&recvs);CHKERRQ(ierr);
997854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&recv_waits);CHKERRQ(ierr);
99889d82c54SBarry Smith   for (i=0; i<nrecvs; i++) {
99932dcc486SBarry Smith     ierr = MPI_Irecv(recvs+2*nmax*i,2*nmax,MPIU_INT,MPI_ANY_SOURCE,tag1,comm,recv_waits+i);CHKERRQ(ierr);
100089d82c54SBarry Smith   }
100189d82c54SBarry Smith 
100289d82c54SBarry Smith   /* pack messages containing lists of local nodes to owners */
1003854ce69bSBarry Smith   ierr      = PetscMalloc1(2*n+1,&sends);CHKERRQ(ierr);
1004854ce69bSBarry Smith   ierr      = PetscMalloc1(size+1,&starts);CHKERRQ(ierr);
100589d82c54SBarry Smith   starts[0] = 0;
1006f6e5521dSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
100789d82c54SBarry Smith   for (i=0; i<n; i++) {
100889d82c54SBarry Smith     sends[starts[owner[i]]++] = lindices[i];
100930dcb7c9SBarry Smith     sends[starts[owner[i]]++] = i;
101089d82c54SBarry Smith   }
101189d82c54SBarry Smith   ierr = PetscFree(owner);CHKERRQ(ierr);
101289d82c54SBarry Smith   starts[0] = 0;
1013f6e5521dSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
101489d82c54SBarry Smith 
101589d82c54SBarry Smith   /* send the messages */
1016854ce69bSBarry Smith   ierr = PetscMalloc1(nsends+1,&send_waits);CHKERRQ(ierr);
1017854ce69bSBarry Smith   ierr = PetscMalloc1(nsends+1,&dest);CHKERRQ(ierr);
101889d82c54SBarry Smith   cnt = 0;
101989d82c54SBarry Smith   for (i=0; i<size; i++) {
102027c402fcSBarry Smith     if (nprocs[2*i]) {
102132dcc486SBarry Smith       ierr      = MPI_Isend(sends+starts[i],2*nprocs[2*i],MPIU_INT,i,tag1,comm,send_waits+cnt);CHKERRQ(ierr);
102230dcb7c9SBarry Smith       dest[cnt] = i;
102389d82c54SBarry Smith       cnt++;
102489d82c54SBarry Smith     }
102589d82c54SBarry Smith   }
102689d82c54SBarry Smith   ierr = PetscFree(starts);CHKERRQ(ierr);
102789d82c54SBarry Smith 
102889d82c54SBarry Smith   /* wait on receives */
1029854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&source);CHKERRQ(ierr);
1030854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&len);CHKERRQ(ierr);
103189d82c54SBarry Smith   cnt  = nrecvs;
1032854ce69bSBarry Smith   ierr = PetscMalloc1(ng+1,&nownedsenders);CHKERRQ(ierr);
103332dcc486SBarry Smith   ierr = PetscMemzero(nownedsenders,ng*sizeof(PetscInt));CHKERRQ(ierr);
103489d82c54SBarry Smith   while (cnt) {
103589d82c54SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
103689d82c54SBarry Smith     /* unpack receives into our local space */
103732dcc486SBarry Smith     ierr          = MPI_Get_count(&recv_status,MPIU_INT,&len[imdex]);CHKERRQ(ierr);
103889d82c54SBarry Smith     source[imdex] = recv_status.MPI_SOURCE;
103930dcb7c9SBarry Smith     len[imdex]    = len[imdex]/2;
1040caba0dd0SBarry Smith     /* count how many local owners for each of my global owned indices */
104130dcb7c9SBarry Smith     for (i=0; i<len[imdex]; i++) nownedsenders[recvs[2*imdex*nmax+2*i]-rstart]++;
104289d82c54SBarry Smith     cnt--;
104389d82c54SBarry Smith   }
104489d82c54SBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
104589d82c54SBarry Smith 
104630dcb7c9SBarry Smith   /* count how many globally owned indices are on an edge multiplied by how many processors own them. */
1047bc8ff85bSBarry Smith   nowned  = 0;
1048bc8ff85bSBarry Smith   nownedm = 0;
1049bc8ff85bSBarry Smith   for (i=0; i<ng; i++) {
1050bc8ff85bSBarry Smith     if (nownedsenders[i] > 1) {nownedm += nownedsenders[i]; nowned++;}
1051bc8ff85bSBarry Smith   }
1052bc8ff85bSBarry Smith 
1053bc8ff85bSBarry Smith   /* create single array to contain rank of all local owners of each globally owned index */
1054854ce69bSBarry Smith   ierr      = PetscMalloc1(nownedm+1,&ownedsenders);CHKERRQ(ierr);
1055854ce69bSBarry Smith   ierr      = PetscMalloc1(ng+1,&starts);CHKERRQ(ierr);
1056bc8ff85bSBarry Smith   starts[0] = 0;
1057bc8ff85bSBarry Smith   for (i=1; i<ng; i++) {
1058bc8ff85bSBarry Smith     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
1059bc8ff85bSBarry Smith     else starts[i] = starts[i-1];
1060bc8ff85bSBarry Smith   }
1061bc8ff85bSBarry Smith 
106230dcb7c9SBarry Smith   /* for each nontrival globally owned node list all arriving processors */
1063bc8ff85bSBarry Smith   for (i=0; i<nrecvs; i++) {
1064bc8ff85bSBarry Smith     for (j=0; j<len[i]; j++) {
106530dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
1066f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) ownedsenders[starts[node]++] = source[i];
1067bc8ff85bSBarry Smith     }
1068bc8ff85bSBarry Smith   }
1069bc8ff85bSBarry Smith 
107007b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
107130dcb7c9SBarry Smith     starts[0] = 0;
107230dcb7c9SBarry Smith     for (i=1; i<ng; i++) {
107330dcb7c9SBarry Smith       if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
107430dcb7c9SBarry Smith       else starts[i] = starts[i-1];
107530dcb7c9SBarry Smith     }
107630dcb7c9SBarry Smith     for (i=0; i<ng; i++) {
107730dcb7c9SBarry Smith       if (nownedsenders[i] > 1) {
10787904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"[%d] global node %D local owner processors: ",rank,i+rstart);CHKERRQ(ierr);
107930dcb7c9SBarry Smith         for (j=0; j<nownedsenders[i]; j++) {
10807904a332SBarry Smith           ierr = PetscSynchronizedPrintf(comm,"%D ",ownedsenders[starts[i]+j]);CHKERRQ(ierr);
108130dcb7c9SBarry Smith         }
108230dcb7c9SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
108330dcb7c9SBarry Smith       }
108430dcb7c9SBarry Smith     }
10850ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
108607b52d57SBarry Smith   } /* -----------------------------------  */
108730dcb7c9SBarry Smith 
10883677ff5aSBarry Smith   /* wait on original sends */
10893a96401aSBarry Smith   if (nsends) {
1090785e854fSJed Brown     ierr = PetscMalloc1(nsends,&send_status);CHKERRQ(ierr);
10913a96401aSBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
10923a96401aSBarry Smith     ierr = PetscFree(send_status);CHKERRQ(ierr);
10933a96401aSBarry Smith   }
109489d82c54SBarry Smith   ierr = PetscFree(send_waits);CHKERRQ(ierr);
10953a96401aSBarry Smith   ierr = PetscFree(sends);CHKERRQ(ierr);
10963677ff5aSBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
10973677ff5aSBarry Smith 
10983677ff5aSBarry Smith   /* pack messages to send back to local owners */
109930dcb7c9SBarry Smith   starts[0] = 0;
110030dcb7c9SBarry Smith   for (i=1; i<ng; i++) {
110130dcb7c9SBarry Smith     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
110230dcb7c9SBarry Smith     else starts[i] = starts[i-1];
110330dcb7c9SBarry Smith   }
110430dcb7c9SBarry Smith   nsends2 = nrecvs;
1105854ce69bSBarry Smith   ierr    = PetscMalloc1(nsends2+1,&nprocs);CHKERRQ(ierr); /* length of each message */
110630dcb7c9SBarry Smith   for (i=0; i<nrecvs; i++) {
110730dcb7c9SBarry Smith     nprocs[i] = 1;
110830dcb7c9SBarry Smith     for (j=0; j<len[i]; j++) {
110930dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
1110f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) nprocs[i] += 2 + nownedsenders[node];
111130dcb7c9SBarry Smith     }
111230dcb7c9SBarry Smith   }
1113f6e5521dSKarl Rupp   nt = 0;
1114f6e5521dSKarl Rupp   for (i=0; i<nsends2; i++) nt += nprocs[i];
1115f6e5521dSKarl Rupp 
1116854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,&sends2);CHKERRQ(ierr);
1117854ce69bSBarry Smith   ierr = PetscMalloc1(nsends2+1,&starts2);CHKERRQ(ierr);
1118f6e5521dSKarl Rupp 
1119f6e5521dSKarl Rupp   starts2[0] = 0;
1120f6e5521dSKarl Rupp   for (i=1; i<nsends2; i++) starts2[i] = starts2[i-1] + nprocs[i-1];
112130dcb7c9SBarry Smith   /*
112230dcb7c9SBarry Smith      Each message is 1 + nprocs[i] long, and consists of
112330dcb7c9SBarry Smith        (0) the number of nodes being sent back
112430dcb7c9SBarry Smith        (1) the local node number,
112530dcb7c9SBarry Smith        (2) the number of processors sharing it,
112630dcb7c9SBarry Smith        (3) the processors sharing it
112730dcb7c9SBarry Smith   */
112830dcb7c9SBarry Smith   for (i=0; i<nsends2; i++) {
112930dcb7c9SBarry Smith     cnt = 1;
113030dcb7c9SBarry Smith     sends2[starts2[i]] = 0;
113130dcb7c9SBarry Smith     for (j=0; j<len[i]; j++) {
113230dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
113330dcb7c9SBarry Smith       if (nownedsenders[node] > 1) {
113430dcb7c9SBarry Smith         sends2[starts2[i]]++;
113530dcb7c9SBarry Smith         sends2[starts2[i]+cnt++] = recvs[2*i*nmax+2*j+1];
113630dcb7c9SBarry Smith         sends2[starts2[i]+cnt++] = nownedsenders[node];
113732dcc486SBarry Smith         ierr = PetscMemcpy(&sends2[starts2[i]+cnt],&ownedsenders[starts[node]],nownedsenders[node]*sizeof(PetscInt));CHKERRQ(ierr);
113830dcb7c9SBarry Smith         cnt += nownedsenders[node];
113930dcb7c9SBarry Smith       }
114030dcb7c9SBarry Smith     }
114130dcb7c9SBarry Smith   }
114230dcb7c9SBarry Smith 
114330dcb7c9SBarry Smith   /* receive the message lengths */
114430dcb7c9SBarry Smith   nrecvs2 = nsends;
1145854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&lens2);CHKERRQ(ierr);
1146854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&starts3);CHKERRQ(ierr);
1147854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&recv_waits);CHKERRQ(ierr);
114830dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
1149d44834fbSBarry Smith     ierr = MPI_Irecv(&lens2[i],1,MPIU_INT,dest[i],tag2,comm,recv_waits+i);CHKERRQ(ierr);
115030dcb7c9SBarry Smith   }
1151d44834fbSBarry Smith 
11528a8e0b3aSBarry Smith   /* send the message lengths */
11538a8e0b3aSBarry Smith   for (i=0; i<nsends2; i++) {
11548a8e0b3aSBarry Smith     ierr = MPI_Send(&nprocs[i],1,MPIU_INT,source[i],tag2,comm);CHKERRQ(ierr);
11558a8e0b3aSBarry Smith   }
11568a8e0b3aSBarry Smith 
1157d44834fbSBarry Smith   /* wait on receives of lens */
11580c468ba9SBarry Smith   if (nrecvs2) {
1159785e854fSJed Brown     ierr = PetscMalloc1(nrecvs2,&recv_statuses);CHKERRQ(ierr);
1160d44834fbSBarry Smith     ierr = MPI_Waitall(nrecvs2,recv_waits,recv_statuses);CHKERRQ(ierr);
1161d44834fbSBarry Smith     ierr = PetscFree(recv_statuses);CHKERRQ(ierr);
11620c468ba9SBarry Smith   }
1163a2ea699eSBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
1164d44834fbSBarry Smith 
116530dcb7c9SBarry Smith   starts3[0] = 0;
1166d44834fbSBarry Smith   nt         = 0;
116730dcb7c9SBarry Smith   for (i=0; i<nrecvs2-1; i++) {
116830dcb7c9SBarry Smith     starts3[i+1] = starts3[i] + lens2[i];
1169d44834fbSBarry Smith     nt          += lens2[i];
117030dcb7c9SBarry Smith   }
117176466f69SStefano Zampini   if (nrecvs2) nt += lens2[nrecvs2-1];
1172d44834fbSBarry Smith 
1173854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,&recvs2);CHKERRQ(ierr);
1174854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs2+1,&recv_waits);CHKERRQ(ierr);
117552b72c4aSBarry Smith   for (i=0; i<nrecvs2; i++) {
117632dcc486SBarry Smith     ierr = MPI_Irecv(recvs2+starts3[i],lens2[i],MPIU_INT,dest[i],tag3,comm,recv_waits+i);CHKERRQ(ierr);
117730dcb7c9SBarry Smith   }
117830dcb7c9SBarry Smith 
117930dcb7c9SBarry Smith   /* send the messages */
1180854ce69bSBarry Smith   ierr = PetscMalloc1(nsends2+1,&send_waits);CHKERRQ(ierr);
118130dcb7c9SBarry Smith   for (i=0; i<nsends2; i++) {
118232dcc486SBarry Smith     ierr = MPI_Isend(sends2+starts2[i],nprocs[i],MPIU_INT,source[i],tag3,comm,send_waits+i);CHKERRQ(ierr);
118330dcb7c9SBarry Smith   }
118430dcb7c9SBarry Smith 
118530dcb7c9SBarry Smith   /* wait on receives */
11860c468ba9SBarry Smith   if (nrecvs2) {
1187785e854fSJed Brown     ierr = PetscMalloc1(nrecvs2,&recv_statuses);CHKERRQ(ierr);
118830dcb7c9SBarry Smith     ierr = MPI_Waitall(nrecvs2,recv_waits,recv_statuses);CHKERRQ(ierr);
118930dcb7c9SBarry Smith     ierr = PetscFree(recv_statuses);CHKERRQ(ierr);
11900c468ba9SBarry Smith   }
119130dcb7c9SBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
119230dcb7c9SBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
119330dcb7c9SBarry Smith 
119407b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
119530dcb7c9SBarry Smith     cnt = 0;
119630dcb7c9SBarry Smith     for (i=0; i<nrecvs2; i++) {
119730dcb7c9SBarry Smith       nt = recvs2[cnt++];
119830dcb7c9SBarry Smith       for (j=0; j<nt; j++) {
11997904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"[%d] local node %D number of subdomains %D: ",rank,recvs2[cnt],recvs2[cnt+1]);CHKERRQ(ierr);
120030dcb7c9SBarry Smith         for (k=0; k<recvs2[cnt+1]; k++) {
12017904a332SBarry Smith           ierr = PetscSynchronizedPrintf(comm,"%D ",recvs2[cnt+2+k]);CHKERRQ(ierr);
120230dcb7c9SBarry Smith         }
120330dcb7c9SBarry Smith         cnt += 2 + recvs2[cnt+1];
120430dcb7c9SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
120530dcb7c9SBarry Smith       }
120630dcb7c9SBarry Smith     }
12070ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
120807b52d57SBarry Smith   } /* -----------------------------------  */
120930dcb7c9SBarry Smith 
121030dcb7c9SBarry Smith   /* count number subdomains for each local node */
1211785e854fSJed Brown   ierr = PetscMalloc1(size,&nprocs);CHKERRQ(ierr);
121232dcc486SBarry Smith   ierr = PetscMemzero(nprocs,size*sizeof(PetscInt));CHKERRQ(ierr);
121330dcb7c9SBarry Smith   cnt  = 0;
121430dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
121530dcb7c9SBarry Smith     nt = recvs2[cnt++];
121630dcb7c9SBarry Smith     for (j=0; j<nt; j++) {
1217f6e5521dSKarl Rupp       for (k=0; k<recvs2[cnt+1]; k++) nprocs[recvs2[cnt+2+k]]++;
121830dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt+1];
121930dcb7c9SBarry Smith     }
122030dcb7c9SBarry Smith   }
122130dcb7c9SBarry Smith   nt = 0; for (i=0; i<size; i++) nt += (nprocs[i] > 0);
122230dcb7c9SBarry Smith   *nproc    = nt;
1223854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,procs);CHKERRQ(ierr);
1224854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,numprocs);CHKERRQ(ierr);
1225854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,indices);CHKERRQ(ierr);
12260298fd71SBarry Smith   for (i=0;i<nt+1;i++) (*indices)[i]=NULL;
1227785e854fSJed Brown   ierr = PetscMalloc1(size,&bprocs);CHKERRQ(ierr);
122830dcb7c9SBarry Smith   cnt  = 0;
122930dcb7c9SBarry Smith   for (i=0; i<size; i++) {
123030dcb7c9SBarry Smith     if (nprocs[i] > 0) {
123130dcb7c9SBarry Smith       bprocs[i]        = cnt;
123230dcb7c9SBarry Smith       (*procs)[cnt]    = i;
123330dcb7c9SBarry Smith       (*numprocs)[cnt] = nprocs[i];
1234785e854fSJed Brown       ierr             = PetscMalloc1(nprocs[i],&(*indices)[cnt]);CHKERRQ(ierr);
123530dcb7c9SBarry Smith       cnt++;
123630dcb7c9SBarry Smith     }
123730dcb7c9SBarry Smith   }
123830dcb7c9SBarry Smith 
123930dcb7c9SBarry Smith   /* make the list of subdomains for each nontrivial local node */
124032dcc486SBarry Smith   ierr = PetscMemzero(*numprocs,nt*sizeof(PetscInt));CHKERRQ(ierr);
124130dcb7c9SBarry Smith   cnt  = 0;
124230dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
124330dcb7c9SBarry Smith     nt = recvs2[cnt++];
124430dcb7c9SBarry Smith     for (j=0; j<nt; j++) {
1245f6e5521dSKarl Rupp       for (k=0; k<recvs2[cnt+1]; k++) (*indices)[bprocs[recvs2[cnt+2+k]]][(*numprocs)[bprocs[recvs2[cnt+2+k]]]++] = recvs2[cnt];
124630dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt+1];
124730dcb7c9SBarry Smith     }
124830dcb7c9SBarry Smith   }
124930dcb7c9SBarry Smith   ierr = PetscFree(bprocs);CHKERRQ(ierr);
125007b52d57SBarry Smith   ierr = PetscFree(recvs2);CHKERRQ(ierr);
125130dcb7c9SBarry Smith 
125207b52d57SBarry Smith   /* sort the node indexing by their global numbers */
125307b52d57SBarry Smith   nt = *nproc;
125407b52d57SBarry Smith   for (i=0; i<nt; i++) {
1255854ce69bSBarry Smith     ierr = PetscMalloc1((*numprocs)[i],&tmp);CHKERRQ(ierr);
1256f6e5521dSKarl Rupp     for (j=0; j<(*numprocs)[i]; j++) tmp[j] = lindices[(*indices)[i][j]];
125707b52d57SBarry Smith     ierr = PetscSortIntWithArray((*numprocs)[i],tmp,(*indices)[i]);CHKERRQ(ierr);
125807b52d57SBarry Smith     ierr = PetscFree(tmp);CHKERRQ(ierr);
125907b52d57SBarry Smith   }
126007b52d57SBarry Smith 
126107b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
126230dcb7c9SBarry Smith     nt = *nproc;
126330dcb7c9SBarry Smith     for (i=0; i<nt; i++) {
12647904a332SBarry Smith       ierr = PetscSynchronizedPrintf(comm,"[%d] subdomain %D number of indices %D: ",rank,(*procs)[i],(*numprocs)[i]);CHKERRQ(ierr);
126530dcb7c9SBarry Smith       for (j=0; j<(*numprocs)[i]; j++) {
12667904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"%D ",(*indices)[i][j]);CHKERRQ(ierr);
126730dcb7c9SBarry Smith       }
126830dcb7c9SBarry Smith       ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
126930dcb7c9SBarry Smith     }
12700ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
127107b52d57SBarry Smith   } /* -----------------------------------  */
127230dcb7c9SBarry Smith 
127330dcb7c9SBarry Smith   /* wait on sends */
127430dcb7c9SBarry Smith   if (nsends2) {
1275785e854fSJed Brown     ierr = PetscMalloc1(nsends2,&send_status);CHKERRQ(ierr);
127630dcb7c9SBarry Smith     ierr = MPI_Waitall(nsends2,send_waits,send_status);CHKERRQ(ierr);
127730dcb7c9SBarry Smith     ierr = PetscFree(send_status);CHKERRQ(ierr);
127830dcb7c9SBarry Smith   }
127930dcb7c9SBarry Smith 
128030dcb7c9SBarry Smith   ierr = PetscFree(starts3);CHKERRQ(ierr);
128130dcb7c9SBarry Smith   ierr = PetscFree(dest);CHKERRQ(ierr);
128230dcb7c9SBarry Smith   ierr = PetscFree(send_waits);CHKERRQ(ierr);
12833677ff5aSBarry Smith 
1284bc8ff85bSBarry Smith   ierr = PetscFree(nownedsenders);CHKERRQ(ierr);
1285bc8ff85bSBarry Smith   ierr = PetscFree(ownedsenders);CHKERRQ(ierr);
1286bc8ff85bSBarry Smith   ierr = PetscFree(starts);CHKERRQ(ierr);
128730dcb7c9SBarry Smith   ierr = PetscFree(starts2);CHKERRQ(ierr);
128830dcb7c9SBarry Smith   ierr = PetscFree(lens2);CHKERRQ(ierr);
128989d82c54SBarry Smith 
129089d82c54SBarry Smith   ierr = PetscFree(source);CHKERRQ(ierr);
129197f1f81fSBarry Smith   ierr = PetscFree(len);CHKERRQ(ierr);
129289d82c54SBarry Smith   ierr = PetscFree(recvs);CHKERRQ(ierr);
12933a96401aSBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
129430dcb7c9SBarry Smith   ierr = PetscFree(sends2);CHKERRQ(ierr);
129524cf384cSBarry Smith 
129624cf384cSBarry Smith   /* put the information about myself as the first entry in the list */
129724cf384cSBarry Smith   first_procs    = (*procs)[0];
129824cf384cSBarry Smith   first_numprocs = (*numprocs)[0];
129924cf384cSBarry Smith   first_indices  = (*indices)[0];
130024cf384cSBarry Smith   for (i=0; i<*nproc; i++) {
130124cf384cSBarry Smith     if ((*procs)[i] == rank) {
130224cf384cSBarry Smith       (*procs)[0]    = (*procs)[i];
130324cf384cSBarry Smith       (*numprocs)[0] = (*numprocs)[i];
130424cf384cSBarry Smith       (*indices)[0]  = (*indices)[i];
130524cf384cSBarry Smith       (*procs)[i]    = first_procs;
130624cf384cSBarry Smith       (*numprocs)[i] = first_numprocs;
130724cf384cSBarry Smith       (*indices)[i]  = first_indices;
130824cf384cSBarry Smith       break;
130924cf384cSBarry Smith     }
131024cf384cSBarry Smith   }
1311268a049cSStefano Zampini 
1312268a049cSStefano Zampini   /* save info for reuse */
1313268a049cSStefano Zampini   mapping->info_nproc = *nproc;
1314268a049cSStefano Zampini   mapping->info_procs = *procs;
1315268a049cSStefano Zampini   mapping->info_numprocs = *numprocs;
1316268a049cSStefano Zampini   mapping->info_indices = *indices;
1317268a049cSStefano Zampini   mapping->info_cached = PETSC_TRUE;
131889d82c54SBarry Smith   PetscFunctionReturn(0);
131989d82c54SBarry Smith }
132089d82c54SBarry Smith 
13216a818285SBarry Smith /*@C
13226a818285SBarry Smith     ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by ISLocalToGlobalMappingGetBlockInfo()
13236a818285SBarry Smith 
13246a818285SBarry Smith     Collective on ISLocalToGlobalMapping
13256a818285SBarry Smith 
13266a818285SBarry Smith     Input Parameters:
13276a818285SBarry Smith .   mapping - the mapping from local to global indexing
13286a818285SBarry Smith 
13296a818285SBarry Smith     Output Parameter:
13306a818285SBarry Smith +   nproc - number of processors that are connected to this one
13316a818285SBarry Smith .   proc - neighboring processors
13326a818285SBarry Smith .   numproc - number of indices for each processor
13336a818285SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
13346a818285SBarry Smith 
13356a818285SBarry Smith     Level: advanced
13366a818285SBarry Smith 
13376a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
13386a818285SBarry Smith           ISLocalToGlobalMappingGetInfo()
13396a818285SBarry Smith @*/
13406a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
13416a818285SBarry Smith {
13426a818285SBarry Smith   PetscErrorCode ierr;
13436a818285SBarry Smith 
13446a818285SBarry Smith   PetscFunctionBegin;
1345cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1346268a049cSStefano Zampini   if (mapping->info_free) {
13476a818285SBarry Smith     ierr = PetscFree(*numprocs);CHKERRQ(ierr);
13486a818285SBarry Smith     if (*indices) {
1349268a049cSStefano Zampini       PetscInt i;
1350268a049cSStefano Zampini 
13516a818285SBarry Smith       ierr = PetscFree((*indices)[0]);CHKERRQ(ierr);
13526a818285SBarry Smith       for (i=1; i<*nproc; i++) {
13536a818285SBarry Smith         ierr = PetscFree((*indices)[i]);CHKERRQ(ierr);
13546a818285SBarry Smith       }
13556a818285SBarry Smith       ierr = PetscFree(*indices);CHKERRQ(ierr);
13566a818285SBarry Smith     }
1357268a049cSStefano Zampini   }
1358268a049cSStefano Zampini   *nproc    = 0;
1359268a049cSStefano Zampini   *procs    = NULL;
1360268a049cSStefano Zampini   *numprocs = NULL;
1361268a049cSStefano Zampini   *indices  = NULL;
13626a818285SBarry Smith   PetscFunctionReturn(0);
13636a818285SBarry Smith }
13646a818285SBarry Smith 
13656a818285SBarry Smith /*@C
13666a818285SBarry Smith     ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each processor and
13676a818285SBarry Smith      each index shared by more than one processor
13686a818285SBarry Smith 
13696a818285SBarry Smith     Collective on ISLocalToGlobalMapping
13706a818285SBarry Smith 
13716a818285SBarry Smith     Input Parameters:
13726a818285SBarry Smith .   mapping - the mapping from local to global indexing
13736a818285SBarry Smith 
13746a818285SBarry Smith     Output Parameter:
13756a818285SBarry Smith +   nproc - number of processors that are connected to this one
13766a818285SBarry Smith .   proc - neighboring processors
13776a818285SBarry Smith .   numproc - number of indices for each subdomain (processor)
13786a818285SBarry Smith -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
13796a818285SBarry Smith 
13806a818285SBarry Smith     Level: advanced
13816a818285SBarry Smith 
13826a818285SBarry Smith     Concepts: mapping^local to global
13836a818285SBarry Smith 
13846a818285SBarry Smith     Fortran Usage:
13856a818285SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
13866a818285SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
13876a818285SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
13886a818285SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
13896a818285SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
13906a818285SBarry Smith 
13916a818285SBarry Smith 
13926a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
13936a818285SBarry Smith           ISLocalToGlobalMappingRestoreInfo()
13946a818285SBarry Smith @*/
13956a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
13966a818285SBarry Smith {
13976a818285SBarry Smith   PetscErrorCode ierr;
1398268a049cSStefano Zampini   PetscInt       **bindices = NULL,*bnumprocs = NULL,bs = mapping->bs,i,j,k;
13996a818285SBarry Smith 
14006a818285SBarry Smith   PetscFunctionBegin;
14016a818285SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1402268a049cSStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockInfo(mapping,nproc,procs,&bnumprocs,&bindices);CHKERRQ(ierr);
1403268a049cSStefano Zampini   if (bs > 1) { /* we need to expand the cached info */
1404732f65e3SBarry Smith     ierr = PetscCalloc1(*nproc,&*indices);CHKERRQ(ierr);
1405268a049cSStefano Zampini     ierr = PetscCalloc1(*nproc,&*numprocs);CHKERRQ(ierr);
14066a818285SBarry Smith     for (i=0; i<*nproc; i++) {
1407268a049cSStefano Zampini       ierr = PetscMalloc1(bs*bnumprocs[i],&(*indices)[i]);CHKERRQ(ierr);
1408268a049cSStefano Zampini       for (j=0; j<bnumprocs[i]; j++) {
14096a818285SBarry Smith         for (k=0; k<bs; k++) {
14106a818285SBarry Smith           (*indices)[i][j*bs+k] = bs*bindices[i][j] + k;
14116a818285SBarry Smith         }
14126a818285SBarry Smith       }
1413268a049cSStefano Zampini       (*numprocs)[i] = bnumprocs[i]*bs;
14146a818285SBarry Smith     }
1415268a049cSStefano Zampini     mapping->info_free = PETSC_TRUE;
1416268a049cSStefano Zampini   } else {
1417268a049cSStefano Zampini     *numprocs = bnumprocs;
1418268a049cSStefano Zampini     *indices  = bindices;
14196a818285SBarry Smith   }
14206a818285SBarry Smith   PetscFunctionReturn(0);
14216a818285SBarry Smith }
14226a818285SBarry Smith 
142307b52d57SBarry Smith /*@C
142407b52d57SBarry Smith     ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by ISLocalToGlobalMappingGetInfo()
142589d82c54SBarry Smith 
142607b52d57SBarry Smith     Collective on ISLocalToGlobalMapping
142707b52d57SBarry Smith 
142807b52d57SBarry Smith     Input Parameters:
142907b52d57SBarry Smith .   mapping - the mapping from local to global indexing
143007b52d57SBarry Smith 
143107b52d57SBarry Smith     Output Parameter:
143207b52d57SBarry Smith +   nproc - number of processors that are connected to this one
143307b52d57SBarry Smith .   proc - neighboring processors
143407b52d57SBarry Smith .   numproc - number of indices for each processor
143507b52d57SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
143607b52d57SBarry Smith 
143707b52d57SBarry Smith     Level: advanced
143807b52d57SBarry Smith 
143907b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
144007b52d57SBarry Smith           ISLocalToGlobalMappingGetInfo()
144107b52d57SBarry Smith @*/
14427087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
144307b52d57SBarry Smith {
14446849ba73SBarry Smith   PetscErrorCode ierr;
144507b52d57SBarry Smith 
144607b52d57SBarry Smith   PetscFunctionBegin;
14476a818285SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockInfo(mapping,nproc,procs,numprocs,indices);CHKERRQ(ierr);
144807b52d57SBarry Smith   PetscFunctionReturn(0);
144907b52d57SBarry Smith }
145086994e45SJed Brown 
145186994e45SJed Brown /*@C
1452107e9a97SBarry Smith    ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped
145386994e45SJed Brown 
145486994e45SJed Brown    Not Collective
145586994e45SJed Brown 
145686994e45SJed Brown    Input Arguments:
145786994e45SJed Brown . ltog - local to global mapping
145886994e45SJed Brown 
145986994e45SJed Brown    Output Arguments:
1460565245c5SBarry Smith . array - array of indices, the length of this array may be obtained with ISLocalToGlobalMappingGetSize()
146186994e45SJed Brown 
146286994e45SJed Brown    Level: advanced
146386994e45SJed Brown 
1464107e9a97SBarry Smith    Notes: ISLocalToGlobalMappingGetSize() returns the length the this array
1465107e9a97SBarry Smith 
1466107e9a97SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreIndices(), ISLocalToGlobalMappingGetBlockIndices(), ISLocalToGlobalMappingRestoreBlockIndices()
146786994e45SJed Brown @*/
14687087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
146986994e45SJed Brown {
147086994e45SJed Brown   PetscFunctionBegin;
147186994e45SJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
147286994e45SJed Brown   PetscValidPointer(array,2);
147345b6f7e9SBarry Smith   if (ltog->bs == 1) {
147486994e45SJed Brown     *array = ltog->indices;
147545b6f7e9SBarry Smith   } else {
147645b6f7e9SBarry Smith     PetscInt       *jj,k,i,j,n = ltog->n, bs = ltog->bs;
147745b6f7e9SBarry Smith     const PetscInt *ii;
147845b6f7e9SBarry Smith     PetscErrorCode ierr;
147945b6f7e9SBarry Smith 
148045b6f7e9SBarry Smith     ierr = PetscMalloc1(bs*n,&jj);CHKERRQ(ierr);
148145b6f7e9SBarry Smith     *array = jj;
148245b6f7e9SBarry Smith     k    = 0;
148345b6f7e9SBarry Smith     ii   = ltog->indices;
148445b6f7e9SBarry Smith     for (i=0; i<n; i++)
148545b6f7e9SBarry Smith       for (j=0; j<bs; j++)
148645b6f7e9SBarry Smith         jj[k++] = bs*ii[i] + j;
148745b6f7e9SBarry Smith   }
148886994e45SJed Brown   PetscFunctionReturn(0);
148986994e45SJed Brown }
149086994e45SJed Brown 
149186994e45SJed Brown /*@C
1492193a2b41SJulian Andrej    ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with ISLocalToGlobalMappingGetIndices()
149386994e45SJed Brown 
149486994e45SJed Brown    Not Collective
149586994e45SJed Brown 
149686994e45SJed Brown    Input Arguments:
149786994e45SJed Brown + ltog - local to global mapping
149886994e45SJed Brown - array - array of indices
149986994e45SJed Brown 
150086994e45SJed Brown    Level: advanced
150186994e45SJed Brown 
150286994e45SJed Brown .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices()
150386994e45SJed Brown @*/
15047087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
150586994e45SJed Brown {
150686994e45SJed Brown   PetscFunctionBegin;
150786994e45SJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
150886994e45SJed Brown   PetscValidPointer(array,2);
150945b6f7e9SBarry Smith   if (ltog->bs == 1 && *array != ltog->indices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
151045b6f7e9SBarry Smith 
151145b6f7e9SBarry Smith   if (ltog->bs > 1) {
151245b6f7e9SBarry Smith     PetscErrorCode ierr;
151345b6f7e9SBarry Smith     ierr = PetscFree(*(void**)array);CHKERRQ(ierr);
151445b6f7e9SBarry Smith   }
151545b6f7e9SBarry Smith   PetscFunctionReturn(0);
151645b6f7e9SBarry Smith }
151745b6f7e9SBarry Smith 
151845b6f7e9SBarry Smith /*@C
151945b6f7e9SBarry Smith    ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block
152045b6f7e9SBarry Smith 
152145b6f7e9SBarry Smith    Not Collective
152245b6f7e9SBarry Smith 
152345b6f7e9SBarry Smith    Input Arguments:
152445b6f7e9SBarry Smith . ltog - local to global mapping
152545b6f7e9SBarry Smith 
152645b6f7e9SBarry Smith    Output Arguments:
152745b6f7e9SBarry Smith . array - array of indices
152845b6f7e9SBarry Smith 
152945b6f7e9SBarry Smith    Level: advanced
153045b6f7e9SBarry Smith 
153145b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreBlockIndices()
153245b6f7e9SBarry Smith @*/
153345b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
153445b6f7e9SBarry Smith {
153545b6f7e9SBarry Smith   PetscFunctionBegin;
153645b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
153745b6f7e9SBarry Smith   PetscValidPointer(array,2);
153845b6f7e9SBarry Smith   *array = ltog->indices;
153945b6f7e9SBarry Smith   PetscFunctionReturn(0);
154045b6f7e9SBarry Smith }
154145b6f7e9SBarry Smith 
154245b6f7e9SBarry Smith /*@C
154345b6f7e9SBarry Smith    ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with ISLocalToGlobalMappingGetBlockIndices()
154445b6f7e9SBarry Smith 
154545b6f7e9SBarry Smith    Not Collective
154645b6f7e9SBarry Smith 
154745b6f7e9SBarry Smith    Input Arguments:
154845b6f7e9SBarry Smith + ltog - local to global mapping
154945b6f7e9SBarry Smith - array - array of indices
155045b6f7e9SBarry Smith 
155145b6f7e9SBarry Smith    Level: advanced
155245b6f7e9SBarry Smith 
155345b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices()
155445b6f7e9SBarry Smith @*/
155545b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
155645b6f7e9SBarry Smith {
155745b6f7e9SBarry Smith   PetscFunctionBegin;
155845b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
155945b6f7e9SBarry Smith   PetscValidPointer(array,2);
156086994e45SJed Brown   if (*array != ltog->indices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
15610298fd71SBarry Smith   *array = NULL;
156286994e45SJed Brown   PetscFunctionReturn(0);
156386994e45SJed Brown }
1564f7efa3c7SJed Brown 
1565f7efa3c7SJed Brown /*@C
1566f7efa3c7SJed Brown    ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings
1567f7efa3c7SJed Brown 
1568f7efa3c7SJed Brown    Not Collective
1569f7efa3c7SJed Brown 
1570f7efa3c7SJed Brown    Input Arguments:
1571f7efa3c7SJed Brown + comm - communicator for the new mapping, must contain the communicator of every mapping to concatenate
1572f7efa3c7SJed Brown . n - number of mappings to concatenate
1573f7efa3c7SJed Brown - ltogs - local to global mappings
1574f7efa3c7SJed Brown 
1575f7efa3c7SJed Brown    Output Arguments:
1576f7efa3c7SJed Brown . ltogcat - new mapping
1577f7efa3c7SJed Brown 
15789d90f715SBarry Smith    Note: this currently always returns a mapping with block size of 1
15799d90f715SBarry Smith 
15809d90f715SBarry Smith    Developer Note: If all the input mapping have the same block size we could easily handle that as a special case
15819d90f715SBarry Smith 
1582f7efa3c7SJed Brown    Level: advanced
1583f7efa3c7SJed Brown 
1584f7efa3c7SJed Brown .seealso: ISLocalToGlobalMappingCreate()
1585f7efa3c7SJed Brown @*/
1586f7efa3c7SJed Brown PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm,PetscInt n,const ISLocalToGlobalMapping ltogs[],ISLocalToGlobalMapping *ltogcat)
1587f7efa3c7SJed Brown {
1588f7efa3c7SJed Brown   PetscInt       i,cnt,m,*idx;
1589f7efa3c7SJed Brown   PetscErrorCode ierr;
1590f7efa3c7SJed Brown 
1591f7efa3c7SJed Brown   PetscFunctionBegin;
1592f7efa3c7SJed Brown   if (n < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have a non-negative number of mappings, given %D",n);
1593f7efa3c7SJed Brown   if (n > 0) PetscValidPointer(ltogs,3);
1594f7efa3c7SJed Brown   for (i=0; i<n; i++) PetscValidHeaderSpecific(ltogs[i],IS_LTOGM_CLASSID,3);
1595f7efa3c7SJed Brown   PetscValidPointer(ltogcat,4);
1596f7efa3c7SJed Brown   for (cnt=0,i=0; i<n; i++) {
1597f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
1598f7efa3c7SJed Brown     cnt += m;
1599f7efa3c7SJed Brown   }
1600785e854fSJed Brown   ierr = PetscMalloc1(cnt,&idx);CHKERRQ(ierr);
1601f7efa3c7SJed Brown   for (cnt=0,i=0; i<n; i++) {
1602f7efa3c7SJed Brown     const PetscInt *subidx;
1603f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
1604f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetIndices(ltogs[i],&subidx);CHKERRQ(ierr);
1605f7efa3c7SJed Brown     ierr = PetscMemcpy(&idx[cnt],subidx,m*sizeof(PetscInt));CHKERRQ(ierr);
1606f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingRestoreIndices(ltogs[i],&subidx);CHKERRQ(ierr);
1607f7efa3c7SJed Brown     cnt += m;
1608f7efa3c7SJed Brown   }
1609f0413b6fSBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,1,cnt,idx,PETSC_OWN_POINTER,ltogcat);CHKERRQ(ierr);
1610f7efa3c7SJed Brown   PetscFunctionReturn(0);
1611f7efa3c7SJed Brown }
161204a59952SBarry Smith 
1613413f72f0SBarry Smith /*MC
1614413f72f0SBarry Smith       ISLOCALTOGLOBALMAPPINGBASIC - basic implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is
1615413f72f0SBarry Smith                                     used this is good for only small and moderate size problems.
1616413f72f0SBarry Smith 
1617413f72f0SBarry Smith    Options Database Keys:
1618413f72f0SBarry Smith +   -islocaltoglobalmapping_type basic - select this method
1619413f72f0SBarry Smith 
1620413f72f0SBarry Smith    Level: beginner
1621413f72f0SBarry Smith 
1622413f72f0SBarry Smith .seealso:  ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType(), ISLOCALTOGLOBALMAPPINGHASH
1623413f72f0SBarry Smith M*/
1624413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Basic(ISLocalToGlobalMapping ltog)
1625413f72f0SBarry Smith {
1626413f72f0SBarry Smith   PetscFunctionBegin;
1627413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Basic;
1628413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Basic;
1629413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Basic;
1630413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Basic;
1631413f72f0SBarry Smith   PetscFunctionReturn(0);
1632413f72f0SBarry Smith }
1633413f72f0SBarry Smith 
1634413f72f0SBarry Smith /*MC
1635413f72f0SBarry Smith       ISLOCALTOGLOBALMAPPINGHASH - hash implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is
1636*ed56e8eaSBarry Smith                                     used this is good for large memory problems.
1637413f72f0SBarry Smith 
1638413f72f0SBarry Smith    Options Database Keys:
1639413f72f0SBarry Smith +   -islocaltoglobalmapping_type hash - select this method
1640413f72f0SBarry Smith 
1641*ed56e8eaSBarry Smith 
1642*ed56e8eaSBarry Smith    Notes: This is selected automatically for large problems if the user does not set the type.
1643*ed56e8eaSBarry Smith 
1644413f72f0SBarry Smith    Level: beginner
1645413f72f0SBarry Smith 
1646413f72f0SBarry Smith .seealso:  ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType(), ISLOCALTOGLOBALMAPPINGHASH
1647413f72f0SBarry Smith M*/
1648413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Hash(ISLocalToGlobalMapping ltog)
1649413f72f0SBarry Smith {
1650413f72f0SBarry Smith   PetscFunctionBegin;
1651413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Hash;
1652413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Hash;
1653413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Hash;
1654413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Hash;
1655413f72f0SBarry Smith   PetscFunctionReturn(0);
1656413f72f0SBarry Smith }
1657413f72f0SBarry Smith 
1658413f72f0SBarry Smith 
1659413f72f0SBarry Smith /*@C
1660413f72f0SBarry Smith     ISLocalToGlobalMappingRegister -  Adds a method for applying a global to local mapping with an ISLocalToGlobalMapping
1661413f72f0SBarry Smith 
1662413f72f0SBarry Smith    Not Collective
1663413f72f0SBarry Smith 
1664413f72f0SBarry Smith    Input Parameters:
1665413f72f0SBarry Smith +  sname - name of a new method
1666413f72f0SBarry Smith -  routine_create - routine to create method context
1667413f72f0SBarry Smith 
1668413f72f0SBarry Smith    Notes:
1669*ed56e8eaSBarry Smith    ISLocalToGlobalMappingRegister() may be called multiple times to add several user-defined mappings.
1670413f72f0SBarry Smith 
1671413f72f0SBarry Smith    Sample usage:
1672413f72f0SBarry Smith .vb
1673413f72f0SBarry Smith    ISLocalToGlobalMappingRegister("my_mapper",MyCreate);
1674413f72f0SBarry Smith .ve
1675413f72f0SBarry Smith 
1676*ed56e8eaSBarry Smith    Then, your mapping can be chosen with the procedural interface via
1677413f72f0SBarry Smith $     ISLocalToGlobalMappingSetType(ltog,"my_mapper")
1678413f72f0SBarry Smith    or at runtime via the option
1679*ed56e8eaSBarry Smith $     -islocaltoglobalmapping_type my_mapper
1680413f72f0SBarry Smith 
1681413f72f0SBarry Smith    Level: advanced
1682413f72f0SBarry Smith 
1683413f72f0SBarry Smith .keywords: ISLocalToGlobalMappingType, register
1684413f72f0SBarry Smith 
1685413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingRegisterAll(), ISLocalToGlobalMappingRegisterDestroy(), ISLOCALTOGLOBALMAPPINGBASIC, ISLOCALTOGLOBALMAPPINGHASH
1686413f72f0SBarry Smith 
1687413f72f0SBarry Smith @*/
1688413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRegister(const char sname[],PetscErrorCode (*function)(ISLocalToGlobalMapping))
1689413f72f0SBarry Smith {
1690413f72f0SBarry Smith   PetscErrorCode ierr;
1691413f72f0SBarry Smith 
1692413f72f0SBarry Smith   PetscFunctionBegin;
1693413f72f0SBarry Smith   ierr = PetscFunctionListAdd(&ISLocalToGlobalMappingList,sname,function);CHKERRQ(ierr);
1694413f72f0SBarry Smith   PetscFunctionReturn(0);
1695413f72f0SBarry Smith }
1696413f72f0SBarry Smith 
1697413f72f0SBarry Smith /*@C
1698*ed56e8eaSBarry Smith    ISLocalToGlobalMappingSetType - Builds ISLocalToGlobalMapping for a particular global to local mapping approach.
1699413f72f0SBarry Smith 
1700413f72f0SBarry Smith    Logically Collective on ISLocalToGlobalMapping
1701413f72f0SBarry Smith 
1702413f72f0SBarry Smith    Input Parameters:
1703413f72f0SBarry Smith +  ltog - the ISLocalToGlobalMapping object
1704413f72f0SBarry Smith -  type - a known method
1705413f72f0SBarry Smith 
1706413f72f0SBarry Smith    Options Database Key:
1707*ed56e8eaSBarry Smith .  -islocaltoglobalmapping_type  <method> - Sets the method; use -help for a list
1708413f72f0SBarry Smith     of available methods (for instance, basic or hash)
1709413f72f0SBarry Smith 
1710413f72f0SBarry Smith    Notes:
1711413f72f0SBarry Smith    See "petsc/include/petscis.h" for available methods
1712413f72f0SBarry Smith 
1713413f72f0SBarry Smith   Normally, it is best to use the ISLocalToGlobalMappingSetFromOptions() command and
1714413f72f0SBarry Smith   then set the ISLocalToGlobalMapping type from the options database rather than by using
1715413f72f0SBarry Smith   this routine.
1716413f72f0SBarry Smith 
1717413f72f0SBarry Smith   Level: intermediate
1718413f72f0SBarry Smith 
1719413f72f0SBarry Smith   Developer Note: ISLocalToGlobalMappingRegister() is used to add new types to ISLocalToGlobalMappingList from which they
1720413f72f0SBarry Smith   are accessed by ISLocalToGlobalMappingSetType().
1721413f72f0SBarry Smith 
1722413f72f0SBarry Smith .keywords: ISLocalToGlobalMapping, set, method
1723413f72f0SBarry Smith 
1724413f72f0SBarry Smith .seealso: PCSetType(), ISLocalToGlobalMappingType, ISLocalToGlobalMappingRegister(), ISLocalToGlobalMappingCreate()
1725413f72f0SBarry Smith 
1726413f72f0SBarry Smith @*/
1727413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType type)
1728413f72f0SBarry Smith {
1729413f72f0SBarry Smith   PetscErrorCode ierr,(*r)(ISLocalToGlobalMapping);
1730413f72f0SBarry Smith   PetscBool      match;
1731413f72f0SBarry Smith 
1732413f72f0SBarry Smith   PetscFunctionBegin;
1733413f72f0SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
1734413f72f0SBarry Smith   PetscValidCharPointer(type,2);
1735413f72f0SBarry Smith 
1736413f72f0SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)ltog,type,&match);CHKERRQ(ierr);
1737413f72f0SBarry Smith   if (match) PetscFunctionReturn(0);
1738413f72f0SBarry Smith 
1739413f72f0SBarry Smith   ierr =  PetscFunctionListFind(ISLocalToGlobalMappingList,type,&r);CHKERRQ(ierr);
1740413f72f0SBarry Smith   if (!r) SETERRQ1(PetscObjectComm((PetscObject)ltog),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested ISLocalToGlobalMapping type %s",type);
1741413f72f0SBarry Smith   /* Destroy the previous private LTOG context */
1742413f72f0SBarry Smith   if (ltog->ops->destroy) {
1743413f72f0SBarry Smith     ierr              = (*ltog->ops->destroy)(ltog);CHKERRQ(ierr);
1744413f72f0SBarry Smith     ltog->ops->destroy = NULL;
1745413f72f0SBarry Smith   }
1746413f72f0SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)ltog,type);CHKERRQ(ierr);
1747413f72f0SBarry Smith   ierr = (*r)(ltog);CHKERRQ(ierr);
1748413f72f0SBarry Smith   PetscFunctionReturn(0);
1749413f72f0SBarry Smith }
1750413f72f0SBarry Smith 
1751413f72f0SBarry Smith PetscBool ISLocalToGlobalMappingRegisterAllCalled = PETSC_FALSE;
1752413f72f0SBarry Smith 
1753413f72f0SBarry Smith /*@C
1754413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAll - Registers all of the local to global mapping components in the IS package.
1755413f72f0SBarry Smith 
1756413f72f0SBarry Smith   Not Collective
1757413f72f0SBarry Smith 
1758413f72f0SBarry Smith   Level: advanced
1759413f72f0SBarry Smith 
1760413f72f0SBarry Smith .keywords: IS, register, all
1761413f72f0SBarry Smith .seealso:  ISRegister(),  ISLocalToGlobalRegister()
1762413f72f0SBarry Smith @*/
1763413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRegisterAll(void)
1764413f72f0SBarry Smith {
1765413f72f0SBarry Smith   PetscErrorCode ierr;
1766413f72f0SBarry Smith 
1767413f72f0SBarry Smith   PetscFunctionBegin;
1768413f72f0SBarry Smith   if (ISLocalToGlobalMappingRegisterAllCalled) PetscFunctionReturn(0);
1769413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAllCalled = PETSC_TRUE;
1770413f72f0SBarry Smith 
1771413f72f0SBarry Smith   ierr = ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGBASIC, ISLocalToGlobalMappingCreate_Basic);CHKERRQ(ierr);
1772413f72f0SBarry Smith   ierr = ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGHASH, ISLocalToGlobalMappingCreate_Hash);CHKERRQ(ierr);
1773413f72f0SBarry Smith   PetscFunctionReturn(0);
1774413f72f0SBarry Smith }
177504a59952SBarry Smith 
1776