xref: /petsc/src/vec/is/utils/isltog.c (revision 2785b3218cc31e82d725b9b492b549e1063b7c2b)
12362add9SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/isimpl.h>    /*I "petscis.h"  I*/
3bbf7bc21SLisandro Dalcin #include <petsc/private/hash.h>
40c312b8eSJed Brown #include <petscsf.h>
5665c2dedSJed Brown #include <petscviewer.h>
62362add9SBarry Smith 
77087cfbeSBarry Smith PetscClassId IS_LTOGM_CLASSID;
8268a049cSStefano Zampini static PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping,PetscInt*,PetscInt**,PetscInt**,PetscInt***);
98e58c17dSMatthew Knepley 
10413f72f0SBarry Smith typedef struct {
11413f72f0SBarry Smith   PetscInt    *globals;
12413f72f0SBarry Smith } ISLocalToGlobalMapping_Basic;
13413f72f0SBarry Smith 
14413f72f0SBarry Smith typedef struct {
15413f72f0SBarry Smith   PetscHashI globalht;
16413f72f0SBarry Smith } ISLocalToGlobalMapping_Hash;
17413f72f0SBarry Smith 
18413f72f0SBarry Smith 
19413f72f0SBarry Smith /* -----------------------------------------------------------------------------------------*/
20413f72f0SBarry Smith 
21413f72f0SBarry Smith /*
22413f72f0SBarry Smith     Creates the global mapping information in the ISLocalToGlobalMapping structure
23413f72f0SBarry Smith 
24413f72f0SBarry Smith     If the user has not selected how to handle the global to local mapping then use HASH for "large" problems
25413f72f0SBarry Smith */
26413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp(ISLocalToGlobalMapping mapping)
27413f72f0SBarry Smith {
28413f72f0SBarry Smith   PetscInt       i,*idx = mapping->indices,n = mapping->n,end,start;
29413f72f0SBarry Smith   PetscErrorCode ierr;
30413f72f0SBarry Smith 
31413f72f0SBarry Smith   PetscFunctionBegin;
32413f72f0SBarry Smith   if (mapping->data) PetscFunctionReturn(0);
33413f72f0SBarry Smith   end   = 0;
34413f72f0SBarry Smith   start = PETSC_MAX_INT;
35413f72f0SBarry Smith 
36413f72f0SBarry Smith   for (i=0; i<n; i++) {
37413f72f0SBarry Smith     if (idx[i] < 0) continue;
38413f72f0SBarry Smith     if (idx[i] < start) start = idx[i];
39413f72f0SBarry Smith     if (idx[i] > end)   end   = idx[i];
40413f72f0SBarry Smith   }
41413f72f0SBarry Smith   if (start > end) {start = 0; end = -1;}
42413f72f0SBarry Smith   mapping->globalstart = start;
43413f72f0SBarry Smith   mapping->globalend   = end;
44413f72f0SBarry Smith   if (!((PetscObject)mapping)->type_name) {
45413f72f0SBarry Smith     if ((end - start) > PetscMax(4*n,1000000)) {
467f79407eSBarry Smith       ierr = ISLocalToGlobalMappingSetType(mapping,ISLOCALTOGLOBALMAPPINGHASH);CHKERRQ(ierr);
47413f72f0SBarry Smith     } else {
487f79407eSBarry Smith       ierr = ISLocalToGlobalMappingSetType(mapping,ISLOCALTOGLOBALMAPPINGBASIC);CHKERRQ(ierr);
49413f72f0SBarry Smith     }
50413f72f0SBarry Smith   }
51413f72f0SBarry Smith   ierr = (*mapping->ops->globaltolocalmappingsetup)(mapping);CHKERRQ(ierr);
52413f72f0SBarry Smith   PetscFunctionReturn(0);
53413f72f0SBarry Smith }
54413f72f0SBarry Smith 
55413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp_Basic(ISLocalToGlobalMapping mapping)
56413f72f0SBarry Smith {
57413f72f0SBarry Smith   PetscErrorCode              ierr;
58413f72f0SBarry Smith   PetscInt                    i,*idx = mapping->indices,n = mapping->n,end,start,*globals;
59413f72f0SBarry Smith   ISLocalToGlobalMapping_Basic *map;
60413f72f0SBarry Smith 
61413f72f0SBarry Smith   PetscFunctionBegin;
62413f72f0SBarry Smith   start            = mapping->globalstart;
63413f72f0SBarry Smith   end              = mapping->globalend;
64413f72f0SBarry Smith   ierr             = PetscNew(&map);CHKERRQ(ierr);
65413f72f0SBarry Smith   ierr             = PetscMalloc1(end-start+2,&globals);CHKERRQ(ierr);
66413f72f0SBarry Smith   map->globals     = globals;
67413f72f0SBarry Smith   for (i=0; i<end-start+1; i++) globals[i] = -1;
68413f72f0SBarry Smith   for (i=0; i<n; i++) {
69413f72f0SBarry Smith     if (idx[i] < 0) continue;
70413f72f0SBarry Smith     globals[idx[i] - start] = i;
71413f72f0SBarry Smith   }
72413f72f0SBarry Smith   mapping->data = (void*)map;
73413f72f0SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)mapping,(end-start+1)*sizeof(PetscInt));CHKERRQ(ierr);
74413f72f0SBarry Smith   PetscFunctionReturn(0);
75413f72f0SBarry Smith }
76413f72f0SBarry Smith 
77413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp_Hash(ISLocalToGlobalMapping mapping)
78413f72f0SBarry Smith {
79413f72f0SBarry Smith   PetscErrorCode              ierr;
80413f72f0SBarry Smith   PetscInt                    i,*idx = mapping->indices,n = mapping->n;
81413f72f0SBarry Smith   ISLocalToGlobalMapping_Hash *map;
82413f72f0SBarry Smith 
83413f72f0SBarry Smith   PetscFunctionBegin;
84413f72f0SBarry Smith   ierr = PetscNew(&map);CHKERRQ(ierr);
85413f72f0SBarry Smith   PetscHashICreate(map->globalht);
86413f72f0SBarry Smith   for (i=0; i<n; i++ ) {
87413f72f0SBarry Smith     if (idx[i] < 0) continue;
88413f72f0SBarry Smith     PetscHashIAdd(map->globalht, idx[i], i);
89413f72f0SBarry Smith   }
90413f72f0SBarry Smith   mapping->data = (void*)map;
91413f72f0SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)mapping,2*n*sizeof(PetscInt));CHKERRQ(ierr);
92413f72f0SBarry Smith   PetscFunctionReturn(0);
93413f72f0SBarry Smith }
94413f72f0SBarry Smith 
95413f72f0SBarry Smith static PetscErrorCode ISLocalToGlobalMappingDestroy_Basic(ISLocalToGlobalMapping mapping)
96413f72f0SBarry Smith {
97413f72f0SBarry Smith   PetscErrorCode              ierr;
98413f72f0SBarry Smith   ISLocalToGlobalMapping_Basic *map  = (ISLocalToGlobalMapping_Basic *)mapping->data;
99413f72f0SBarry Smith 
100413f72f0SBarry Smith   PetscFunctionBegin;
101413f72f0SBarry Smith   if (!map) PetscFunctionReturn(0);
102413f72f0SBarry Smith   ierr = PetscFree(map->globals);CHKERRQ(ierr);
103413f72f0SBarry Smith   ierr = PetscFree(mapping->data);CHKERRQ(ierr);
104413f72f0SBarry Smith   PetscFunctionReturn(0);
105413f72f0SBarry Smith }
106413f72f0SBarry Smith 
107413f72f0SBarry Smith static PetscErrorCode ISLocalToGlobalMappingDestroy_Hash(ISLocalToGlobalMapping mapping)
108413f72f0SBarry Smith {
109413f72f0SBarry Smith   PetscErrorCode              ierr;
110413f72f0SBarry Smith   ISLocalToGlobalMapping_Hash *map  = (ISLocalToGlobalMapping_Hash*)mapping->data;
111413f72f0SBarry Smith 
112413f72f0SBarry Smith   PetscFunctionBegin;
113413f72f0SBarry Smith   if (!map) PetscFunctionReturn(0);
114413f72f0SBarry Smith   PetscHashIDestroy(map->globalht);
115413f72f0SBarry Smith   ierr = PetscFree(mapping->data);CHKERRQ(ierr);
116413f72f0SBarry Smith   PetscFunctionReturn(0);
117413f72f0SBarry Smith }
118413f72f0SBarry Smith 
119413f72f0SBarry Smith #define GTOLTYPE _Basic
120413f72f0SBarry Smith #define GTOLNAME _Basic
121541bf97eSAdrian Croucher #define GTOLBS mapping->bs
122413f72f0SBarry Smith #define GTOL(g, local) do {                  \
123413f72f0SBarry Smith     local = map->globals[g/bs - start];      \
124413f72f0SBarry Smith     local = bs*local + (g % bs);             \
125413f72f0SBarry Smith   } while (0)
126541bf97eSAdrian Croucher 
127413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
128413f72f0SBarry Smith 
129413f72f0SBarry Smith #define GTOLTYPE _Basic
130413f72f0SBarry Smith #define GTOLNAME Block_Basic
131541bf97eSAdrian Croucher #define GTOLBS 1
132413f72f0SBarry Smith #define GTOL(g, local) do {                  \
133413f72f0SBarry Smith     local = map->globals[g - start];         \
134413f72f0SBarry Smith   } while (0)
135413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
136413f72f0SBarry Smith 
137413f72f0SBarry Smith #define GTOLTYPE _Hash
138413f72f0SBarry Smith #define GTOLNAME _Hash
139541bf97eSAdrian Croucher #define GTOLBS mapping->bs
140413f72f0SBarry Smith #define GTOL(g, local) do {                  \
141413f72f0SBarry Smith     PetscHashIMap(map->globalht,g/bs,local); \
142413f72f0SBarry Smith     local = bs*local + (g % bs);             \
143413f72f0SBarry Smith    } while (0)
144413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
145413f72f0SBarry Smith 
146413f72f0SBarry Smith #define GTOLTYPE _Hash
147413f72f0SBarry Smith #define GTOLNAME Block_Hash
148541bf97eSAdrian Croucher #define GTOLBS 1
149413f72f0SBarry Smith #define GTOL(g, local) do {                  \
150413f72f0SBarry Smith     PetscHashIMap(map->globalht,g,local);    \
151413f72f0SBarry Smith   } while (0)
152413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
153413f72f0SBarry Smith 
1546658fb44Sstefano_zampini /*@
1556658fb44Sstefano_zampini     ISLocalToGlobalMappingDuplicate - Duplicates the local to global mapping object
1566658fb44Sstefano_zampini 
1576658fb44Sstefano_zampini     Not Collective
1586658fb44Sstefano_zampini 
1596658fb44Sstefano_zampini     Input Parameter:
1606658fb44Sstefano_zampini .   ltog - local to global mapping
1616658fb44Sstefano_zampini 
1626658fb44Sstefano_zampini     Output Parameter:
1636658fb44Sstefano_zampini .   nltog - the duplicated local to global mapping
1646658fb44Sstefano_zampini 
1656658fb44Sstefano_zampini     Level: advanced
1666658fb44Sstefano_zampini 
1676658fb44Sstefano_zampini     Concepts: mapping^local to global
1686658fb44Sstefano_zampini 
1696658fb44Sstefano_zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
1706658fb44Sstefano_zampini @*/
1716658fb44Sstefano_zampini PetscErrorCode  ISLocalToGlobalMappingDuplicate(ISLocalToGlobalMapping ltog,ISLocalToGlobalMapping* nltog)
1726658fb44Sstefano_zampini {
1736658fb44Sstefano_zampini   PetscErrorCode ierr;
1746658fb44Sstefano_zampini 
1756658fb44Sstefano_zampini   PetscFunctionBegin;
1766658fb44Sstefano_zampini   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
1776658fb44Sstefano_zampini   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)ltog),ltog->bs,ltog->n,ltog->indices,PETSC_COPY_VALUES,nltog);CHKERRQ(ierr);
1786658fb44Sstefano_zampini   PetscFunctionReturn(0);
1796658fb44Sstefano_zampini }
1806658fb44Sstefano_zampini 
181565245c5SBarry Smith /*@
182107e9a97SBarry Smith     ISLocalToGlobalMappingGetSize - Gets the local size of a local to global mapping
1833b9aefa3SBarry Smith 
1843b9aefa3SBarry Smith     Not Collective
1853b9aefa3SBarry Smith 
1863b9aefa3SBarry Smith     Input Parameter:
1873b9aefa3SBarry Smith .   ltog - local to global mapping
1883b9aefa3SBarry Smith 
1893b9aefa3SBarry Smith     Output Parameter:
190107e9a97SBarry Smith .   n - the number of entries in the local mapping, ISLocalToGlobalMappingGetIndices() returns an array of this length
1913b9aefa3SBarry Smith 
1923b9aefa3SBarry Smith     Level: advanced
1933b9aefa3SBarry Smith 
194273d9f13SBarry Smith     Concepts: mapping^local to global
1953b9aefa3SBarry Smith 
1963b9aefa3SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
1973b9aefa3SBarry Smith @*/
1987087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping mapping,PetscInt *n)
1993b9aefa3SBarry Smith {
2003b9aefa3SBarry Smith   PetscFunctionBegin;
2010700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
2024482741eSBarry Smith   PetscValidIntPointer(n,2);
203107e9a97SBarry Smith   *n = mapping->bs*mapping->n;
2043b9aefa3SBarry Smith   PetscFunctionReturn(0);
2053b9aefa3SBarry Smith }
2063b9aefa3SBarry Smith 
2075a5d4f66SBarry Smith /*@C
2085a5d4f66SBarry Smith     ISLocalToGlobalMappingView - View a local to global mapping
2095a5d4f66SBarry Smith 
210b9cd556bSLois Curfman McInnes     Not Collective
211b9cd556bSLois Curfman McInnes 
2125a5d4f66SBarry Smith     Input Parameters:
2133b9aefa3SBarry Smith +   ltog - local to global mapping
2143b9aefa3SBarry Smith -   viewer - viewer
2155a5d4f66SBarry Smith 
216a997ad1aSLois Curfman McInnes     Level: advanced
217a997ad1aSLois Curfman McInnes 
218273d9f13SBarry Smith     Concepts: mapping^local to global
2195a5d4f66SBarry Smith 
2205a5d4f66SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
2215a5d4f66SBarry Smith @*/
2227087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingView(ISLocalToGlobalMapping mapping,PetscViewer viewer)
2235a5d4f66SBarry Smith {
22432dcc486SBarry Smith   PetscInt       i;
22532dcc486SBarry Smith   PetscMPIInt    rank;
226ace3abfcSBarry Smith   PetscBool      iascii;
2276849ba73SBarry Smith   PetscErrorCode ierr;
2285a5d4f66SBarry Smith 
2295a5d4f66SBarry Smith   PetscFunctionBegin;
2300700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
2313050cee2SBarry Smith   if (!viewer) {
232ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mapping),&viewer);CHKERRQ(ierr);
2333050cee2SBarry Smith   }
2340700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
2355a5d4f66SBarry Smith 
236ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mapping),&rank);CHKERRQ(ierr);
237251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
23832077d6dSBarry Smith   if (iascii) {
23998c3331eSBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)mapping,viewer);CHKERRQ(ierr);
2401575c14dSBarry Smith     ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
2415a5d4f66SBarry Smith     for (i=0; i<mapping->n; i++) {
2427904a332SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D %D\n",rank,i,mapping->indices[i]);CHKERRQ(ierr);
2436831982aSBarry Smith     }
244b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
2451575c14dSBarry Smith     ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
2461575c14dSBarry Smith   }
2475a5d4f66SBarry Smith   PetscFunctionReturn(0);
2485a5d4f66SBarry Smith }
2495a5d4f66SBarry Smith 
2501f428162SBarry Smith /*@
2512bdab257SBarry Smith     ISLocalToGlobalMappingCreateIS - Creates a mapping between a local (0 to n)
2522bdab257SBarry Smith     ordering and a global parallel ordering.
2532bdab257SBarry Smith 
2540f5bd95cSBarry Smith     Not collective
255b9cd556bSLois Curfman McInnes 
256a997ad1aSLois Curfman McInnes     Input Parameter:
2578c03b21aSDmitry Karpeev .   is - index set containing the global numbers for each local number
2582bdab257SBarry Smith 
259a997ad1aSLois Curfman McInnes     Output Parameter:
2602bdab257SBarry Smith .   mapping - new mapping data structure
2612bdab257SBarry Smith 
262f0413b6fSBarry Smith     Notes: the block size of the IS determines the block size of the mapping
263a997ad1aSLois Curfman McInnes     Level: advanced
264a997ad1aSLois Curfman McInnes 
265273d9f13SBarry Smith     Concepts: mapping^local to global
2662bdab257SBarry Smith 
2677e99dc12SLawrence Mitchell .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetFromOptions()
2682bdab257SBarry Smith @*/
2697087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingCreateIS(IS is,ISLocalToGlobalMapping *mapping)
2702bdab257SBarry Smith {
2716849ba73SBarry Smith   PetscErrorCode ierr;
2723bbf0e92SBarry Smith   PetscInt       n,bs;
2735d0c19d7SBarry Smith   const PetscInt *indices;
2742bdab257SBarry Smith   MPI_Comm       comm;
2753bbf0e92SBarry Smith   PetscBool      isblock;
2763a40ed3dSBarry Smith 
2773a40ed3dSBarry Smith   PetscFunctionBegin;
2780700a824SBarry Smith   PetscValidHeaderSpecific(is,IS_CLASSID,1);
2794482741eSBarry Smith   PetscValidPointer(mapping,2);
2802bdab257SBarry Smith 
2812bdab257SBarry Smith   ierr = PetscObjectGetComm((PetscObject)is,&comm);CHKERRQ(ierr);
2823b9aefa3SBarry Smith   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
2833bbf0e92SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)is,ISBLOCK,&isblock);CHKERRQ(ierr);
2846006e8d2SBarry Smith   if (!isblock) {
285f0413b6fSBarry Smith     ierr = ISGetIndices(is,&indices);CHKERRQ(ierr);
286f0413b6fSBarry Smith     ierr = ISLocalToGlobalMappingCreate(comm,1,n,indices,PETSC_COPY_VALUES,mapping);CHKERRQ(ierr);
2872bdab257SBarry Smith     ierr = ISRestoreIndices(is,&indices);CHKERRQ(ierr);
2886006e8d2SBarry Smith   } else {
2896006e8d2SBarry Smith     ierr = ISGetBlockSize(is,&bs);CHKERRQ(ierr);
290f0413b6fSBarry Smith     ierr = ISBlockGetIndices(is,&indices);CHKERRQ(ierr);
29128bc9809SBarry Smith     ierr = ISLocalToGlobalMappingCreate(comm,bs,n/bs,indices,PETSC_COPY_VALUES,mapping);CHKERRQ(ierr);
292f0413b6fSBarry Smith     ierr = ISBlockRestoreIndices(is,&indices);CHKERRQ(ierr);
2936006e8d2SBarry Smith   }
2943a40ed3dSBarry Smith   PetscFunctionReturn(0);
2952bdab257SBarry Smith }
2965a5d4f66SBarry Smith 
297a4d96a55SJed Brown /*@C
298a4d96a55SJed Brown     ISLocalToGlobalMappingCreateSF - Creates a mapping between a local (0 to n)
299a4d96a55SJed Brown     ordering and a global parallel ordering.
300a4d96a55SJed Brown 
301a4d96a55SJed Brown     Collective
302a4d96a55SJed Brown 
303a4d96a55SJed Brown     Input Parameter:
304a4d96a55SJed Brown +   sf - star forest mapping contiguous local indices to (rank, offset)
305a4d96a55SJed Brown -   start - first global index on this process
306a4d96a55SJed Brown 
307a4d96a55SJed Brown     Output Parameter:
308a4d96a55SJed Brown .   mapping - new mapping data structure
309a4d96a55SJed Brown 
310a4d96a55SJed Brown     Level: advanced
311a4d96a55SJed Brown 
312a4d96a55SJed Brown     Concepts: mapping^local to global
313a4d96a55SJed Brown 
3147e99dc12SLawrence Mitchell .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingSetFromOptions()
315a4d96a55SJed Brown @*/
316a4d96a55SJed Brown PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF sf,PetscInt start,ISLocalToGlobalMapping *mapping)
317a4d96a55SJed Brown {
318a4d96a55SJed Brown   PetscErrorCode ierr;
319a4d96a55SJed Brown   PetscInt       i,maxlocal,nroots,nleaves,*globals,*ltog;
320a4d96a55SJed Brown   const PetscInt *ilocal;
321a4d96a55SJed Brown   MPI_Comm       comm;
322a4d96a55SJed Brown 
323a4d96a55SJed Brown   PetscFunctionBegin;
324a4d96a55SJed Brown   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
325a4d96a55SJed Brown   PetscValidPointer(mapping,3);
326a4d96a55SJed Brown 
327a4d96a55SJed Brown   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
3280298fd71SBarry Smith   ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr);
329f6e5521dSKarl Rupp   if (ilocal) {
330f6e5521dSKarl Rupp     for (i=0,maxlocal=0; i<nleaves; i++) maxlocal = PetscMax(maxlocal,ilocal[i]+1);
331f6e5521dSKarl Rupp   }
332a4d96a55SJed Brown   else maxlocal = nleaves;
333785e854fSJed Brown   ierr = PetscMalloc1(nroots,&globals);CHKERRQ(ierr);
334785e854fSJed Brown   ierr = PetscMalloc1(maxlocal,&ltog);CHKERRQ(ierr);
335a4d96a55SJed Brown   for (i=0; i<nroots; i++) globals[i] = start + i;
336a4d96a55SJed Brown   for (i=0; i<maxlocal; i++) ltog[i] = -1;
337a4d96a55SJed Brown   ierr = PetscSFBcastBegin(sf,MPIU_INT,globals,ltog);CHKERRQ(ierr);
338a4d96a55SJed Brown   ierr = PetscSFBcastEnd(sf,MPIU_INT,globals,ltog);CHKERRQ(ierr);
339f0413b6fSBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,1,maxlocal,ltog,PETSC_OWN_POINTER,mapping);CHKERRQ(ierr);
340a4d96a55SJed Brown   ierr = PetscFree(globals);CHKERRQ(ierr);
341a4d96a55SJed Brown   PetscFunctionReturn(0);
342a4d96a55SJed Brown }
343b46b645bSBarry Smith 
34463fa5c83Sstefano_zampini /*@
34563fa5c83Sstefano_zampini     ISLocalToGlobalMappingSetBlockSize - Sets the blocksize of the mapping
34663fa5c83Sstefano_zampini 
34763fa5c83Sstefano_zampini     Not collective
34863fa5c83Sstefano_zampini 
34963fa5c83Sstefano_zampini     Input Parameters:
35063fa5c83Sstefano_zampini .   mapping - mapping data structure
35163fa5c83Sstefano_zampini .   bs - the blocksize
35263fa5c83Sstefano_zampini 
35363fa5c83Sstefano_zampini     Level: advanced
35463fa5c83Sstefano_zampini 
35563fa5c83Sstefano_zampini     Concepts: mapping^local to global
35663fa5c83Sstefano_zampini 
35763fa5c83Sstefano_zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS()
35863fa5c83Sstefano_zampini @*/
35963fa5c83Sstefano_zampini PetscErrorCode  ISLocalToGlobalMappingSetBlockSize(ISLocalToGlobalMapping mapping,PetscInt bs)
36063fa5c83Sstefano_zampini {
361a59f3c4dSstefano_zampini   PetscInt       *nid;
362a59f3c4dSstefano_zampini   const PetscInt *oid;
363a59f3c4dSstefano_zampini   PetscInt       i,cn,on,obs,nn;
36463fa5c83Sstefano_zampini   PetscErrorCode ierr;
36563fa5c83Sstefano_zampini 
36663fa5c83Sstefano_zampini   PetscFunctionBegin;
36763fa5c83Sstefano_zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
36863fa5c83Sstefano_zampini   if (bs < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid block size %D",bs);
36963fa5c83Sstefano_zampini   if (bs == mapping->bs) PetscFunctionReturn(0);
37063fa5c83Sstefano_zampini   on  = mapping->n;
37163fa5c83Sstefano_zampini   obs = mapping->bs;
37263fa5c83Sstefano_zampini   oid = mapping->indices;
37363fa5c83Sstefano_zampini   nn  = (on*obs)/bs;
37463fa5c83Sstefano_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);
375a59f3c4dSstefano_zampini 
37663fa5c83Sstefano_zampini   ierr = PetscMalloc1(nn,&nid);CHKERRQ(ierr);
377a59f3c4dSstefano_zampini   ierr = ISLocalToGlobalMappingGetIndices(mapping,&oid);CHKERRQ(ierr);
378a59f3c4dSstefano_zampini   for (i=0;i<nn;i++) {
379a59f3c4dSstefano_zampini     PetscInt j;
380a59f3c4dSstefano_zampini     for (j=0,cn=0;j<bs-1;j++) {
381a59f3c4dSstefano_zampini       if (oid[i*bs+j] < 0) { cn++; continue; }
382a59f3c4dSstefano_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]);
383a59f3c4dSstefano_zampini     }
384a59f3c4dSstefano_zampini     if (oid[i*bs+j] < 0) cn++;
3858b7cb0e6Sstefano_zampini     if (cn) {
386a59f3c4dSstefano_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);
387a59f3c4dSstefano_zampini       nid[i] = -1;
3888b7cb0e6Sstefano_zampini     } else {
389a59f3c4dSstefano_zampini       nid[i] = oid[i*bs]/bs;
39063fa5c83Sstefano_zampini     }
39163fa5c83Sstefano_zampini   }
392a59f3c4dSstefano_zampini   ierr = ISLocalToGlobalMappingRestoreIndices(mapping,&oid);CHKERRQ(ierr);
393a59f3c4dSstefano_zampini 
39463fa5c83Sstefano_zampini   mapping->n           = nn;
39563fa5c83Sstefano_zampini   mapping->bs          = bs;
39663fa5c83Sstefano_zampini   ierr                 = PetscFree(mapping->indices);CHKERRQ(ierr);
39763fa5c83Sstefano_zampini   mapping->indices     = nid;
398c9345713Sstefano_zampini   mapping->globalstart = 0;
399c9345713Sstefano_zampini   mapping->globalend   = 0;
400413f72f0SBarry Smith   if (mapping->ops->destroy) {
401413f72f0SBarry Smith     ierr = (*mapping->ops->destroy)(mapping);CHKERRQ(ierr);
402413f72f0SBarry Smith   }
40363fa5c83Sstefano_zampini   PetscFunctionReturn(0);
40463fa5c83Sstefano_zampini }
40563fa5c83Sstefano_zampini 
40645b6f7e9SBarry Smith /*@
40745b6f7e9SBarry Smith     ISLocalToGlobalMappingGetBlockSize - Gets the blocksize of the mapping
40845b6f7e9SBarry Smith     ordering and a global parallel ordering.
40945b6f7e9SBarry Smith 
41045b6f7e9SBarry Smith     Not Collective
41145b6f7e9SBarry Smith 
41245b6f7e9SBarry Smith     Input Parameters:
41345b6f7e9SBarry Smith .   mapping - mapping data structure
41445b6f7e9SBarry Smith 
41545b6f7e9SBarry Smith     Output Parameter:
41645b6f7e9SBarry Smith .   bs - the blocksize
41745b6f7e9SBarry Smith 
41845b6f7e9SBarry Smith     Level: advanced
41945b6f7e9SBarry Smith 
42045b6f7e9SBarry Smith     Concepts: mapping^local to global
42145b6f7e9SBarry Smith 
42245b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS()
42345b6f7e9SBarry Smith @*/
42445b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping mapping,PetscInt *bs)
42545b6f7e9SBarry Smith {
42645b6f7e9SBarry Smith   PetscFunctionBegin;
427cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
42845b6f7e9SBarry Smith   *bs = mapping->bs;
42945b6f7e9SBarry Smith   PetscFunctionReturn(0);
43045b6f7e9SBarry Smith }
43145b6f7e9SBarry Smith 
432ba5bb76aSSatish Balay /*@
43390f02eecSBarry Smith     ISLocalToGlobalMappingCreate - Creates a mapping between a local (0 to n)
43490f02eecSBarry Smith     ordering and a global parallel ordering.
4352362add9SBarry Smith 
43689d82c54SBarry Smith     Not Collective, but communicator may have more than one process
437b9cd556bSLois Curfman McInnes 
4382362add9SBarry Smith     Input Parameters:
43989d82c54SBarry Smith +   comm - MPI communicator
440f0413b6fSBarry Smith .   bs - the block size
44128bc9809SBarry Smith .   n - the number of local elements divided by the block size, or equivalently the number of block indices
44228bc9809SBarry 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
443d5ad8652SBarry Smith -   mode - see PetscCopyMode
4442362add9SBarry Smith 
445a997ad1aSLois Curfman McInnes     Output Parameter:
44690f02eecSBarry Smith .   mapping - new mapping data structure
4472362add9SBarry Smith 
448f0413b6fSBarry 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
449413f72f0SBarry Smith 
4509a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
451413f72f0SBarry 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.
452413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
453413f72f0SBarry Smith 
454a997ad1aSLois Curfman McInnes     Level: advanced
455a997ad1aSLois Curfman McInnes 
456273d9f13SBarry Smith     Concepts: mapping^local to global
4572362add9SBarry Smith 
458413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingSetFromOptions(), ISLOCALTOGLOBALMAPPINGBASIC, ISLOCALTOGLOBALMAPPINGHASH
459413f72f0SBarry Smith           ISLocalToGlobalMappingSetType(), ISLocalToGlobalMappingType
4602362add9SBarry Smith @*/
46160c7cefcSBarry Smith PetscErrorCode  ISLocalToGlobalMappingCreate(MPI_Comm comm,PetscInt bs,PetscInt n,const PetscInt indices[],PetscCopyMode mode,ISLocalToGlobalMapping *mapping)
4622362add9SBarry Smith {
4636849ba73SBarry Smith   PetscErrorCode ierr;
46432dcc486SBarry Smith   PetscInt       *in;
465b46b645bSBarry Smith 
466b46b645bSBarry Smith   PetscFunctionBegin;
46773911063SBarry Smith   if (n) PetscValidIntPointer(indices,3);
4684482741eSBarry Smith   PetscValidPointer(mapping,4);
469b46b645bSBarry Smith 
4700298fd71SBarry Smith   *mapping = NULL;
471607a6623SBarry Smith   ierr = ISInitializePackage();CHKERRQ(ierr);
4722362add9SBarry Smith 
47373107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(*mapping,IS_LTOGM_CLASSID,"ISLocalToGlobalMapping","Local to global mapping","IS",
47460c7cefcSBarry Smith                            comm,ISLocalToGlobalMappingDestroy,ISLocalToGlobalMappingView);CHKERRQ(ierr);
475d4bb536fSBarry Smith   (*mapping)->n             = n;
476f0413b6fSBarry Smith   (*mapping)->bs            = bs;
477268a049cSStefano Zampini   (*mapping)->info_cached   = PETSC_FALSE;
478268a049cSStefano Zampini   (*mapping)->info_free     = PETSC_FALSE;
479268a049cSStefano Zampini   (*mapping)->info_procs    = NULL;
480268a049cSStefano Zampini   (*mapping)->info_numprocs = NULL;
481268a049cSStefano Zampini   (*mapping)->info_indices  = NULL;
482413f72f0SBarry Smith 
483413f72f0SBarry Smith   (*mapping)->ops->globaltolocalmappingapply      = NULL;
484413f72f0SBarry Smith   (*mapping)->ops->globaltolocalmappingapplyblock = NULL;
485413f72f0SBarry Smith   (*mapping)->ops->destroy                        = NULL;
486d5ad8652SBarry Smith   if (mode == PETSC_COPY_VALUES) {
487785e854fSJed Brown     ierr = PetscMalloc1(n,&in);CHKERRQ(ierr);
488d5ad8652SBarry Smith     ierr = PetscMemcpy(in,indices,n*sizeof(PetscInt));CHKERRQ(ierr);
489d5ad8652SBarry Smith     (*mapping)->indices = in;
4906389a1a1SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt));CHKERRQ(ierr);
4916389a1a1SBarry Smith   } else if (mode == PETSC_OWN_POINTER) {
4926389a1a1SBarry Smith     (*mapping)->indices = (PetscInt*)indices;
4936389a1a1SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt));CHKERRQ(ierr);
4946389a1a1SBarry Smith   }
49560c7cefcSBarry Smith   else SETERRQ(comm,PETSC_ERR_SUP,"Cannot currently use PETSC_USE_POINTER");
4963a40ed3dSBarry Smith   PetscFunctionReturn(0);
4972362add9SBarry Smith }
4982362add9SBarry Smith 
499413f72f0SBarry Smith PetscFunctionList ISLocalToGlobalMappingList = NULL;
500413f72f0SBarry Smith 
501413f72f0SBarry Smith 
50290f02eecSBarry Smith /*@
5037e99dc12SLawrence Mitchell    ISLocalToGlobalMappingSetFromOptions - Set mapping options from the options database.
5047e99dc12SLawrence Mitchell 
5057e99dc12SLawrence Mitchell    Not collective
5067e99dc12SLawrence Mitchell 
5077e99dc12SLawrence Mitchell    Input Parameters:
5087e99dc12SLawrence Mitchell .  mapping - mapping data structure
5097e99dc12SLawrence Mitchell 
5107e99dc12SLawrence Mitchell    Level: advanced
5117e99dc12SLawrence Mitchell 
5127e99dc12SLawrence Mitchell @*/
5137e99dc12SLawrence Mitchell PetscErrorCode ISLocalToGlobalMappingSetFromOptions(ISLocalToGlobalMapping mapping)
5147e99dc12SLawrence Mitchell {
5157e99dc12SLawrence Mitchell   PetscErrorCode             ierr;
516413f72f0SBarry Smith   char                       type[256];
517413f72f0SBarry Smith   ISLocalToGlobalMappingType defaulttype = "Not set";
5187e99dc12SLawrence Mitchell   PetscBool                  flg;
5197e99dc12SLawrence Mitchell 
5207e99dc12SLawrence Mitchell   PetscFunctionBegin;
5217e99dc12SLawrence Mitchell   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
522413f72f0SBarry Smith   ierr = ISLocalToGlobalMappingRegisterAll();CHKERRQ(ierr);
5237e99dc12SLawrence Mitchell   ierr = PetscObjectOptionsBegin((PetscObject)mapping);CHKERRQ(ierr);
524413f72f0SBarry Smith   ierr = PetscOptionsFList("-islocaltoglobalmapping_type","ISLocalToGlobalMapping method","ISLocalToGlobalMappingSetType",ISLocalToGlobalMappingList,(char*)(((PetscObject)mapping)->type_name) ? ((PetscObject)mapping)->type_name : defaulttype,type,256,&flg);CHKERRQ(ierr);
525413f72f0SBarry Smith   if (flg) {
526413f72f0SBarry Smith     ierr = ISLocalToGlobalMappingSetType(mapping,type);CHKERRQ(ierr);
527413f72f0SBarry Smith   }
5287e99dc12SLawrence Mitchell   ierr = PetscOptionsEnd();CHKERRQ(ierr);
5297e99dc12SLawrence Mitchell   PetscFunctionReturn(0);
5307e99dc12SLawrence Mitchell }
5317e99dc12SLawrence Mitchell 
5327e99dc12SLawrence Mitchell /*@
53390f02eecSBarry Smith    ISLocalToGlobalMappingDestroy - Destroys a mapping between a local (0 to n)
53490f02eecSBarry Smith    ordering and a global parallel ordering.
53590f02eecSBarry Smith 
5360f5bd95cSBarry Smith    Note Collective
537b9cd556bSLois Curfman McInnes 
53890f02eecSBarry Smith    Input Parameters:
53990f02eecSBarry Smith .  mapping - mapping data structure
54090f02eecSBarry Smith 
541a997ad1aSLois Curfman McInnes    Level: advanced
542a997ad1aSLois Curfman McInnes 
5433acfe500SLois Curfman McInnes .seealso: ISLocalToGlobalMappingCreate()
54490f02eecSBarry Smith @*/
5456bf464f9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping *mapping)
54690f02eecSBarry Smith {
547dfbe8321SBarry Smith   PetscErrorCode ierr;
5485fd66863SKarl Rupp 
5493a40ed3dSBarry Smith   PetscFunctionBegin;
5506bf464f9SBarry Smith   if (!*mapping) PetscFunctionReturn(0);
5516bf464f9SBarry Smith   PetscValidHeaderSpecific((*mapping),IS_LTOGM_CLASSID,1);
552997056adSBarry Smith   if (--((PetscObject)(*mapping))->refct > 0) {*mapping = 0;PetscFunctionReturn(0);}
5536bf464f9SBarry Smith   ierr = PetscFree((*mapping)->indices);CHKERRQ(ierr);
554268a049cSStefano Zampini   ierr = PetscFree((*mapping)->info_procs);CHKERRQ(ierr);
555268a049cSStefano Zampini   ierr = PetscFree((*mapping)->info_numprocs);CHKERRQ(ierr);
556268a049cSStefano Zampini   if ((*mapping)->info_indices) {
557268a049cSStefano Zampini     PetscInt i;
558268a049cSStefano Zampini 
559268a049cSStefano Zampini     ierr = PetscFree(((*mapping)->info_indices)[0]);CHKERRQ(ierr);
560268a049cSStefano Zampini     for (i=1; i<(*mapping)->info_nproc; i++) {
561268a049cSStefano Zampini       ierr = PetscFree(((*mapping)->info_indices)[i]);CHKERRQ(ierr);
562268a049cSStefano Zampini     }
563268a049cSStefano Zampini     ierr = PetscFree((*mapping)->info_indices);CHKERRQ(ierr);
564268a049cSStefano Zampini   }
565413f72f0SBarry Smith   if ((*mapping)->ops->destroy) {
566413f72f0SBarry Smith     ierr = (*(*mapping)->ops->destroy)(*mapping);CHKERRQ(ierr);
567413f72f0SBarry Smith   }
568d38fa0fbSBarry Smith   ierr     = PetscHeaderDestroy(mapping);CHKERRQ(ierr);
569992144d0SBarry Smith   *mapping = 0;
5703a40ed3dSBarry Smith   PetscFunctionReturn(0);
57190f02eecSBarry Smith }
57290f02eecSBarry Smith 
57390f02eecSBarry Smith /*@
5743acfe500SLois Curfman McInnes     ISLocalToGlobalMappingApplyIS - Creates from an IS in the local numbering
5753acfe500SLois Curfman McInnes     a new index set using the global numbering defined in an ISLocalToGlobalMapping
5763acfe500SLois Curfman McInnes     context.
57790f02eecSBarry Smith 
578b9cd556bSLois Curfman McInnes     Not collective
579b9cd556bSLois Curfman McInnes 
58090f02eecSBarry Smith     Input Parameters:
581b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
582b9cd556bSLois Curfman McInnes -   is - index set in local numbering
58390f02eecSBarry Smith 
58490f02eecSBarry Smith     Output Parameters:
58590f02eecSBarry Smith .   newis - index set in global numbering
58690f02eecSBarry Smith 
587a997ad1aSLois Curfman McInnes     Level: advanced
588a997ad1aSLois Curfman McInnes 
589273d9f13SBarry Smith     Concepts: mapping^local to global
5903acfe500SLois Curfman McInnes 
59190f02eecSBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),
592d4bb536fSBarry Smith           ISLocalToGlobalMappingDestroy(), ISGlobalToLocalMappingApply()
59390f02eecSBarry Smith @*/
5947087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping,IS is,IS *newis)
59590f02eecSBarry Smith {
5966849ba73SBarry Smith   PetscErrorCode ierr;
597e24637baSBarry Smith   PetscInt       n,*idxout;
5985d0c19d7SBarry Smith   const PetscInt *idxin;
5993a40ed3dSBarry Smith 
6003a40ed3dSBarry Smith   PetscFunctionBegin;
6010700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
6020700a824SBarry Smith   PetscValidHeaderSpecific(is,IS_CLASSID,2);
6034482741eSBarry Smith   PetscValidPointer(newis,3);
60490f02eecSBarry Smith 
6053b9aefa3SBarry Smith   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
60690f02eecSBarry Smith   ierr = ISGetIndices(is,&idxin);CHKERRQ(ierr);
607785e854fSJed Brown   ierr = PetscMalloc1(n,&idxout);CHKERRQ(ierr);
608e24637baSBarry Smith   ierr = ISLocalToGlobalMappingApply(mapping,n,idxin,idxout);CHKERRQ(ierr);
6093b9aefa3SBarry Smith   ierr = ISRestoreIndices(is,&idxin);CHKERRQ(ierr);
610543f3098SMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),n,idxout,PETSC_OWN_POINTER,newis);CHKERRQ(ierr);
6113a40ed3dSBarry Smith   PetscFunctionReturn(0);
61290f02eecSBarry Smith }
61390f02eecSBarry Smith 
614b89cb25eSSatish Balay /*@
6153acfe500SLois Curfman McInnes    ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering
6163acfe500SLois Curfman McInnes    and converts them to the global numbering.
61790f02eecSBarry Smith 
618b9cd556bSLois Curfman McInnes    Not collective
619b9cd556bSLois Curfman McInnes 
620bb25748dSBarry Smith    Input Parameters:
621b9cd556bSLois Curfman McInnes +  mapping - the local to global mapping context
622bb25748dSBarry Smith .  N - number of integers
623b9cd556bSLois Curfman McInnes -  in - input indices in local numbering
624bb25748dSBarry Smith 
625bb25748dSBarry Smith    Output Parameter:
626bb25748dSBarry Smith .  out - indices in global numbering
627bb25748dSBarry Smith 
628b9cd556bSLois Curfman McInnes    Notes:
629b9cd556bSLois Curfman McInnes    The in and out array parameters may be identical.
630d4bb536fSBarry Smith 
631a997ad1aSLois Curfman McInnes    Level: advanced
632a997ad1aSLois Curfman McInnes 
63345b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
6340752156aSBarry Smith           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
635d4bb536fSBarry Smith           AOPetscToApplication(), ISGlobalToLocalMappingApply()
636bb25748dSBarry Smith 
637273d9f13SBarry Smith     Concepts: mapping^local to global
638afcb2eb5SJed Brown @*/
639afcb2eb5SJed Brown PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
640afcb2eb5SJed Brown {
641cbc1caf0SMatthew G. Knepley   PetscInt i,bs,Nmax;
64245b6f7e9SBarry Smith 
64345b6f7e9SBarry Smith   PetscFunctionBegin;
644cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
645cbc1caf0SMatthew G. Knepley   bs   = mapping->bs;
646cbc1caf0SMatthew G. Knepley   Nmax = bs*mapping->n;
64745b6f7e9SBarry Smith   if (bs == 1) {
648cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
64945b6f7e9SBarry Smith     for (i=0; i<N; i++) {
65045b6f7e9SBarry Smith       if (in[i] < 0) {
65145b6f7e9SBarry Smith         out[i] = in[i];
65245b6f7e9SBarry Smith         continue;
65345b6f7e9SBarry Smith       }
654e24637baSBarry 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);
65545b6f7e9SBarry Smith       out[i] = idx[in[i]];
65645b6f7e9SBarry Smith     }
65745b6f7e9SBarry Smith   } else {
658cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
65945b6f7e9SBarry Smith     for (i=0; i<N; i++) {
66045b6f7e9SBarry Smith       if (in[i] < 0) {
66145b6f7e9SBarry Smith         out[i] = in[i];
66245b6f7e9SBarry Smith         continue;
66345b6f7e9SBarry Smith       }
664e24637baSBarry 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);
66545b6f7e9SBarry Smith       out[i] = idx[in[i]/bs]*bs + (in[i] % bs);
66645b6f7e9SBarry Smith     }
66745b6f7e9SBarry Smith   }
66845b6f7e9SBarry Smith   PetscFunctionReturn(0);
66945b6f7e9SBarry Smith }
67045b6f7e9SBarry Smith 
67145b6f7e9SBarry Smith /*@
6726006e8d2SBarry Smith    ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local block numbering and converts them to the global block numbering
67345b6f7e9SBarry Smith 
67445b6f7e9SBarry Smith    Not collective
67545b6f7e9SBarry Smith 
67645b6f7e9SBarry Smith    Input Parameters:
67745b6f7e9SBarry Smith +  mapping - the local to global mapping context
67845b6f7e9SBarry Smith .  N - number of integers
6796006e8d2SBarry Smith -  in - input indices in local block numbering
68045b6f7e9SBarry Smith 
68145b6f7e9SBarry Smith    Output Parameter:
6826006e8d2SBarry Smith .  out - indices in global block numbering
68345b6f7e9SBarry Smith 
68445b6f7e9SBarry Smith    Notes:
68545b6f7e9SBarry Smith    The in and out array parameters may be identical.
68645b6f7e9SBarry Smith 
6876006e8d2SBarry Smith    Example:
6886006e8d2SBarry 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
6896006e8d2SBarry Smith      (the first block) would produce 0 and the mapping applied to 1 (the second block) would produce 3.
6906006e8d2SBarry Smith 
69145b6f7e9SBarry Smith    Level: advanced
69245b6f7e9SBarry Smith 
69345b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
69445b6f7e9SBarry Smith           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
69545b6f7e9SBarry Smith           AOPetscToApplication(), ISGlobalToLocalMappingApply()
69645b6f7e9SBarry Smith 
69745b6f7e9SBarry Smith     Concepts: mapping^local to global
69845b6f7e9SBarry Smith @*/
69945b6f7e9SBarry Smith PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
70045b6f7e9SBarry Smith {
701cbc1caf0SMatthew G. Knepley 
702cbc1caf0SMatthew G. Knepley   PetscFunctionBegin;
703cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
704cbc1caf0SMatthew G. Knepley   {
705afcb2eb5SJed Brown     PetscInt i,Nmax = mapping->n;
706afcb2eb5SJed Brown     const PetscInt *idx = mapping->indices;
707d4bb536fSBarry Smith 
708afcb2eb5SJed Brown     for (i=0; i<N; i++) {
709afcb2eb5SJed Brown       if (in[i] < 0) {
710afcb2eb5SJed Brown         out[i] = in[i];
711afcb2eb5SJed Brown         continue;
712afcb2eb5SJed Brown       }
713e24637baSBarry 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);
714afcb2eb5SJed Brown       out[i] = idx[in[i]];
715afcb2eb5SJed Brown     }
716cbc1caf0SMatthew G. Knepley   }
717afcb2eb5SJed Brown   PetscFunctionReturn(0);
718afcb2eb5SJed Brown }
719d4bb536fSBarry Smith 
7207e99dc12SLawrence Mitchell /*@
721a997ad1aSLois Curfman McInnes     ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers
722a997ad1aSLois Curfman McInnes     specified with a global numbering.
723d4bb536fSBarry Smith 
724b9cd556bSLois Curfman McInnes     Not collective
725b9cd556bSLois Curfman McInnes 
726d4bb536fSBarry Smith     Input Parameters:
727b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
728d4bb536fSBarry Smith .   type - IS_GTOLM_MASK - replaces global indices with no local value with -1
729d4bb536fSBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
730d4bb536fSBarry Smith .   n - number of global indices to map
731b9cd556bSLois Curfman McInnes -   idx - global indices to map
732d4bb536fSBarry Smith 
733d4bb536fSBarry Smith     Output Parameters:
734b9cd556bSLois Curfman McInnes +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
735b9cd556bSLois Curfman McInnes -   idxout - local index of each global index, one must pass in an array long enough
736e182c471SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApply() with
7370298fd71SBarry Smith              idxout == NULL to determine the required length (returned in nout)
738e182c471SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApply()
739e182c471SBarry Smith              a second time to set the values.
740d4bb536fSBarry Smith 
741b9cd556bSLois Curfman McInnes     Notes:
7420298fd71SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
743d4bb536fSBarry Smith 
7449a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
745413f72f0SBarry 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.
746413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
7470f5bd95cSBarry Smith 
748a997ad1aSLois Curfman McInnes     Level: advanced
749a997ad1aSLois Curfman McInnes 
75032fd6b96SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
75132fd6b96SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
75232fd6b96SBarry Smith 
753273d9f13SBarry Smith     Concepts: mapping^global to local
754d4bb536fSBarry Smith 
7559d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),
756413f72f0SBarry Smith           ISLocalToGlobalMappingDestroy()
757d4bb536fSBarry Smith @*/
758413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
759d4bb536fSBarry Smith {
7609d90f715SBarry Smith   PetscErrorCode ierr;
7619d90f715SBarry Smith 
7629d90f715SBarry Smith   PetscFunctionBegin;
7639d90f715SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
764413f72f0SBarry Smith   if (!mapping->data) {
765413f72f0SBarry Smith     ierr = ISGlobalToLocalMappingSetUp(mapping);CHKERRQ(ierr);
7669d90f715SBarry Smith   }
767413f72f0SBarry Smith   ierr = (*mapping->ops->globaltolocalmappingapply)(mapping,type,n,idx,nout,idxout);CHKERRQ(ierr);
7689d90f715SBarry Smith   PetscFunctionReturn(0);
7699d90f715SBarry Smith }
7709d90f715SBarry Smith 
771d4fe737eSStefano Zampini /*@
772d4fe737eSStefano Zampini     ISGlobalToLocalMappingApplyIS - Creates from an IS in the global numbering
773d4fe737eSStefano Zampini     a new index set using the local numbering defined in an ISLocalToGlobalMapping
774d4fe737eSStefano Zampini     context.
775d4fe737eSStefano Zampini 
776d4fe737eSStefano Zampini     Not collective
777d4fe737eSStefano Zampini 
778d4fe737eSStefano Zampini     Input Parameters:
779d4fe737eSStefano Zampini +   mapping - mapping between local and global numbering
780*2785b321SStefano Zampini .   type - IS_GTOLM_MASK - replaces global indices with no local value with -1
781*2785b321SStefano Zampini            IS_GTOLM_DROP - drops the indices with no local value from the output list
782d4fe737eSStefano Zampini -   is - index set in global numbering
783d4fe737eSStefano Zampini 
784d4fe737eSStefano Zampini     Output Parameters:
785d4fe737eSStefano Zampini .   newis - index set in local numbering
786d4fe737eSStefano Zampini 
787d4fe737eSStefano Zampini     Level: advanced
788d4fe737eSStefano Zampini 
789d4fe737eSStefano Zampini     Concepts: mapping^local to global
790d4fe737eSStefano Zampini 
791d4fe737eSStefano Zampini .seealso: ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(),
792d4fe737eSStefano Zampini           ISLocalToGlobalMappingDestroy()
793d4fe737eSStefano Zampini @*/
794413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,IS is,IS *newis)
795d4fe737eSStefano Zampini {
796d4fe737eSStefano Zampini   PetscErrorCode ierr;
797d4fe737eSStefano Zampini   PetscInt       n,nout,*idxout;
798d4fe737eSStefano Zampini   const PetscInt *idxin;
799d4fe737eSStefano Zampini 
800d4fe737eSStefano Zampini   PetscFunctionBegin;
801d4fe737eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
802d4fe737eSStefano Zampini   PetscValidHeaderSpecific(is,IS_CLASSID,3);
803d4fe737eSStefano Zampini   PetscValidPointer(newis,4);
804d4fe737eSStefano Zampini 
805d4fe737eSStefano Zampini   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
806d4fe737eSStefano Zampini   ierr = ISGetIndices(is,&idxin);CHKERRQ(ierr);
807d4fe737eSStefano Zampini   if (type == IS_GTOLM_MASK) {
808d4fe737eSStefano Zampini     ierr = PetscMalloc1(n,&idxout);CHKERRQ(ierr);
809d4fe737eSStefano Zampini   } else {
810d4fe737eSStefano Zampini     ierr = ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,NULL);CHKERRQ(ierr);
811d4fe737eSStefano Zampini     ierr = PetscMalloc1(nout,&idxout);CHKERRQ(ierr);
812d4fe737eSStefano Zampini   }
813d4fe737eSStefano Zampini   ierr = ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,idxout);CHKERRQ(ierr);
814d4fe737eSStefano Zampini   ierr = ISRestoreIndices(is,&idxin);CHKERRQ(ierr);
815d4fe737eSStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),nout,idxout,PETSC_OWN_POINTER,newis);CHKERRQ(ierr);
816d4fe737eSStefano Zampini   PetscFunctionReturn(0);
817d4fe737eSStefano Zampini }
818d4fe737eSStefano Zampini 
8199d90f715SBarry Smith /*@
8209d90f715SBarry Smith     ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers
8219d90f715SBarry Smith     specified with a block global numbering.
8229d90f715SBarry Smith 
8239d90f715SBarry Smith     Not collective
8249d90f715SBarry Smith 
8259d90f715SBarry Smith     Input Parameters:
8269d90f715SBarry Smith +   mapping - mapping between local and global numbering
8279d90f715SBarry Smith .   type - IS_GTOLM_MASK - replaces global indices with no local value with -1
8289d90f715SBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
8299d90f715SBarry Smith .   n - number of global indices to map
8309d90f715SBarry Smith -   idx - global indices to map
8319d90f715SBarry Smith 
8329d90f715SBarry Smith     Output Parameters:
8339d90f715SBarry Smith +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
8349d90f715SBarry Smith -   idxout - local index of each global index, one must pass in an array long enough
8359d90f715SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApplyBlock() with
8369d90f715SBarry Smith              idxout == NULL to determine the required length (returned in nout)
8379d90f715SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApplyBlock()
8389d90f715SBarry Smith              a second time to set the values.
8399d90f715SBarry Smith 
8409d90f715SBarry Smith     Notes:
8419d90f715SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
8429d90f715SBarry Smith 
8439a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
844413f72f0SBarry 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.
845413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
8469d90f715SBarry Smith 
8479d90f715SBarry Smith     Level: advanced
8489d90f715SBarry Smith 
8499d90f715SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
8509d90f715SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
8519d90f715SBarry Smith 
8529d90f715SBarry Smith     Concepts: mapping^global to local
8539d90f715SBarry Smith 
8549d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(),
855413f72f0SBarry Smith           ISLocalToGlobalMappingDestroy()
8569d90f715SBarry Smith @*/
857413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,
8589d90f715SBarry Smith                                                  PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
8599d90f715SBarry Smith {
8606849ba73SBarry Smith   PetscErrorCode ierr;
861d4bb536fSBarry Smith 
8623a40ed3dSBarry Smith   PetscFunctionBegin;
8630700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
864413f72f0SBarry Smith   if (!mapping->data) {
865413f72f0SBarry Smith     ierr = ISGlobalToLocalMappingSetUp(mapping);CHKERRQ(ierr);
866d4bb536fSBarry Smith   }
867413f72f0SBarry Smith   ierr = (*mapping->ops->globaltolocalmappingapplyblock)(mapping,type,n,idx,nout,idxout);CHKERRQ(ierr);
8683a40ed3dSBarry Smith   PetscFunctionReturn(0);
869d4bb536fSBarry Smith }
87090f02eecSBarry Smith 
871413f72f0SBarry Smith 
87289d82c54SBarry Smith /*@C
8736a818285SBarry Smith     ISLocalToGlobalMappingGetBlockInfo - Gets the neighbor information for each processor and
87489d82c54SBarry Smith      each index shared by more than one processor
87589d82c54SBarry Smith 
87689d82c54SBarry Smith     Collective on ISLocalToGlobalMapping
87789d82c54SBarry Smith 
87889d82c54SBarry Smith     Input Parameters:
87989d82c54SBarry Smith .   mapping - the mapping from local to global indexing
88089d82c54SBarry Smith 
88189d82c54SBarry Smith     Output Parameter:
88289d82c54SBarry Smith +   nproc - number of processors that are connected to this one
88389d82c54SBarry Smith .   proc - neighboring processors
88407b52d57SBarry Smith .   numproc - number of indices for each subdomain (processor)
8853463a7baSJed Brown -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
88689d82c54SBarry Smith 
88789d82c54SBarry Smith     Level: advanced
88889d82c54SBarry Smith 
889273d9f13SBarry Smith     Concepts: mapping^local to global
89089d82c54SBarry Smith 
8912cfcea29SBarry Smith     Fortran Usage:
8922cfcea29SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
8932cfcea29SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
8942cfcea29SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
8952cfcea29SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
8962cfcea29SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
8972cfcea29SBarry Smith 
8982cfcea29SBarry Smith 
89907b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
90007b52d57SBarry Smith           ISLocalToGlobalMappingRestoreInfo()
90189d82c54SBarry Smith @*/
9026a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
90389d82c54SBarry Smith {
9046849ba73SBarry Smith   PetscErrorCode ierr;
905268a049cSStefano Zampini 
906268a049cSStefano Zampini   PetscFunctionBegin;
907268a049cSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
908268a049cSStefano Zampini   if (mapping->info_cached) {
909268a049cSStefano Zampini     *nproc    = mapping->info_nproc;
910268a049cSStefano Zampini     *procs    = mapping->info_procs;
911268a049cSStefano Zampini     *numprocs = mapping->info_numprocs;
912268a049cSStefano Zampini     *indices  = mapping->info_indices;
913268a049cSStefano Zampini   } else {
914268a049cSStefano Zampini     ierr = ISLocalToGlobalMappingGetBlockInfo_Private(mapping,nproc,procs,numprocs,indices);CHKERRQ(ierr);
915268a049cSStefano Zampini   }
916268a049cSStefano Zampini   PetscFunctionReturn(0);
917268a049cSStefano Zampini }
918268a049cSStefano Zampini 
919268a049cSStefano Zampini static PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
920268a049cSStefano Zampini {
921268a049cSStefano Zampini   PetscErrorCode ierr;
92297f1f81fSBarry Smith   PetscMPIInt    size,rank,tag1,tag2,tag3,*len,*source,imdex;
92332dcc486SBarry Smith   PetscInt       i,n = mapping->n,Ng,ng,max = 0,*lindices = mapping->indices;
92432dcc486SBarry Smith   PetscInt       *nprocs,*owner,nsends,*sends,j,*starts,nmax,nrecvs,*recvs,proc;
92597f1f81fSBarry Smith   PetscInt       cnt,scale,*ownedsenders,*nownedsenders,rstart,nowned;
92632dcc486SBarry Smith   PetscInt       node,nownedm,nt,*sends2,nsends2,*starts2,*lens2,*dest,nrecvs2,*starts3,*recvs2,k,*bprocs,*tmp;
92732dcc486SBarry Smith   PetscInt       first_procs,first_numprocs,*first_indices;
92889d82c54SBarry Smith   MPI_Request    *recv_waits,*send_waits;
92930dcb7c9SBarry Smith   MPI_Status     recv_status,*send_status,*recv_statuses;
930ce94432eSBarry Smith   MPI_Comm       comm;
931ace3abfcSBarry Smith   PetscBool      debug = PETSC_FALSE;
93289d82c54SBarry Smith 
93389d82c54SBarry Smith   PetscFunctionBegin;
934ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)mapping,&comm);CHKERRQ(ierr);
93524cf384cSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
93624cf384cSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
93724cf384cSBarry Smith   if (size == 1) {
93824cf384cSBarry Smith     *nproc         = 0;
9390298fd71SBarry Smith     *procs         = NULL;
94095dccacaSBarry Smith     ierr           = PetscNew(numprocs);CHKERRQ(ierr);
9411e2105dcSBarry Smith     (*numprocs)[0] = 0;
94295dccacaSBarry Smith     ierr           = PetscNew(indices);CHKERRQ(ierr);
9430298fd71SBarry Smith     (*indices)[0]  = NULL;
944268a049cSStefano Zampini     /* save info for reuse */
945268a049cSStefano Zampini     mapping->info_nproc = *nproc;
946268a049cSStefano Zampini     mapping->info_procs = *procs;
947268a049cSStefano Zampini     mapping->info_numprocs = *numprocs;
948268a049cSStefano Zampini     mapping->info_indices = *indices;
949268a049cSStefano Zampini     mapping->info_cached = PETSC_TRUE;
95024cf384cSBarry Smith     PetscFunctionReturn(0);
95124cf384cSBarry Smith   }
95224cf384cSBarry Smith 
953c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)mapping)->options,NULL,"-islocaltoglobalmappinggetinfo_debug",&debug,NULL);CHKERRQ(ierr);
95407b52d57SBarry Smith 
9553677ff5aSBarry Smith   /*
9566a818285SBarry Smith     Notes on ISLocalToGlobalMappingGetBlockInfo
9573677ff5aSBarry Smith 
9583677ff5aSBarry Smith     globally owned node - the nodes that have been assigned to this processor in global
9593677ff5aSBarry Smith            numbering, just for this routine.
9603677ff5aSBarry Smith 
9613677ff5aSBarry Smith     nontrivial globally owned node - node assigned to this processor that is on a subdomain
9623677ff5aSBarry Smith            boundary (i.e. is has more than one local owner)
9633677ff5aSBarry Smith 
9643677ff5aSBarry Smith     locally owned node - node that exists on this processors subdomain
9653677ff5aSBarry Smith 
9663677ff5aSBarry Smith     nontrivial locally owned node - node that is not in the interior (i.e. has more than one
9673677ff5aSBarry Smith            local subdomain
9683677ff5aSBarry Smith   */
96924cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag1);CHKERRQ(ierr);
97024cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag2);CHKERRQ(ierr);
97124cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag3);CHKERRQ(ierr);
97289d82c54SBarry Smith 
97389d82c54SBarry Smith   for (i=0; i<n; i++) {
97489d82c54SBarry Smith     if (lindices[i] > max) max = lindices[i];
97589d82c54SBarry Smith   }
976b2566f29SBarry Smith   ierr   = MPIU_Allreduce(&max,&Ng,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr);
97778058e43SBarry Smith   Ng++;
97889d82c54SBarry Smith   ierr   = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
97989d82c54SBarry Smith   ierr   = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
980bc8ff85bSBarry Smith   scale  = Ng/size + 1;
981a2e34c3dSBarry Smith   ng     = scale; if (rank == size-1) ng = Ng - scale*(size-1); ng = PetscMax(1,ng);
982caba0dd0SBarry Smith   rstart = scale*rank;
98389d82c54SBarry Smith 
98489d82c54SBarry Smith   /* determine ownership ranges of global indices */
985785e854fSJed Brown   ierr = PetscMalloc1(2*size,&nprocs);CHKERRQ(ierr);
98632dcc486SBarry Smith   ierr = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr);
98789d82c54SBarry Smith 
98889d82c54SBarry Smith   /* determine owners of each local node  */
989785e854fSJed Brown   ierr = PetscMalloc1(n,&owner);CHKERRQ(ierr);
99089d82c54SBarry Smith   for (i=0; i<n; i++) {
9913677ff5aSBarry Smith     proc             = lindices[i]/scale; /* processor that globally owns this index */
99227c402fcSBarry Smith     nprocs[2*proc+1] = 1;                 /* processor globally owns at least one of ours */
9933677ff5aSBarry Smith     owner[i]         = proc;
99427c402fcSBarry Smith     nprocs[2*proc]++;                     /* count of how many that processor globally owns of ours */
99589d82c54SBarry Smith   }
99627c402fcSBarry Smith   nsends = 0; for (i=0; i<size; i++) nsends += nprocs[2*i+1];
9977904a332SBarry Smith   ierr = PetscInfo1(mapping,"Number of global owners for my local data %D\n",nsends);CHKERRQ(ierr);
99889d82c54SBarry Smith 
99989d82c54SBarry Smith   /* inform other processors of number of messages and max length*/
100027c402fcSBarry Smith   ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr);
10017904a332SBarry Smith   ierr = PetscInfo1(mapping,"Number of local owners for my global data %D\n",nrecvs);CHKERRQ(ierr);
100289d82c54SBarry Smith 
100389d82c54SBarry Smith   /* post receives for owned rows */
1004785e854fSJed Brown   ierr = PetscMalloc1((2*nrecvs+1)*(nmax+1),&recvs);CHKERRQ(ierr);
1005854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&recv_waits);CHKERRQ(ierr);
100689d82c54SBarry Smith   for (i=0; i<nrecvs; i++) {
100732dcc486SBarry Smith     ierr = MPI_Irecv(recvs+2*nmax*i,2*nmax,MPIU_INT,MPI_ANY_SOURCE,tag1,comm,recv_waits+i);CHKERRQ(ierr);
100889d82c54SBarry Smith   }
100989d82c54SBarry Smith 
101089d82c54SBarry Smith   /* pack messages containing lists of local nodes to owners */
1011854ce69bSBarry Smith   ierr      = PetscMalloc1(2*n+1,&sends);CHKERRQ(ierr);
1012854ce69bSBarry Smith   ierr      = PetscMalloc1(size+1,&starts);CHKERRQ(ierr);
101389d82c54SBarry Smith   starts[0] = 0;
1014f6e5521dSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
101589d82c54SBarry Smith   for (i=0; i<n; i++) {
101689d82c54SBarry Smith     sends[starts[owner[i]]++] = lindices[i];
101730dcb7c9SBarry Smith     sends[starts[owner[i]]++] = i;
101889d82c54SBarry Smith   }
101989d82c54SBarry Smith   ierr = PetscFree(owner);CHKERRQ(ierr);
102089d82c54SBarry Smith   starts[0] = 0;
1021f6e5521dSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
102289d82c54SBarry Smith 
102389d82c54SBarry Smith   /* send the messages */
1024854ce69bSBarry Smith   ierr = PetscMalloc1(nsends+1,&send_waits);CHKERRQ(ierr);
1025854ce69bSBarry Smith   ierr = PetscMalloc1(nsends+1,&dest);CHKERRQ(ierr);
102689d82c54SBarry Smith   cnt = 0;
102789d82c54SBarry Smith   for (i=0; i<size; i++) {
102827c402fcSBarry Smith     if (nprocs[2*i]) {
102932dcc486SBarry Smith       ierr      = MPI_Isend(sends+starts[i],2*nprocs[2*i],MPIU_INT,i,tag1,comm,send_waits+cnt);CHKERRQ(ierr);
103030dcb7c9SBarry Smith       dest[cnt] = i;
103189d82c54SBarry Smith       cnt++;
103289d82c54SBarry Smith     }
103389d82c54SBarry Smith   }
103489d82c54SBarry Smith   ierr = PetscFree(starts);CHKERRQ(ierr);
103589d82c54SBarry Smith 
103689d82c54SBarry Smith   /* wait on receives */
1037854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&source);CHKERRQ(ierr);
1038854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&len);CHKERRQ(ierr);
103989d82c54SBarry Smith   cnt  = nrecvs;
1040854ce69bSBarry Smith   ierr = PetscMalloc1(ng+1,&nownedsenders);CHKERRQ(ierr);
104132dcc486SBarry Smith   ierr = PetscMemzero(nownedsenders,ng*sizeof(PetscInt));CHKERRQ(ierr);
104289d82c54SBarry Smith   while (cnt) {
104389d82c54SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
104489d82c54SBarry Smith     /* unpack receives into our local space */
104532dcc486SBarry Smith     ierr          = MPI_Get_count(&recv_status,MPIU_INT,&len[imdex]);CHKERRQ(ierr);
104689d82c54SBarry Smith     source[imdex] = recv_status.MPI_SOURCE;
104730dcb7c9SBarry Smith     len[imdex]    = len[imdex]/2;
1048caba0dd0SBarry Smith     /* count how many local owners for each of my global owned indices */
104930dcb7c9SBarry Smith     for (i=0; i<len[imdex]; i++) nownedsenders[recvs[2*imdex*nmax+2*i]-rstart]++;
105089d82c54SBarry Smith     cnt--;
105189d82c54SBarry Smith   }
105289d82c54SBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
105389d82c54SBarry Smith 
105430dcb7c9SBarry Smith   /* count how many globally owned indices are on an edge multiplied by how many processors own them. */
1055bc8ff85bSBarry Smith   nowned  = 0;
1056bc8ff85bSBarry Smith   nownedm = 0;
1057bc8ff85bSBarry Smith   for (i=0; i<ng; i++) {
1058bc8ff85bSBarry Smith     if (nownedsenders[i] > 1) {nownedm += nownedsenders[i]; nowned++;}
1059bc8ff85bSBarry Smith   }
1060bc8ff85bSBarry Smith 
1061bc8ff85bSBarry Smith   /* create single array to contain rank of all local owners of each globally owned index */
1062854ce69bSBarry Smith   ierr      = PetscMalloc1(nownedm+1,&ownedsenders);CHKERRQ(ierr);
1063854ce69bSBarry Smith   ierr      = PetscMalloc1(ng+1,&starts);CHKERRQ(ierr);
1064bc8ff85bSBarry Smith   starts[0] = 0;
1065bc8ff85bSBarry Smith   for (i=1; i<ng; i++) {
1066bc8ff85bSBarry Smith     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
1067bc8ff85bSBarry Smith     else starts[i] = starts[i-1];
1068bc8ff85bSBarry Smith   }
1069bc8ff85bSBarry Smith 
107030dcb7c9SBarry Smith   /* for each nontrival globally owned node list all arriving processors */
1071bc8ff85bSBarry Smith   for (i=0; i<nrecvs; i++) {
1072bc8ff85bSBarry Smith     for (j=0; j<len[i]; j++) {
107330dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
1074f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) ownedsenders[starts[node]++] = source[i];
1075bc8ff85bSBarry Smith     }
1076bc8ff85bSBarry Smith   }
1077bc8ff85bSBarry Smith 
107807b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
107930dcb7c9SBarry Smith     starts[0] = 0;
108030dcb7c9SBarry Smith     for (i=1; i<ng; i++) {
108130dcb7c9SBarry Smith       if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
108230dcb7c9SBarry Smith       else starts[i] = starts[i-1];
108330dcb7c9SBarry Smith     }
108430dcb7c9SBarry Smith     for (i=0; i<ng; i++) {
108530dcb7c9SBarry Smith       if (nownedsenders[i] > 1) {
10867904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"[%d] global node %D local owner processors: ",rank,i+rstart);CHKERRQ(ierr);
108730dcb7c9SBarry Smith         for (j=0; j<nownedsenders[i]; j++) {
10887904a332SBarry Smith           ierr = PetscSynchronizedPrintf(comm,"%D ",ownedsenders[starts[i]+j]);CHKERRQ(ierr);
108930dcb7c9SBarry Smith         }
109030dcb7c9SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
109130dcb7c9SBarry Smith       }
109230dcb7c9SBarry Smith     }
10930ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
109407b52d57SBarry Smith   } /* -----------------------------------  */
109530dcb7c9SBarry Smith 
10963677ff5aSBarry Smith   /* wait on original sends */
10973a96401aSBarry Smith   if (nsends) {
1098785e854fSJed Brown     ierr = PetscMalloc1(nsends,&send_status);CHKERRQ(ierr);
10993a96401aSBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
11003a96401aSBarry Smith     ierr = PetscFree(send_status);CHKERRQ(ierr);
11013a96401aSBarry Smith   }
110289d82c54SBarry Smith   ierr = PetscFree(send_waits);CHKERRQ(ierr);
11033a96401aSBarry Smith   ierr = PetscFree(sends);CHKERRQ(ierr);
11043677ff5aSBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
11053677ff5aSBarry Smith 
11063677ff5aSBarry Smith   /* pack messages to send back to local owners */
110730dcb7c9SBarry Smith   starts[0] = 0;
110830dcb7c9SBarry Smith   for (i=1; i<ng; i++) {
110930dcb7c9SBarry Smith     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
111030dcb7c9SBarry Smith     else starts[i] = starts[i-1];
111130dcb7c9SBarry Smith   }
111230dcb7c9SBarry Smith   nsends2 = nrecvs;
1113854ce69bSBarry Smith   ierr    = PetscMalloc1(nsends2+1,&nprocs);CHKERRQ(ierr); /* length of each message */
111430dcb7c9SBarry Smith   for (i=0; i<nrecvs; i++) {
111530dcb7c9SBarry Smith     nprocs[i] = 1;
111630dcb7c9SBarry Smith     for (j=0; j<len[i]; j++) {
111730dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
1118f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) nprocs[i] += 2 + nownedsenders[node];
111930dcb7c9SBarry Smith     }
112030dcb7c9SBarry Smith   }
1121f6e5521dSKarl Rupp   nt = 0;
1122f6e5521dSKarl Rupp   for (i=0; i<nsends2; i++) nt += nprocs[i];
1123f6e5521dSKarl Rupp 
1124854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,&sends2);CHKERRQ(ierr);
1125854ce69bSBarry Smith   ierr = PetscMalloc1(nsends2+1,&starts2);CHKERRQ(ierr);
1126f6e5521dSKarl Rupp 
1127f6e5521dSKarl Rupp   starts2[0] = 0;
1128f6e5521dSKarl Rupp   for (i=1; i<nsends2; i++) starts2[i] = starts2[i-1] + nprocs[i-1];
112930dcb7c9SBarry Smith   /*
113030dcb7c9SBarry Smith      Each message is 1 + nprocs[i] long, and consists of
113130dcb7c9SBarry Smith        (0) the number of nodes being sent back
113230dcb7c9SBarry Smith        (1) the local node number,
113330dcb7c9SBarry Smith        (2) the number of processors sharing it,
113430dcb7c9SBarry Smith        (3) the processors sharing it
113530dcb7c9SBarry Smith   */
113630dcb7c9SBarry Smith   for (i=0; i<nsends2; i++) {
113730dcb7c9SBarry Smith     cnt = 1;
113830dcb7c9SBarry Smith     sends2[starts2[i]] = 0;
113930dcb7c9SBarry Smith     for (j=0; j<len[i]; j++) {
114030dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
114130dcb7c9SBarry Smith       if (nownedsenders[node] > 1) {
114230dcb7c9SBarry Smith         sends2[starts2[i]]++;
114330dcb7c9SBarry Smith         sends2[starts2[i]+cnt++] = recvs[2*i*nmax+2*j+1];
114430dcb7c9SBarry Smith         sends2[starts2[i]+cnt++] = nownedsenders[node];
114532dcc486SBarry Smith         ierr = PetscMemcpy(&sends2[starts2[i]+cnt],&ownedsenders[starts[node]],nownedsenders[node]*sizeof(PetscInt));CHKERRQ(ierr);
114630dcb7c9SBarry Smith         cnt += nownedsenders[node];
114730dcb7c9SBarry Smith       }
114830dcb7c9SBarry Smith     }
114930dcb7c9SBarry Smith   }
115030dcb7c9SBarry Smith 
115130dcb7c9SBarry Smith   /* receive the message lengths */
115230dcb7c9SBarry Smith   nrecvs2 = nsends;
1153854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&lens2);CHKERRQ(ierr);
1154854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&starts3);CHKERRQ(ierr);
1155854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&recv_waits);CHKERRQ(ierr);
115630dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
1157d44834fbSBarry Smith     ierr = MPI_Irecv(&lens2[i],1,MPIU_INT,dest[i],tag2,comm,recv_waits+i);CHKERRQ(ierr);
115830dcb7c9SBarry Smith   }
1159d44834fbSBarry Smith 
11608a8e0b3aSBarry Smith   /* send the message lengths */
11618a8e0b3aSBarry Smith   for (i=0; i<nsends2; i++) {
11628a8e0b3aSBarry Smith     ierr = MPI_Send(&nprocs[i],1,MPIU_INT,source[i],tag2,comm);CHKERRQ(ierr);
11638a8e0b3aSBarry Smith   }
11648a8e0b3aSBarry Smith 
1165d44834fbSBarry Smith   /* wait on receives of lens */
11660c468ba9SBarry Smith   if (nrecvs2) {
1167785e854fSJed Brown     ierr = PetscMalloc1(nrecvs2,&recv_statuses);CHKERRQ(ierr);
1168d44834fbSBarry Smith     ierr = MPI_Waitall(nrecvs2,recv_waits,recv_statuses);CHKERRQ(ierr);
1169d44834fbSBarry Smith     ierr = PetscFree(recv_statuses);CHKERRQ(ierr);
11700c468ba9SBarry Smith   }
1171a2ea699eSBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
1172d44834fbSBarry Smith 
117330dcb7c9SBarry Smith   starts3[0] = 0;
1174d44834fbSBarry Smith   nt         = 0;
117530dcb7c9SBarry Smith   for (i=0; i<nrecvs2-1; i++) {
117630dcb7c9SBarry Smith     starts3[i+1] = starts3[i] + lens2[i];
1177d44834fbSBarry Smith     nt          += lens2[i];
117830dcb7c9SBarry Smith   }
117976466f69SStefano Zampini   if (nrecvs2) nt += lens2[nrecvs2-1];
1180d44834fbSBarry Smith 
1181854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,&recvs2);CHKERRQ(ierr);
1182854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs2+1,&recv_waits);CHKERRQ(ierr);
118352b72c4aSBarry Smith   for (i=0; i<nrecvs2; i++) {
118432dcc486SBarry Smith     ierr = MPI_Irecv(recvs2+starts3[i],lens2[i],MPIU_INT,dest[i],tag3,comm,recv_waits+i);CHKERRQ(ierr);
118530dcb7c9SBarry Smith   }
118630dcb7c9SBarry Smith 
118730dcb7c9SBarry Smith   /* send the messages */
1188854ce69bSBarry Smith   ierr = PetscMalloc1(nsends2+1,&send_waits);CHKERRQ(ierr);
118930dcb7c9SBarry Smith   for (i=0; i<nsends2; i++) {
119032dcc486SBarry Smith     ierr = MPI_Isend(sends2+starts2[i],nprocs[i],MPIU_INT,source[i],tag3,comm,send_waits+i);CHKERRQ(ierr);
119130dcb7c9SBarry Smith   }
119230dcb7c9SBarry Smith 
119330dcb7c9SBarry Smith   /* wait on receives */
11940c468ba9SBarry Smith   if (nrecvs2) {
1195785e854fSJed Brown     ierr = PetscMalloc1(nrecvs2,&recv_statuses);CHKERRQ(ierr);
119630dcb7c9SBarry Smith     ierr = MPI_Waitall(nrecvs2,recv_waits,recv_statuses);CHKERRQ(ierr);
119730dcb7c9SBarry Smith     ierr = PetscFree(recv_statuses);CHKERRQ(ierr);
11980c468ba9SBarry Smith   }
119930dcb7c9SBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
120030dcb7c9SBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
120130dcb7c9SBarry Smith 
120207b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
120330dcb7c9SBarry Smith     cnt = 0;
120430dcb7c9SBarry Smith     for (i=0; i<nrecvs2; i++) {
120530dcb7c9SBarry Smith       nt = recvs2[cnt++];
120630dcb7c9SBarry Smith       for (j=0; j<nt; j++) {
12077904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"[%d] local node %D number of subdomains %D: ",rank,recvs2[cnt],recvs2[cnt+1]);CHKERRQ(ierr);
120830dcb7c9SBarry Smith         for (k=0; k<recvs2[cnt+1]; k++) {
12097904a332SBarry Smith           ierr = PetscSynchronizedPrintf(comm,"%D ",recvs2[cnt+2+k]);CHKERRQ(ierr);
121030dcb7c9SBarry Smith         }
121130dcb7c9SBarry Smith         cnt += 2 + recvs2[cnt+1];
121230dcb7c9SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
121330dcb7c9SBarry Smith       }
121430dcb7c9SBarry Smith     }
12150ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
121607b52d57SBarry Smith   } /* -----------------------------------  */
121730dcb7c9SBarry Smith 
121830dcb7c9SBarry Smith   /* count number subdomains for each local node */
1219785e854fSJed Brown   ierr = PetscMalloc1(size,&nprocs);CHKERRQ(ierr);
122032dcc486SBarry Smith   ierr = PetscMemzero(nprocs,size*sizeof(PetscInt));CHKERRQ(ierr);
122130dcb7c9SBarry Smith   cnt  = 0;
122230dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
122330dcb7c9SBarry Smith     nt = recvs2[cnt++];
122430dcb7c9SBarry Smith     for (j=0; j<nt; j++) {
1225f6e5521dSKarl Rupp       for (k=0; k<recvs2[cnt+1]; k++) nprocs[recvs2[cnt+2+k]]++;
122630dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt+1];
122730dcb7c9SBarry Smith     }
122830dcb7c9SBarry Smith   }
122930dcb7c9SBarry Smith   nt = 0; for (i=0; i<size; i++) nt += (nprocs[i] > 0);
123030dcb7c9SBarry Smith   *nproc    = nt;
1231854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,procs);CHKERRQ(ierr);
1232854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,numprocs);CHKERRQ(ierr);
1233854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,indices);CHKERRQ(ierr);
12340298fd71SBarry Smith   for (i=0;i<nt+1;i++) (*indices)[i]=NULL;
1235785e854fSJed Brown   ierr = PetscMalloc1(size,&bprocs);CHKERRQ(ierr);
123630dcb7c9SBarry Smith   cnt  = 0;
123730dcb7c9SBarry Smith   for (i=0; i<size; i++) {
123830dcb7c9SBarry Smith     if (nprocs[i] > 0) {
123930dcb7c9SBarry Smith       bprocs[i]        = cnt;
124030dcb7c9SBarry Smith       (*procs)[cnt]    = i;
124130dcb7c9SBarry Smith       (*numprocs)[cnt] = nprocs[i];
1242785e854fSJed Brown       ierr             = PetscMalloc1(nprocs[i],&(*indices)[cnt]);CHKERRQ(ierr);
124330dcb7c9SBarry Smith       cnt++;
124430dcb7c9SBarry Smith     }
124530dcb7c9SBarry Smith   }
124630dcb7c9SBarry Smith 
124730dcb7c9SBarry Smith   /* make the list of subdomains for each nontrivial local node */
124832dcc486SBarry Smith   ierr = PetscMemzero(*numprocs,nt*sizeof(PetscInt));CHKERRQ(ierr);
124930dcb7c9SBarry Smith   cnt  = 0;
125030dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
125130dcb7c9SBarry Smith     nt = recvs2[cnt++];
125230dcb7c9SBarry Smith     for (j=0; j<nt; j++) {
1253f6e5521dSKarl Rupp       for (k=0; k<recvs2[cnt+1]; k++) (*indices)[bprocs[recvs2[cnt+2+k]]][(*numprocs)[bprocs[recvs2[cnt+2+k]]]++] = recvs2[cnt];
125430dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt+1];
125530dcb7c9SBarry Smith     }
125630dcb7c9SBarry Smith   }
125730dcb7c9SBarry Smith   ierr = PetscFree(bprocs);CHKERRQ(ierr);
125807b52d57SBarry Smith   ierr = PetscFree(recvs2);CHKERRQ(ierr);
125930dcb7c9SBarry Smith 
126007b52d57SBarry Smith   /* sort the node indexing by their global numbers */
126107b52d57SBarry Smith   nt = *nproc;
126207b52d57SBarry Smith   for (i=0; i<nt; i++) {
1263854ce69bSBarry Smith     ierr = PetscMalloc1((*numprocs)[i],&tmp);CHKERRQ(ierr);
1264f6e5521dSKarl Rupp     for (j=0; j<(*numprocs)[i]; j++) tmp[j] = lindices[(*indices)[i][j]];
126507b52d57SBarry Smith     ierr = PetscSortIntWithArray((*numprocs)[i],tmp,(*indices)[i]);CHKERRQ(ierr);
126607b52d57SBarry Smith     ierr = PetscFree(tmp);CHKERRQ(ierr);
126707b52d57SBarry Smith   }
126807b52d57SBarry Smith 
126907b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
127030dcb7c9SBarry Smith     nt = *nproc;
127130dcb7c9SBarry Smith     for (i=0; i<nt; i++) {
12727904a332SBarry Smith       ierr = PetscSynchronizedPrintf(comm,"[%d] subdomain %D number of indices %D: ",rank,(*procs)[i],(*numprocs)[i]);CHKERRQ(ierr);
127330dcb7c9SBarry Smith       for (j=0; j<(*numprocs)[i]; j++) {
12747904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"%D ",(*indices)[i][j]);CHKERRQ(ierr);
127530dcb7c9SBarry Smith       }
127630dcb7c9SBarry Smith       ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
127730dcb7c9SBarry Smith     }
12780ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
127907b52d57SBarry Smith   } /* -----------------------------------  */
128030dcb7c9SBarry Smith 
128130dcb7c9SBarry Smith   /* wait on sends */
128230dcb7c9SBarry Smith   if (nsends2) {
1283785e854fSJed Brown     ierr = PetscMalloc1(nsends2,&send_status);CHKERRQ(ierr);
128430dcb7c9SBarry Smith     ierr = MPI_Waitall(nsends2,send_waits,send_status);CHKERRQ(ierr);
128530dcb7c9SBarry Smith     ierr = PetscFree(send_status);CHKERRQ(ierr);
128630dcb7c9SBarry Smith   }
128730dcb7c9SBarry Smith 
128830dcb7c9SBarry Smith   ierr = PetscFree(starts3);CHKERRQ(ierr);
128930dcb7c9SBarry Smith   ierr = PetscFree(dest);CHKERRQ(ierr);
129030dcb7c9SBarry Smith   ierr = PetscFree(send_waits);CHKERRQ(ierr);
12913677ff5aSBarry Smith 
1292bc8ff85bSBarry Smith   ierr = PetscFree(nownedsenders);CHKERRQ(ierr);
1293bc8ff85bSBarry Smith   ierr = PetscFree(ownedsenders);CHKERRQ(ierr);
1294bc8ff85bSBarry Smith   ierr = PetscFree(starts);CHKERRQ(ierr);
129530dcb7c9SBarry Smith   ierr = PetscFree(starts2);CHKERRQ(ierr);
129630dcb7c9SBarry Smith   ierr = PetscFree(lens2);CHKERRQ(ierr);
129789d82c54SBarry Smith 
129889d82c54SBarry Smith   ierr = PetscFree(source);CHKERRQ(ierr);
129997f1f81fSBarry Smith   ierr = PetscFree(len);CHKERRQ(ierr);
130089d82c54SBarry Smith   ierr = PetscFree(recvs);CHKERRQ(ierr);
13013a96401aSBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
130230dcb7c9SBarry Smith   ierr = PetscFree(sends2);CHKERRQ(ierr);
130324cf384cSBarry Smith 
130424cf384cSBarry Smith   /* put the information about myself as the first entry in the list */
130524cf384cSBarry Smith   first_procs    = (*procs)[0];
130624cf384cSBarry Smith   first_numprocs = (*numprocs)[0];
130724cf384cSBarry Smith   first_indices  = (*indices)[0];
130824cf384cSBarry Smith   for (i=0; i<*nproc; i++) {
130924cf384cSBarry Smith     if ((*procs)[i] == rank) {
131024cf384cSBarry Smith       (*procs)[0]    = (*procs)[i];
131124cf384cSBarry Smith       (*numprocs)[0] = (*numprocs)[i];
131224cf384cSBarry Smith       (*indices)[0]  = (*indices)[i];
131324cf384cSBarry Smith       (*procs)[i]    = first_procs;
131424cf384cSBarry Smith       (*numprocs)[i] = first_numprocs;
131524cf384cSBarry Smith       (*indices)[i]  = first_indices;
131624cf384cSBarry Smith       break;
131724cf384cSBarry Smith     }
131824cf384cSBarry Smith   }
1319268a049cSStefano Zampini 
1320268a049cSStefano Zampini   /* save info for reuse */
1321268a049cSStefano Zampini   mapping->info_nproc = *nproc;
1322268a049cSStefano Zampini   mapping->info_procs = *procs;
1323268a049cSStefano Zampini   mapping->info_numprocs = *numprocs;
1324268a049cSStefano Zampini   mapping->info_indices = *indices;
1325268a049cSStefano Zampini   mapping->info_cached = PETSC_TRUE;
132689d82c54SBarry Smith   PetscFunctionReturn(0);
132789d82c54SBarry Smith }
132889d82c54SBarry Smith 
13296a818285SBarry Smith /*@C
13306a818285SBarry Smith     ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by ISLocalToGlobalMappingGetBlockInfo()
13316a818285SBarry Smith 
13326a818285SBarry Smith     Collective on ISLocalToGlobalMapping
13336a818285SBarry Smith 
13346a818285SBarry Smith     Input Parameters:
13356a818285SBarry Smith .   mapping - the mapping from local to global indexing
13366a818285SBarry Smith 
13376a818285SBarry Smith     Output Parameter:
13386a818285SBarry Smith +   nproc - number of processors that are connected to this one
13396a818285SBarry Smith .   proc - neighboring processors
13406a818285SBarry Smith .   numproc - number of indices for each processor
13416a818285SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
13426a818285SBarry Smith 
13436a818285SBarry Smith     Level: advanced
13446a818285SBarry Smith 
13456a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
13466a818285SBarry Smith           ISLocalToGlobalMappingGetInfo()
13476a818285SBarry Smith @*/
13486a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
13496a818285SBarry Smith {
13506a818285SBarry Smith   PetscErrorCode ierr;
13516a818285SBarry Smith 
13526a818285SBarry Smith   PetscFunctionBegin;
1353cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1354268a049cSStefano Zampini   if (mapping->info_free) {
13556a818285SBarry Smith     ierr = PetscFree(*numprocs);CHKERRQ(ierr);
13566a818285SBarry Smith     if (*indices) {
1357268a049cSStefano Zampini       PetscInt i;
1358268a049cSStefano Zampini 
13596a818285SBarry Smith       ierr = PetscFree((*indices)[0]);CHKERRQ(ierr);
13606a818285SBarry Smith       for (i=1; i<*nproc; i++) {
13616a818285SBarry Smith         ierr = PetscFree((*indices)[i]);CHKERRQ(ierr);
13626a818285SBarry Smith       }
13636a818285SBarry Smith       ierr = PetscFree(*indices);CHKERRQ(ierr);
13646a818285SBarry Smith     }
1365268a049cSStefano Zampini   }
1366268a049cSStefano Zampini   *nproc    = 0;
1367268a049cSStefano Zampini   *procs    = NULL;
1368268a049cSStefano Zampini   *numprocs = NULL;
1369268a049cSStefano Zampini   *indices  = NULL;
13706a818285SBarry Smith   PetscFunctionReturn(0);
13716a818285SBarry Smith }
13726a818285SBarry Smith 
13736a818285SBarry Smith /*@C
13746a818285SBarry Smith     ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each processor and
13756a818285SBarry Smith      each index shared by more than one processor
13766a818285SBarry Smith 
13776a818285SBarry Smith     Collective on ISLocalToGlobalMapping
13786a818285SBarry Smith 
13796a818285SBarry Smith     Input Parameters:
13806a818285SBarry Smith .   mapping - the mapping from local to global indexing
13816a818285SBarry Smith 
13826a818285SBarry Smith     Output Parameter:
13836a818285SBarry Smith +   nproc - number of processors that are connected to this one
13846a818285SBarry Smith .   proc - neighboring processors
13856a818285SBarry Smith .   numproc - number of indices for each subdomain (processor)
13866a818285SBarry Smith -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
13876a818285SBarry Smith 
13886a818285SBarry Smith     Level: advanced
13896a818285SBarry Smith 
13906a818285SBarry Smith     Concepts: mapping^local to global
13916a818285SBarry Smith 
13926a818285SBarry Smith     Fortran Usage:
13936a818285SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
13946a818285SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
13956a818285SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
13966a818285SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
13976a818285SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
13986a818285SBarry Smith 
13996a818285SBarry Smith 
14006a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
14016a818285SBarry Smith           ISLocalToGlobalMappingRestoreInfo()
14026a818285SBarry Smith @*/
14036a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
14046a818285SBarry Smith {
14056a818285SBarry Smith   PetscErrorCode ierr;
1406268a049cSStefano Zampini   PetscInt       **bindices = NULL,*bnumprocs = NULL,bs = mapping->bs,i,j,k;
14076a818285SBarry Smith 
14086a818285SBarry Smith   PetscFunctionBegin;
14096a818285SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1410268a049cSStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockInfo(mapping,nproc,procs,&bnumprocs,&bindices);CHKERRQ(ierr);
1411268a049cSStefano Zampini   if (bs > 1) { /* we need to expand the cached info */
1412732f65e3SBarry Smith     ierr = PetscCalloc1(*nproc,&*indices);CHKERRQ(ierr);
1413268a049cSStefano Zampini     ierr = PetscCalloc1(*nproc,&*numprocs);CHKERRQ(ierr);
14146a818285SBarry Smith     for (i=0; i<*nproc; i++) {
1415268a049cSStefano Zampini       ierr = PetscMalloc1(bs*bnumprocs[i],&(*indices)[i]);CHKERRQ(ierr);
1416268a049cSStefano Zampini       for (j=0; j<bnumprocs[i]; j++) {
14176a818285SBarry Smith         for (k=0; k<bs; k++) {
14186a818285SBarry Smith           (*indices)[i][j*bs+k] = bs*bindices[i][j] + k;
14196a818285SBarry Smith         }
14206a818285SBarry Smith       }
1421268a049cSStefano Zampini       (*numprocs)[i] = bnumprocs[i]*bs;
14226a818285SBarry Smith     }
1423268a049cSStefano Zampini     mapping->info_free = PETSC_TRUE;
1424268a049cSStefano Zampini   } else {
1425268a049cSStefano Zampini     *numprocs = bnumprocs;
1426268a049cSStefano Zampini     *indices  = bindices;
14276a818285SBarry Smith   }
14286a818285SBarry Smith   PetscFunctionReturn(0);
14296a818285SBarry Smith }
14306a818285SBarry Smith 
143107b52d57SBarry Smith /*@C
143207b52d57SBarry Smith     ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by ISLocalToGlobalMappingGetInfo()
143389d82c54SBarry Smith 
143407b52d57SBarry Smith     Collective on ISLocalToGlobalMapping
143507b52d57SBarry Smith 
143607b52d57SBarry Smith     Input Parameters:
143707b52d57SBarry Smith .   mapping - the mapping from local to global indexing
143807b52d57SBarry Smith 
143907b52d57SBarry Smith     Output Parameter:
144007b52d57SBarry Smith +   nproc - number of processors that are connected to this one
144107b52d57SBarry Smith .   proc - neighboring processors
144207b52d57SBarry Smith .   numproc - number of indices for each processor
144307b52d57SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
144407b52d57SBarry Smith 
144507b52d57SBarry Smith     Level: advanced
144607b52d57SBarry Smith 
144707b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
144807b52d57SBarry Smith           ISLocalToGlobalMappingGetInfo()
144907b52d57SBarry Smith @*/
14507087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
145107b52d57SBarry Smith {
14526849ba73SBarry Smith   PetscErrorCode ierr;
145307b52d57SBarry Smith 
145407b52d57SBarry Smith   PetscFunctionBegin;
14556a818285SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockInfo(mapping,nproc,procs,numprocs,indices);CHKERRQ(ierr);
145607b52d57SBarry Smith   PetscFunctionReturn(0);
145707b52d57SBarry Smith }
145886994e45SJed Brown 
145986994e45SJed Brown /*@C
1460107e9a97SBarry Smith    ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped
146186994e45SJed Brown 
146286994e45SJed Brown    Not Collective
146386994e45SJed Brown 
146486994e45SJed Brown    Input Arguments:
146586994e45SJed Brown . ltog - local to global mapping
146686994e45SJed Brown 
146786994e45SJed Brown    Output Arguments:
1468565245c5SBarry Smith . array - array of indices, the length of this array may be obtained with ISLocalToGlobalMappingGetSize()
146986994e45SJed Brown 
147086994e45SJed Brown    Level: advanced
147186994e45SJed Brown 
1472107e9a97SBarry Smith    Notes: ISLocalToGlobalMappingGetSize() returns the length the this array
1473107e9a97SBarry Smith 
1474107e9a97SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreIndices(), ISLocalToGlobalMappingGetBlockIndices(), ISLocalToGlobalMappingRestoreBlockIndices()
147586994e45SJed Brown @*/
14767087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
147786994e45SJed Brown {
147886994e45SJed Brown   PetscFunctionBegin;
147986994e45SJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
148086994e45SJed Brown   PetscValidPointer(array,2);
148145b6f7e9SBarry Smith   if (ltog->bs == 1) {
148286994e45SJed Brown     *array = ltog->indices;
148345b6f7e9SBarry Smith   } else {
148445b6f7e9SBarry Smith     PetscInt       *jj,k,i,j,n = ltog->n, bs = ltog->bs;
148545b6f7e9SBarry Smith     const PetscInt *ii;
148645b6f7e9SBarry Smith     PetscErrorCode ierr;
148745b6f7e9SBarry Smith 
148845b6f7e9SBarry Smith     ierr = PetscMalloc1(bs*n,&jj);CHKERRQ(ierr);
148945b6f7e9SBarry Smith     *array = jj;
149045b6f7e9SBarry Smith     k    = 0;
149145b6f7e9SBarry Smith     ii   = ltog->indices;
149245b6f7e9SBarry Smith     for (i=0; i<n; i++)
149345b6f7e9SBarry Smith       for (j=0; j<bs; j++)
149445b6f7e9SBarry Smith         jj[k++] = bs*ii[i] + j;
149545b6f7e9SBarry Smith   }
149686994e45SJed Brown   PetscFunctionReturn(0);
149786994e45SJed Brown }
149886994e45SJed Brown 
149986994e45SJed Brown /*@C
1500193a2b41SJulian Andrej    ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with ISLocalToGlobalMappingGetIndices()
150186994e45SJed Brown 
150286994e45SJed Brown    Not Collective
150386994e45SJed Brown 
150486994e45SJed Brown    Input Arguments:
150586994e45SJed Brown + ltog - local to global mapping
150686994e45SJed Brown - array - array of indices
150786994e45SJed Brown 
150886994e45SJed Brown    Level: advanced
150986994e45SJed Brown 
151086994e45SJed Brown .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices()
151186994e45SJed Brown @*/
15127087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
151386994e45SJed Brown {
151486994e45SJed Brown   PetscFunctionBegin;
151586994e45SJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
151686994e45SJed Brown   PetscValidPointer(array,2);
151745b6f7e9SBarry Smith   if (ltog->bs == 1 && *array != ltog->indices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
151845b6f7e9SBarry Smith 
151945b6f7e9SBarry Smith   if (ltog->bs > 1) {
152045b6f7e9SBarry Smith     PetscErrorCode ierr;
152145b6f7e9SBarry Smith     ierr = PetscFree(*(void**)array);CHKERRQ(ierr);
152245b6f7e9SBarry Smith   }
152345b6f7e9SBarry Smith   PetscFunctionReturn(0);
152445b6f7e9SBarry Smith }
152545b6f7e9SBarry Smith 
152645b6f7e9SBarry Smith /*@C
152745b6f7e9SBarry Smith    ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block
152845b6f7e9SBarry Smith 
152945b6f7e9SBarry Smith    Not Collective
153045b6f7e9SBarry Smith 
153145b6f7e9SBarry Smith    Input Arguments:
153245b6f7e9SBarry Smith . ltog - local to global mapping
153345b6f7e9SBarry Smith 
153445b6f7e9SBarry Smith    Output Arguments:
153545b6f7e9SBarry Smith . array - array of indices
153645b6f7e9SBarry Smith 
153745b6f7e9SBarry Smith    Level: advanced
153845b6f7e9SBarry Smith 
153945b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreBlockIndices()
154045b6f7e9SBarry Smith @*/
154145b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
154245b6f7e9SBarry Smith {
154345b6f7e9SBarry Smith   PetscFunctionBegin;
154445b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
154545b6f7e9SBarry Smith   PetscValidPointer(array,2);
154645b6f7e9SBarry Smith   *array = ltog->indices;
154745b6f7e9SBarry Smith   PetscFunctionReturn(0);
154845b6f7e9SBarry Smith }
154945b6f7e9SBarry Smith 
155045b6f7e9SBarry Smith /*@C
155145b6f7e9SBarry Smith    ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with ISLocalToGlobalMappingGetBlockIndices()
155245b6f7e9SBarry Smith 
155345b6f7e9SBarry Smith    Not Collective
155445b6f7e9SBarry Smith 
155545b6f7e9SBarry Smith    Input Arguments:
155645b6f7e9SBarry Smith + ltog - local to global mapping
155745b6f7e9SBarry Smith - array - array of indices
155845b6f7e9SBarry Smith 
155945b6f7e9SBarry Smith    Level: advanced
156045b6f7e9SBarry Smith 
156145b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices()
156245b6f7e9SBarry Smith @*/
156345b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
156445b6f7e9SBarry Smith {
156545b6f7e9SBarry Smith   PetscFunctionBegin;
156645b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
156745b6f7e9SBarry Smith   PetscValidPointer(array,2);
156886994e45SJed Brown   if (*array != ltog->indices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
15690298fd71SBarry Smith   *array = NULL;
157086994e45SJed Brown   PetscFunctionReturn(0);
157186994e45SJed Brown }
1572f7efa3c7SJed Brown 
1573f7efa3c7SJed Brown /*@C
1574f7efa3c7SJed Brown    ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings
1575f7efa3c7SJed Brown 
1576f7efa3c7SJed Brown    Not Collective
1577f7efa3c7SJed Brown 
1578f7efa3c7SJed Brown    Input Arguments:
1579f7efa3c7SJed Brown + comm - communicator for the new mapping, must contain the communicator of every mapping to concatenate
1580f7efa3c7SJed Brown . n - number of mappings to concatenate
1581f7efa3c7SJed Brown - ltogs - local to global mappings
1582f7efa3c7SJed Brown 
1583f7efa3c7SJed Brown    Output Arguments:
1584f7efa3c7SJed Brown . ltogcat - new mapping
1585f7efa3c7SJed Brown 
15869d90f715SBarry Smith    Note: this currently always returns a mapping with block size of 1
15879d90f715SBarry Smith 
15889d90f715SBarry Smith    Developer Note: If all the input mapping have the same block size we could easily handle that as a special case
15899d90f715SBarry Smith 
1590f7efa3c7SJed Brown    Level: advanced
1591f7efa3c7SJed Brown 
1592f7efa3c7SJed Brown .seealso: ISLocalToGlobalMappingCreate()
1593f7efa3c7SJed Brown @*/
1594f7efa3c7SJed Brown PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm,PetscInt n,const ISLocalToGlobalMapping ltogs[],ISLocalToGlobalMapping *ltogcat)
1595f7efa3c7SJed Brown {
1596f7efa3c7SJed Brown   PetscInt       i,cnt,m,*idx;
1597f7efa3c7SJed Brown   PetscErrorCode ierr;
1598f7efa3c7SJed Brown 
1599f7efa3c7SJed Brown   PetscFunctionBegin;
1600f7efa3c7SJed Brown   if (n < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have a non-negative number of mappings, given %D",n);
1601f7efa3c7SJed Brown   if (n > 0) PetscValidPointer(ltogs,3);
1602f7efa3c7SJed Brown   for (i=0; i<n; i++) PetscValidHeaderSpecific(ltogs[i],IS_LTOGM_CLASSID,3);
1603f7efa3c7SJed Brown   PetscValidPointer(ltogcat,4);
1604f7efa3c7SJed Brown   for (cnt=0,i=0; i<n; i++) {
1605f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
1606f7efa3c7SJed Brown     cnt += m;
1607f7efa3c7SJed Brown   }
1608785e854fSJed Brown   ierr = PetscMalloc1(cnt,&idx);CHKERRQ(ierr);
1609f7efa3c7SJed Brown   for (cnt=0,i=0; i<n; i++) {
1610f7efa3c7SJed Brown     const PetscInt *subidx;
1611f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
1612f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetIndices(ltogs[i],&subidx);CHKERRQ(ierr);
1613f7efa3c7SJed Brown     ierr = PetscMemcpy(&idx[cnt],subidx,m*sizeof(PetscInt));CHKERRQ(ierr);
1614f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingRestoreIndices(ltogs[i],&subidx);CHKERRQ(ierr);
1615f7efa3c7SJed Brown     cnt += m;
1616f7efa3c7SJed Brown   }
1617f0413b6fSBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,1,cnt,idx,PETSC_OWN_POINTER,ltogcat);CHKERRQ(ierr);
1618f7efa3c7SJed Brown   PetscFunctionReturn(0);
1619f7efa3c7SJed Brown }
162004a59952SBarry Smith 
1621413f72f0SBarry Smith /*MC
1622413f72f0SBarry Smith       ISLOCALTOGLOBALMAPPINGBASIC - basic implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is
1623413f72f0SBarry Smith                                     used this is good for only small and moderate size problems.
1624413f72f0SBarry Smith 
1625413f72f0SBarry Smith    Options Database Keys:
1626413f72f0SBarry Smith +   -islocaltoglobalmapping_type basic - select this method
1627413f72f0SBarry Smith 
1628413f72f0SBarry Smith    Level: beginner
1629413f72f0SBarry Smith 
1630413f72f0SBarry Smith .seealso:  ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType(), ISLOCALTOGLOBALMAPPINGHASH
1631413f72f0SBarry Smith M*/
1632413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Basic(ISLocalToGlobalMapping ltog)
1633413f72f0SBarry Smith {
1634413f72f0SBarry Smith   PetscFunctionBegin;
1635413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Basic;
1636413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Basic;
1637413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Basic;
1638413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Basic;
1639413f72f0SBarry Smith   PetscFunctionReturn(0);
1640413f72f0SBarry Smith }
1641413f72f0SBarry Smith 
1642413f72f0SBarry Smith /*MC
1643413f72f0SBarry Smith       ISLOCALTOGLOBALMAPPINGHASH - hash implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is
1644ed56e8eaSBarry Smith                                     used this is good for large memory problems.
1645413f72f0SBarry Smith 
1646413f72f0SBarry Smith    Options Database Keys:
1647413f72f0SBarry Smith +   -islocaltoglobalmapping_type hash - select this method
1648413f72f0SBarry Smith 
1649ed56e8eaSBarry Smith 
1650ed56e8eaSBarry Smith    Notes: This is selected automatically for large problems if the user does not set the type.
1651ed56e8eaSBarry Smith 
1652413f72f0SBarry Smith    Level: beginner
1653413f72f0SBarry Smith 
1654413f72f0SBarry Smith .seealso:  ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType(), ISLOCALTOGLOBALMAPPINGHASH
1655413f72f0SBarry Smith M*/
1656413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Hash(ISLocalToGlobalMapping ltog)
1657413f72f0SBarry Smith {
1658413f72f0SBarry Smith   PetscFunctionBegin;
1659413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Hash;
1660413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Hash;
1661413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Hash;
1662413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Hash;
1663413f72f0SBarry Smith   PetscFunctionReturn(0);
1664413f72f0SBarry Smith }
1665413f72f0SBarry Smith 
1666413f72f0SBarry Smith 
1667413f72f0SBarry Smith /*@C
1668413f72f0SBarry Smith     ISLocalToGlobalMappingRegister -  Adds a method for applying a global to local mapping with an ISLocalToGlobalMapping
1669413f72f0SBarry Smith 
1670413f72f0SBarry Smith    Not Collective
1671413f72f0SBarry Smith 
1672413f72f0SBarry Smith    Input Parameters:
1673413f72f0SBarry Smith +  sname - name of a new method
1674413f72f0SBarry Smith -  routine_create - routine to create method context
1675413f72f0SBarry Smith 
1676413f72f0SBarry Smith    Notes:
1677ed56e8eaSBarry Smith    ISLocalToGlobalMappingRegister() may be called multiple times to add several user-defined mappings.
1678413f72f0SBarry Smith 
1679413f72f0SBarry Smith    Sample usage:
1680413f72f0SBarry Smith .vb
1681413f72f0SBarry Smith    ISLocalToGlobalMappingRegister("my_mapper",MyCreate);
1682413f72f0SBarry Smith .ve
1683413f72f0SBarry Smith 
1684ed56e8eaSBarry Smith    Then, your mapping can be chosen with the procedural interface via
1685413f72f0SBarry Smith $     ISLocalToGlobalMappingSetType(ltog,"my_mapper")
1686413f72f0SBarry Smith    or at runtime via the option
1687ed56e8eaSBarry Smith $     -islocaltoglobalmapping_type my_mapper
1688413f72f0SBarry Smith 
1689413f72f0SBarry Smith    Level: advanced
1690413f72f0SBarry Smith 
1691413f72f0SBarry Smith .keywords: ISLocalToGlobalMappingType, register
1692413f72f0SBarry Smith 
1693413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingRegisterAll(), ISLocalToGlobalMappingRegisterDestroy(), ISLOCALTOGLOBALMAPPINGBASIC, ISLOCALTOGLOBALMAPPINGHASH
1694413f72f0SBarry Smith 
1695413f72f0SBarry Smith @*/
1696413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRegister(const char sname[],PetscErrorCode (*function)(ISLocalToGlobalMapping))
1697413f72f0SBarry Smith {
1698413f72f0SBarry Smith   PetscErrorCode ierr;
1699413f72f0SBarry Smith 
1700413f72f0SBarry Smith   PetscFunctionBegin;
1701413f72f0SBarry Smith   ierr = PetscFunctionListAdd(&ISLocalToGlobalMappingList,sname,function);CHKERRQ(ierr);
1702413f72f0SBarry Smith   PetscFunctionReturn(0);
1703413f72f0SBarry Smith }
1704413f72f0SBarry Smith 
1705413f72f0SBarry Smith /*@C
1706ed56e8eaSBarry Smith    ISLocalToGlobalMappingSetType - Builds ISLocalToGlobalMapping for a particular global to local mapping approach.
1707413f72f0SBarry Smith 
1708413f72f0SBarry Smith    Logically Collective on ISLocalToGlobalMapping
1709413f72f0SBarry Smith 
1710413f72f0SBarry Smith    Input Parameters:
1711413f72f0SBarry Smith +  ltog - the ISLocalToGlobalMapping object
1712413f72f0SBarry Smith -  type - a known method
1713413f72f0SBarry Smith 
1714413f72f0SBarry Smith    Options Database Key:
1715ed56e8eaSBarry Smith .  -islocaltoglobalmapping_type  <method> - Sets the method; use -help for a list
1716413f72f0SBarry Smith     of available methods (for instance, basic or hash)
1717413f72f0SBarry Smith 
1718413f72f0SBarry Smith    Notes:
1719413f72f0SBarry Smith    See "petsc/include/petscis.h" for available methods
1720413f72f0SBarry Smith 
1721413f72f0SBarry Smith   Normally, it is best to use the ISLocalToGlobalMappingSetFromOptions() command and
1722413f72f0SBarry Smith   then set the ISLocalToGlobalMapping type from the options database rather than by using
1723413f72f0SBarry Smith   this routine.
1724413f72f0SBarry Smith 
1725413f72f0SBarry Smith   Level: intermediate
1726413f72f0SBarry Smith 
1727413f72f0SBarry Smith   Developer Note: ISLocalToGlobalMappingRegister() is used to add new types to ISLocalToGlobalMappingList from which they
1728413f72f0SBarry Smith   are accessed by ISLocalToGlobalMappingSetType().
1729413f72f0SBarry Smith 
1730413f72f0SBarry Smith .keywords: ISLocalToGlobalMapping, set, method
1731413f72f0SBarry Smith 
1732413f72f0SBarry Smith .seealso: PCSetType(), ISLocalToGlobalMappingType, ISLocalToGlobalMappingRegister(), ISLocalToGlobalMappingCreate()
1733413f72f0SBarry Smith 
1734413f72f0SBarry Smith @*/
1735413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType type)
1736413f72f0SBarry Smith {
1737413f72f0SBarry Smith   PetscErrorCode ierr,(*r)(ISLocalToGlobalMapping);
1738413f72f0SBarry Smith   PetscBool      match;
1739413f72f0SBarry Smith 
1740413f72f0SBarry Smith   PetscFunctionBegin;
1741413f72f0SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
1742413f72f0SBarry Smith   PetscValidCharPointer(type,2);
1743413f72f0SBarry Smith 
1744413f72f0SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)ltog,type,&match);CHKERRQ(ierr);
1745413f72f0SBarry Smith   if (match) PetscFunctionReturn(0);
1746413f72f0SBarry Smith 
1747413f72f0SBarry Smith   ierr =  PetscFunctionListFind(ISLocalToGlobalMappingList,type,&r);CHKERRQ(ierr);
1748413f72f0SBarry Smith   if (!r) SETERRQ1(PetscObjectComm((PetscObject)ltog),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested ISLocalToGlobalMapping type %s",type);
1749413f72f0SBarry Smith   /* Destroy the previous private LTOG context */
1750413f72f0SBarry Smith   if (ltog->ops->destroy) {
1751413f72f0SBarry Smith     ierr              = (*ltog->ops->destroy)(ltog);CHKERRQ(ierr);
1752413f72f0SBarry Smith     ltog->ops->destroy = NULL;
1753413f72f0SBarry Smith   }
1754413f72f0SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)ltog,type);CHKERRQ(ierr);
1755413f72f0SBarry Smith   ierr = (*r)(ltog);CHKERRQ(ierr);
1756413f72f0SBarry Smith   PetscFunctionReturn(0);
1757413f72f0SBarry Smith }
1758413f72f0SBarry Smith 
1759413f72f0SBarry Smith PetscBool ISLocalToGlobalMappingRegisterAllCalled = PETSC_FALSE;
1760413f72f0SBarry Smith 
1761413f72f0SBarry Smith /*@C
1762413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAll - Registers all of the local to global mapping components in the IS package.
1763413f72f0SBarry Smith 
1764413f72f0SBarry Smith   Not Collective
1765413f72f0SBarry Smith 
1766413f72f0SBarry Smith   Level: advanced
1767413f72f0SBarry Smith 
1768413f72f0SBarry Smith .keywords: IS, register, all
1769413f72f0SBarry Smith .seealso:  ISRegister(),  ISLocalToGlobalRegister()
1770413f72f0SBarry Smith @*/
1771413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRegisterAll(void)
1772413f72f0SBarry Smith {
1773413f72f0SBarry Smith   PetscErrorCode ierr;
1774413f72f0SBarry Smith 
1775413f72f0SBarry Smith   PetscFunctionBegin;
1776413f72f0SBarry Smith   if (ISLocalToGlobalMappingRegisterAllCalled) PetscFunctionReturn(0);
1777413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAllCalled = PETSC_TRUE;
1778413f72f0SBarry Smith 
1779413f72f0SBarry Smith   ierr = ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGBASIC, ISLocalToGlobalMappingCreate_Basic);CHKERRQ(ierr);
1780413f72f0SBarry Smith   ierr = ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGHASH, ISLocalToGlobalMappingCreate_Hash);CHKERRQ(ierr);
1781413f72f0SBarry Smith   PetscFunctionReturn(0);
1782413f72f0SBarry Smith }
178304a59952SBarry Smith 
1784