xref: /petsc/src/vec/is/utils/isltog.c (revision 95452b02e12c0ee11232c7ff2b24b568a8e07e43)
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 
262*95452b02SPatrick Sanan     Notes:
263*95452b02SPatrick Sanan     the block size of the IS determines the block size of the mapping
264a997ad1aSLois Curfman McInnes     Level: advanced
265a997ad1aSLois Curfman McInnes 
266273d9f13SBarry Smith     Concepts: mapping^local to global
2672bdab257SBarry Smith 
2687e99dc12SLawrence Mitchell .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetFromOptions()
2692bdab257SBarry Smith @*/
2707087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingCreateIS(IS is,ISLocalToGlobalMapping *mapping)
2712bdab257SBarry Smith {
2726849ba73SBarry Smith   PetscErrorCode ierr;
2733bbf0e92SBarry Smith   PetscInt       n,bs;
2745d0c19d7SBarry Smith   const PetscInt *indices;
2752bdab257SBarry Smith   MPI_Comm       comm;
2763bbf0e92SBarry Smith   PetscBool      isblock;
2773a40ed3dSBarry Smith 
2783a40ed3dSBarry Smith   PetscFunctionBegin;
2790700a824SBarry Smith   PetscValidHeaderSpecific(is,IS_CLASSID,1);
2804482741eSBarry Smith   PetscValidPointer(mapping,2);
2812bdab257SBarry Smith 
2822bdab257SBarry Smith   ierr = PetscObjectGetComm((PetscObject)is,&comm);CHKERRQ(ierr);
2833b9aefa3SBarry Smith   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
2843bbf0e92SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)is,ISBLOCK,&isblock);CHKERRQ(ierr);
2856006e8d2SBarry Smith   if (!isblock) {
286f0413b6fSBarry Smith     ierr = ISGetIndices(is,&indices);CHKERRQ(ierr);
287f0413b6fSBarry Smith     ierr = ISLocalToGlobalMappingCreate(comm,1,n,indices,PETSC_COPY_VALUES,mapping);CHKERRQ(ierr);
2882bdab257SBarry Smith     ierr = ISRestoreIndices(is,&indices);CHKERRQ(ierr);
2896006e8d2SBarry Smith   } else {
2906006e8d2SBarry Smith     ierr = ISGetBlockSize(is,&bs);CHKERRQ(ierr);
291f0413b6fSBarry Smith     ierr = ISBlockGetIndices(is,&indices);CHKERRQ(ierr);
29228bc9809SBarry Smith     ierr = ISLocalToGlobalMappingCreate(comm,bs,n/bs,indices,PETSC_COPY_VALUES,mapping);CHKERRQ(ierr);
293f0413b6fSBarry Smith     ierr = ISBlockRestoreIndices(is,&indices);CHKERRQ(ierr);
2946006e8d2SBarry Smith   }
2953a40ed3dSBarry Smith   PetscFunctionReturn(0);
2962bdab257SBarry Smith }
2975a5d4f66SBarry Smith 
298a4d96a55SJed Brown /*@C
299a4d96a55SJed Brown     ISLocalToGlobalMappingCreateSF - Creates a mapping between a local (0 to n)
300a4d96a55SJed Brown     ordering and a global parallel ordering.
301a4d96a55SJed Brown 
302a4d96a55SJed Brown     Collective
303a4d96a55SJed Brown 
304a4d96a55SJed Brown     Input Parameter:
305a4d96a55SJed Brown +   sf - star forest mapping contiguous local indices to (rank, offset)
306a4d96a55SJed Brown -   start - first global index on this process
307a4d96a55SJed Brown 
308a4d96a55SJed Brown     Output Parameter:
309a4d96a55SJed Brown .   mapping - new mapping data structure
310a4d96a55SJed Brown 
311a4d96a55SJed Brown     Level: advanced
312a4d96a55SJed Brown 
313a4d96a55SJed Brown     Concepts: mapping^local to global
314a4d96a55SJed Brown 
3157e99dc12SLawrence Mitchell .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingSetFromOptions()
316a4d96a55SJed Brown @*/
317a4d96a55SJed Brown PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF sf,PetscInt start,ISLocalToGlobalMapping *mapping)
318a4d96a55SJed Brown {
319a4d96a55SJed Brown   PetscErrorCode ierr;
320a4d96a55SJed Brown   PetscInt       i,maxlocal,nroots,nleaves,*globals,*ltog;
321a4d96a55SJed Brown   const PetscInt *ilocal;
322a4d96a55SJed Brown   MPI_Comm       comm;
323a4d96a55SJed Brown 
324a4d96a55SJed Brown   PetscFunctionBegin;
325a4d96a55SJed Brown   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
326a4d96a55SJed Brown   PetscValidPointer(mapping,3);
327a4d96a55SJed Brown 
328a4d96a55SJed Brown   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
3290298fd71SBarry Smith   ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr);
330f6e5521dSKarl Rupp   if (ilocal) {
331f6e5521dSKarl Rupp     for (i=0,maxlocal=0; i<nleaves; i++) maxlocal = PetscMax(maxlocal,ilocal[i]+1);
332f6e5521dSKarl Rupp   }
333a4d96a55SJed Brown   else maxlocal = nleaves;
334785e854fSJed Brown   ierr = PetscMalloc1(nroots,&globals);CHKERRQ(ierr);
335785e854fSJed Brown   ierr = PetscMalloc1(maxlocal,&ltog);CHKERRQ(ierr);
336a4d96a55SJed Brown   for (i=0; i<nroots; i++) globals[i] = start + i;
337a4d96a55SJed Brown   for (i=0; i<maxlocal; i++) ltog[i] = -1;
338a4d96a55SJed Brown   ierr = PetscSFBcastBegin(sf,MPIU_INT,globals,ltog);CHKERRQ(ierr);
339a4d96a55SJed Brown   ierr = PetscSFBcastEnd(sf,MPIU_INT,globals,ltog);CHKERRQ(ierr);
340f0413b6fSBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,1,maxlocal,ltog,PETSC_OWN_POINTER,mapping);CHKERRQ(ierr);
341a4d96a55SJed Brown   ierr = PetscFree(globals);CHKERRQ(ierr);
342a4d96a55SJed Brown   PetscFunctionReturn(0);
343a4d96a55SJed Brown }
344b46b645bSBarry Smith 
34563fa5c83Sstefano_zampini /*@
34663fa5c83Sstefano_zampini     ISLocalToGlobalMappingSetBlockSize - Sets the blocksize of the mapping
34763fa5c83Sstefano_zampini 
34863fa5c83Sstefano_zampini     Not collective
34963fa5c83Sstefano_zampini 
35063fa5c83Sstefano_zampini     Input Parameters:
35163fa5c83Sstefano_zampini .   mapping - mapping data structure
35263fa5c83Sstefano_zampini .   bs - the blocksize
35363fa5c83Sstefano_zampini 
35463fa5c83Sstefano_zampini     Level: advanced
35563fa5c83Sstefano_zampini 
35663fa5c83Sstefano_zampini     Concepts: mapping^local to global
35763fa5c83Sstefano_zampini 
35863fa5c83Sstefano_zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS()
35963fa5c83Sstefano_zampini @*/
36063fa5c83Sstefano_zampini PetscErrorCode  ISLocalToGlobalMappingSetBlockSize(ISLocalToGlobalMapping mapping,PetscInt bs)
36163fa5c83Sstefano_zampini {
362a59f3c4dSstefano_zampini   PetscInt       *nid;
363a59f3c4dSstefano_zampini   const PetscInt *oid;
364a59f3c4dSstefano_zampini   PetscInt       i,cn,on,obs,nn;
36563fa5c83Sstefano_zampini   PetscErrorCode ierr;
36663fa5c83Sstefano_zampini 
36763fa5c83Sstefano_zampini   PetscFunctionBegin;
36863fa5c83Sstefano_zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
36963fa5c83Sstefano_zampini   if (bs < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid block size %D",bs);
37063fa5c83Sstefano_zampini   if (bs == mapping->bs) PetscFunctionReturn(0);
37163fa5c83Sstefano_zampini   on  = mapping->n;
37263fa5c83Sstefano_zampini   obs = mapping->bs;
37363fa5c83Sstefano_zampini   oid = mapping->indices;
37463fa5c83Sstefano_zampini   nn  = (on*obs)/bs;
37563fa5c83Sstefano_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);
376a59f3c4dSstefano_zampini 
37763fa5c83Sstefano_zampini   ierr = PetscMalloc1(nn,&nid);CHKERRQ(ierr);
378a59f3c4dSstefano_zampini   ierr = ISLocalToGlobalMappingGetIndices(mapping,&oid);CHKERRQ(ierr);
379a59f3c4dSstefano_zampini   for (i=0;i<nn;i++) {
380a59f3c4dSstefano_zampini     PetscInt j;
381a59f3c4dSstefano_zampini     for (j=0,cn=0;j<bs-1;j++) {
382a59f3c4dSstefano_zampini       if (oid[i*bs+j] < 0) { cn++; continue; }
383a59f3c4dSstefano_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]);
384a59f3c4dSstefano_zampini     }
385a59f3c4dSstefano_zampini     if (oid[i*bs+j] < 0) cn++;
3868b7cb0e6Sstefano_zampini     if (cn) {
387a59f3c4dSstefano_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);
388a59f3c4dSstefano_zampini       nid[i] = -1;
3898b7cb0e6Sstefano_zampini     } else {
390a59f3c4dSstefano_zampini       nid[i] = oid[i*bs]/bs;
39163fa5c83Sstefano_zampini     }
39263fa5c83Sstefano_zampini   }
393a59f3c4dSstefano_zampini   ierr = ISLocalToGlobalMappingRestoreIndices(mapping,&oid);CHKERRQ(ierr);
394a59f3c4dSstefano_zampini 
39563fa5c83Sstefano_zampini   mapping->n           = nn;
39663fa5c83Sstefano_zampini   mapping->bs          = bs;
39763fa5c83Sstefano_zampini   ierr                 = PetscFree(mapping->indices);CHKERRQ(ierr);
39863fa5c83Sstefano_zampini   mapping->indices     = nid;
399c9345713Sstefano_zampini   mapping->globalstart = 0;
400c9345713Sstefano_zampini   mapping->globalend   = 0;
401413f72f0SBarry Smith   if (mapping->ops->destroy) {
402413f72f0SBarry Smith     ierr = (*mapping->ops->destroy)(mapping);CHKERRQ(ierr);
403413f72f0SBarry Smith   }
40463fa5c83Sstefano_zampini   PetscFunctionReturn(0);
40563fa5c83Sstefano_zampini }
40663fa5c83Sstefano_zampini 
40745b6f7e9SBarry Smith /*@
40845b6f7e9SBarry Smith     ISLocalToGlobalMappingGetBlockSize - Gets the blocksize of the mapping
40945b6f7e9SBarry Smith     ordering and a global parallel ordering.
41045b6f7e9SBarry Smith 
41145b6f7e9SBarry Smith     Not Collective
41245b6f7e9SBarry Smith 
41345b6f7e9SBarry Smith     Input Parameters:
41445b6f7e9SBarry Smith .   mapping - mapping data structure
41545b6f7e9SBarry Smith 
41645b6f7e9SBarry Smith     Output Parameter:
41745b6f7e9SBarry Smith .   bs - the blocksize
41845b6f7e9SBarry Smith 
41945b6f7e9SBarry Smith     Level: advanced
42045b6f7e9SBarry Smith 
42145b6f7e9SBarry Smith     Concepts: mapping^local to global
42245b6f7e9SBarry Smith 
42345b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS()
42445b6f7e9SBarry Smith @*/
42545b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping mapping,PetscInt *bs)
42645b6f7e9SBarry Smith {
42745b6f7e9SBarry Smith   PetscFunctionBegin;
428cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
42945b6f7e9SBarry Smith   *bs = mapping->bs;
43045b6f7e9SBarry Smith   PetscFunctionReturn(0);
43145b6f7e9SBarry Smith }
43245b6f7e9SBarry Smith 
433ba5bb76aSSatish Balay /*@
43490f02eecSBarry Smith     ISLocalToGlobalMappingCreate - Creates a mapping between a local (0 to n)
43590f02eecSBarry Smith     ordering and a global parallel ordering.
4362362add9SBarry Smith 
43789d82c54SBarry Smith     Not Collective, but communicator may have more than one process
438b9cd556bSLois Curfman McInnes 
4392362add9SBarry Smith     Input Parameters:
44089d82c54SBarry Smith +   comm - MPI communicator
441f0413b6fSBarry Smith .   bs - the block size
44228bc9809SBarry Smith .   n - the number of local elements divided by the block size, or equivalently the number of block indices
44328bc9809SBarry 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
444d5ad8652SBarry Smith -   mode - see PetscCopyMode
4452362add9SBarry Smith 
446a997ad1aSLois Curfman McInnes     Output Parameter:
44790f02eecSBarry Smith .   mapping - new mapping data structure
4482362add9SBarry Smith 
449*95452b02SPatrick Sanan     Notes:
450*95452b02SPatrick Sanan     There is one integer value in indices per block and it represents the actual indices bs*idx + j, where j=0,..,bs-1
451413f72f0SBarry Smith 
4529a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
453413f72f0SBarry 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.
454413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
455413f72f0SBarry Smith 
456a997ad1aSLois Curfman McInnes     Level: advanced
457a997ad1aSLois Curfman McInnes 
458273d9f13SBarry Smith     Concepts: mapping^local to global
4592362add9SBarry Smith 
460413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingSetFromOptions(), ISLOCALTOGLOBALMAPPINGBASIC, ISLOCALTOGLOBALMAPPINGHASH
461413f72f0SBarry Smith           ISLocalToGlobalMappingSetType(), ISLocalToGlobalMappingType
4622362add9SBarry Smith @*/
46360c7cefcSBarry Smith PetscErrorCode  ISLocalToGlobalMappingCreate(MPI_Comm comm,PetscInt bs,PetscInt n,const PetscInt indices[],PetscCopyMode mode,ISLocalToGlobalMapping *mapping)
4642362add9SBarry Smith {
4656849ba73SBarry Smith   PetscErrorCode ierr;
46632dcc486SBarry Smith   PetscInt       *in;
467b46b645bSBarry Smith 
468b46b645bSBarry Smith   PetscFunctionBegin;
46973911063SBarry Smith   if (n) PetscValidIntPointer(indices,3);
4704482741eSBarry Smith   PetscValidPointer(mapping,4);
471b46b645bSBarry Smith 
4720298fd71SBarry Smith   *mapping = NULL;
473607a6623SBarry Smith   ierr = ISInitializePackage();CHKERRQ(ierr);
4742362add9SBarry Smith 
47573107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(*mapping,IS_LTOGM_CLASSID,"ISLocalToGlobalMapping","Local to global mapping","IS",
47660c7cefcSBarry Smith                            comm,ISLocalToGlobalMappingDestroy,ISLocalToGlobalMappingView);CHKERRQ(ierr);
477d4bb536fSBarry Smith   (*mapping)->n             = n;
478f0413b6fSBarry Smith   (*mapping)->bs            = bs;
479268a049cSStefano Zampini   (*mapping)->info_cached   = PETSC_FALSE;
480268a049cSStefano Zampini   (*mapping)->info_free     = PETSC_FALSE;
481268a049cSStefano Zampini   (*mapping)->info_procs    = NULL;
482268a049cSStefano Zampini   (*mapping)->info_numprocs = NULL;
483268a049cSStefano Zampini   (*mapping)->info_indices  = NULL;
484413f72f0SBarry Smith 
485413f72f0SBarry Smith   (*mapping)->ops->globaltolocalmappingapply      = NULL;
486413f72f0SBarry Smith   (*mapping)->ops->globaltolocalmappingapplyblock = NULL;
487413f72f0SBarry Smith   (*mapping)->ops->destroy                        = NULL;
488d5ad8652SBarry Smith   if (mode == PETSC_COPY_VALUES) {
489785e854fSJed Brown     ierr = PetscMalloc1(n,&in);CHKERRQ(ierr);
490d5ad8652SBarry Smith     ierr = PetscMemcpy(in,indices,n*sizeof(PetscInt));CHKERRQ(ierr);
491d5ad8652SBarry Smith     (*mapping)->indices = in;
4926389a1a1SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt));CHKERRQ(ierr);
4936389a1a1SBarry Smith   } else if (mode == PETSC_OWN_POINTER) {
4946389a1a1SBarry Smith     (*mapping)->indices = (PetscInt*)indices;
4956389a1a1SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt));CHKERRQ(ierr);
4966389a1a1SBarry Smith   }
49760c7cefcSBarry Smith   else SETERRQ(comm,PETSC_ERR_SUP,"Cannot currently use PETSC_USE_POINTER");
4983a40ed3dSBarry Smith   PetscFunctionReturn(0);
4992362add9SBarry Smith }
5002362add9SBarry Smith 
501413f72f0SBarry Smith PetscFunctionList ISLocalToGlobalMappingList = NULL;
502413f72f0SBarry Smith 
503413f72f0SBarry Smith 
50490f02eecSBarry Smith /*@
5057e99dc12SLawrence Mitchell    ISLocalToGlobalMappingSetFromOptions - Set mapping options from the options database.
5067e99dc12SLawrence Mitchell 
5077e99dc12SLawrence Mitchell    Not collective
5087e99dc12SLawrence Mitchell 
5097e99dc12SLawrence Mitchell    Input Parameters:
5107e99dc12SLawrence Mitchell .  mapping - mapping data structure
5117e99dc12SLawrence Mitchell 
5127e99dc12SLawrence Mitchell    Level: advanced
5137e99dc12SLawrence Mitchell 
5147e99dc12SLawrence Mitchell @*/
5157e99dc12SLawrence Mitchell PetscErrorCode ISLocalToGlobalMappingSetFromOptions(ISLocalToGlobalMapping mapping)
5167e99dc12SLawrence Mitchell {
5177e99dc12SLawrence Mitchell   PetscErrorCode             ierr;
518413f72f0SBarry Smith   char                       type[256];
519413f72f0SBarry Smith   ISLocalToGlobalMappingType defaulttype = "Not set";
5207e99dc12SLawrence Mitchell   PetscBool                  flg;
5217e99dc12SLawrence Mitchell 
5227e99dc12SLawrence Mitchell   PetscFunctionBegin;
5237e99dc12SLawrence Mitchell   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
524413f72f0SBarry Smith   ierr = ISLocalToGlobalMappingRegisterAll();CHKERRQ(ierr);
5257e99dc12SLawrence Mitchell   ierr = PetscObjectOptionsBegin((PetscObject)mapping);CHKERRQ(ierr);
526413f72f0SBarry Smith   ierr = PetscOptionsFList("-islocaltoglobalmapping_type","ISLocalToGlobalMapping method","ISLocalToGlobalMappingSetType",ISLocalToGlobalMappingList,(char*)(((PetscObject)mapping)->type_name) ? ((PetscObject)mapping)->type_name : defaulttype,type,256,&flg);CHKERRQ(ierr);
527413f72f0SBarry Smith   if (flg) {
528413f72f0SBarry Smith     ierr = ISLocalToGlobalMappingSetType(mapping,type);CHKERRQ(ierr);
529413f72f0SBarry Smith   }
5307e99dc12SLawrence Mitchell   ierr = PetscOptionsEnd();CHKERRQ(ierr);
5317e99dc12SLawrence Mitchell   PetscFunctionReturn(0);
5327e99dc12SLawrence Mitchell }
5337e99dc12SLawrence Mitchell 
5347e99dc12SLawrence Mitchell /*@
53590f02eecSBarry Smith    ISLocalToGlobalMappingDestroy - Destroys a mapping between a local (0 to n)
53690f02eecSBarry Smith    ordering and a global parallel ordering.
53790f02eecSBarry Smith 
5380f5bd95cSBarry Smith    Note Collective
539b9cd556bSLois Curfman McInnes 
54090f02eecSBarry Smith    Input Parameters:
54190f02eecSBarry Smith .  mapping - mapping data structure
54290f02eecSBarry Smith 
543a997ad1aSLois Curfman McInnes    Level: advanced
544a997ad1aSLois Curfman McInnes 
5453acfe500SLois Curfman McInnes .seealso: ISLocalToGlobalMappingCreate()
54690f02eecSBarry Smith @*/
5476bf464f9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping *mapping)
54890f02eecSBarry Smith {
549dfbe8321SBarry Smith   PetscErrorCode ierr;
5505fd66863SKarl Rupp 
5513a40ed3dSBarry Smith   PetscFunctionBegin;
5526bf464f9SBarry Smith   if (!*mapping) PetscFunctionReturn(0);
5536bf464f9SBarry Smith   PetscValidHeaderSpecific((*mapping),IS_LTOGM_CLASSID,1);
554997056adSBarry Smith   if (--((PetscObject)(*mapping))->refct > 0) {*mapping = 0;PetscFunctionReturn(0);}
5556bf464f9SBarry Smith   ierr = PetscFree((*mapping)->indices);CHKERRQ(ierr);
556268a049cSStefano Zampini   ierr = PetscFree((*mapping)->info_procs);CHKERRQ(ierr);
557268a049cSStefano Zampini   ierr = PetscFree((*mapping)->info_numprocs);CHKERRQ(ierr);
558268a049cSStefano Zampini   if ((*mapping)->info_indices) {
559268a049cSStefano Zampini     PetscInt i;
560268a049cSStefano Zampini 
561268a049cSStefano Zampini     ierr = PetscFree(((*mapping)->info_indices)[0]);CHKERRQ(ierr);
562268a049cSStefano Zampini     for (i=1; i<(*mapping)->info_nproc; i++) {
563268a049cSStefano Zampini       ierr = PetscFree(((*mapping)->info_indices)[i]);CHKERRQ(ierr);
564268a049cSStefano Zampini     }
565268a049cSStefano Zampini     ierr = PetscFree((*mapping)->info_indices);CHKERRQ(ierr);
566268a049cSStefano Zampini   }
567413f72f0SBarry Smith   if ((*mapping)->ops->destroy) {
568413f72f0SBarry Smith     ierr = (*(*mapping)->ops->destroy)(*mapping);CHKERRQ(ierr);
569413f72f0SBarry Smith   }
570d38fa0fbSBarry Smith   ierr     = PetscHeaderDestroy(mapping);CHKERRQ(ierr);
571992144d0SBarry Smith   *mapping = 0;
5723a40ed3dSBarry Smith   PetscFunctionReturn(0);
57390f02eecSBarry Smith }
57490f02eecSBarry Smith 
57590f02eecSBarry Smith /*@
5763acfe500SLois Curfman McInnes     ISLocalToGlobalMappingApplyIS - Creates from an IS in the local numbering
5773acfe500SLois Curfman McInnes     a new index set using the global numbering defined in an ISLocalToGlobalMapping
5783acfe500SLois Curfman McInnes     context.
57990f02eecSBarry Smith 
580b9cd556bSLois Curfman McInnes     Not collective
581b9cd556bSLois Curfman McInnes 
58290f02eecSBarry Smith     Input Parameters:
583b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
584b9cd556bSLois Curfman McInnes -   is - index set in local numbering
58590f02eecSBarry Smith 
58690f02eecSBarry Smith     Output Parameters:
58790f02eecSBarry Smith .   newis - index set in global numbering
58890f02eecSBarry Smith 
589a997ad1aSLois Curfman McInnes     Level: advanced
590a997ad1aSLois Curfman McInnes 
591273d9f13SBarry Smith     Concepts: mapping^local to global
5923acfe500SLois Curfman McInnes 
59390f02eecSBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),
594d4bb536fSBarry Smith           ISLocalToGlobalMappingDestroy(), ISGlobalToLocalMappingApply()
59590f02eecSBarry Smith @*/
5967087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping,IS is,IS *newis)
59790f02eecSBarry Smith {
5986849ba73SBarry Smith   PetscErrorCode ierr;
599e24637baSBarry Smith   PetscInt       n,*idxout;
6005d0c19d7SBarry Smith   const PetscInt *idxin;
6013a40ed3dSBarry Smith 
6023a40ed3dSBarry Smith   PetscFunctionBegin;
6030700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
6040700a824SBarry Smith   PetscValidHeaderSpecific(is,IS_CLASSID,2);
6054482741eSBarry Smith   PetscValidPointer(newis,3);
60690f02eecSBarry Smith 
6073b9aefa3SBarry Smith   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
60890f02eecSBarry Smith   ierr = ISGetIndices(is,&idxin);CHKERRQ(ierr);
609785e854fSJed Brown   ierr = PetscMalloc1(n,&idxout);CHKERRQ(ierr);
610e24637baSBarry Smith   ierr = ISLocalToGlobalMappingApply(mapping,n,idxin,idxout);CHKERRQ(ierr);
6113b9aefa3SBarry Smith   ierr = ISRestoreIndices(is,&idxin);CHKERRQ(ierr);
612543f3098SMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),n,idxout,PETSC_OWN_POINTER,newis);CHKERRQ(ierr);
6133a40ed3dSBarry Smith   PetscFunctionReturn(0);
61490f02eecSBarry Smith }
61590f02eecSBarry Smith 
616b89cb25eSSatish Balay /*@
6173acfe500SLois Curfman McInnes    ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering
6183acfe500SLois Curfman McInnes    and converts them to the global numbering.
61990f02eecSBarry Smith 
620b9cd556bSLois Curfman McInnes    Not collective
621b9cd556bSLois Curfman McInnes 
622bb25748dSBarry Smith    Input Parameters:
623b9cd556bSLois Curfman McInnes +  mapping - the local to global mapping context
624bb25748dSBarry Smith .  N - number of integers
625b9cd556bSLois Curfman McInnes -  in - input indices in local numbering
626bb25748dSBarry Smith 
627bb25748dSBarry Smith    Output Parameter:
628bb25748dSBarry Smith .  out - indices in global numbering
629bb25748dSBarry Smith 
630b9cd556bSLois Curfman McInnes    Notes:
631b9cd556bSLois Curfman McInnes    The in and out array parameters may be identical.
632d4bb536fSBarry Smith 
633a997ad1aSLois Curfman McInnes    Level: advanced
634a997ad1aSLois Curfman McInnes 
63545b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
6360752156aSBarry Smith           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
637d4bb536fSBarry Smith           AOPetscToApplication(), ISGlobalToLocalMappingApply()
638bb25748dSBarry Smith 
639273d9f13SBarry Smith     Concepts: mapping^local to global
640afcb2eb5SJed Brown @*/
641afcb2eb5SJed Brown PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
642afcb2eb5SJed Brown {
643cbc1caf0SMatthew G. Knepley   PetscInt i,bs,Nmax;
64445b6f7e9SBarry Smith 
64545b6f7e9SBarry Smith   PetscFunctionBegin;
646cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
647cbc1caf0SMatthew G. Knepley   bs   = mapping->bs;
648cbc1caf0SMatthew G. Knepley   Nmax = bs*mapping->n;
64945b6f7e9SBarry Smith   if (bs == 1) {
650cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
65145b6f7e9SBarry Smith     for (i=0; i<N; i++) {
65245b6f7e9SBarry Smith       if (in[i] < 0) {
65345b6f7e9SBarry Smith         out[i] = in[i];
65445b6f7e9SBarry Smith         continue;
65545b6f7e9SBarry Smith       }
656e24637baSBarry 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);
65745b6f7e9SBarry Smith       out[i] = idx[in[i]];
65845b6f7e9SBarry Smith     }
65945b6f7e9SBarry Smith   } else {
660cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
66145b6f7e9SBarry Smith     for (i=0; i<N; i++) {
66245b6f7e9SBarry Smith       if (in[i] < 0) {
66345b6f7e9SBarry Smith         out[i] = in[i];
66445b6f7e9SBarry Smith         continue;
66545b6f7e9SBarry Smith       }
666e24637baSBarry 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);
66745b6f7e9SBarry Smith       out[i] = idx[in[i]/bs]*bs + (in[i] % bs);
66845b6f7e9SBarry Smith     }
66945b6f7e9SBarry Smith   }
67045b6f7e9SBarry Smith   PetscFunctionReturn(0);
67145b6f7e9SBarry Smith }
67245b6f7e9SBarry Smith 
67345b6f7e9SBarry Smith /*@
6746006e8d2SBarry Smith    ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local block numbering and converts them to the global block numbering
67545b6f7e9SBarry Smith 
67645b6f7e9SBarry Smith    Not collective
67745b6f7e9SBarry Smith 
67845b6f7e9SBarry Smith    Input Parameters:
67945b6f7e9SBarry Smith +  mapping - the local to global mapping context
68045b6f7e9SBarry Smith .  N - number of integers
6816006e8d2SBarry Smith -  in - input indices in local block numbering
68245b6f7e9SBarry Smith 
68345b6f7e9SBarry Smith    Output Parameter:
6846006e8d2SBarry Smith .  out - indices in global block numbering
68545b6f7e9SBarry Smith 
68645b6f7e9SBarry Smith    Notes:
68745b6f7e9SBarry Smith    The in and out array parameters may be identical.
68845b6f7e9SBarry Smith 
6896006e8d2SBarry Smith    Example:
6906006e8d2SBarry 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
6916006e8d2SBarry Smith      (the first block) would produce 0 and the mapping applied to 1 (the second block) would produce 3.
6926006e8d2SBarry Smith 
69345b6f7e9SBarry Smith    Level: advanced
69445b6f7e9SBarry Smith 
69545b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
69645b6f7e9SBarry Smith           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
69745b6f7e9SBarry Smith           AOPetscToApplication(), ISGlobalToLocalMappingApply()
69845b6f7e9SBarry Smith 
69945b6f7e9SBarry Smith     Concepts: mapping^local to global
70045b6f7e9SBarry Smith @*/
70145b6f7e9SBarry Smith PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
70245b6f7e9SBarry Smith {
703cbc1caf0SMatthew G. Knepley 
704cbc1caf0SMatthew G. Knepley   PetscFunctionBegin;
705cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
706cbc1caf0SMatthew G. Knepley   {
707afcb2eb5SJed Brown     PetscInt i,Nmax = mapping->n;
708afcb2eb5SJed Brown     const PetscInt *idx = mapping->indices;
709d4bb536fSBarry Smith 
710afcb2eb5SJed Brown     for (i=0; i<N; i++) {
711afcb2eb5SJed Brown       if (in[i] < 0) {
712afcb2eb5SJed Brown         out[i] = in[i];
713afcb2eb5SJed Brown         continue;
714afcb2eb5SJed Brown       }
715e24637baSBarry 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);
716afcb2eb5SJed Brown       out[i] = idx[in[i]];
717afcb2eb5SJed Brown     }
718cbc1caf0SMatthew G. Knepley   }
719afcb2eb5SJed Brown   PetscFunctionReturn(0);
720afcb2eb5SJed Brown }
721d4bb536fSBarry Smith 
7227e99dc12SLawrence Mitchell /*@
723a997ad1aSLois Curfman McInnes     ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers
724a997ad1aSLois Curfman McInnes     specified with a global numbering.
725d4bb536fSBarry Smith 
726b9cd556bSLois Curfman McInnes     Not collective
727b9cd556bSLois Curfman McInnes 
728d4bb536fSBarry Smith     Input Parameters:
729b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
730d4bb536fSBarry Smith .   type - IS_GTOLM_MASK - replaces global indices with no local value with -1
731d4bb536fSBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
732d4bb536fSBarry Smith .   n - number of global indices to map
733b9cd556bSLois Curfman McInnes -   idx - global indices to map
734d4bb536fSBarry Smith 
735d4bb536fSBarry Smith     Output Parameters:
736b9cd556bSLois Curfman McInnes +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
737b9cd556bSLois Curfman McInnes -   idxout - local index of each global index, one must pass in an array long enough
738e182c471SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApply() with
7390298fd71SBarry Smith              idxout == NULL to determine the required length (returned in nout)
740e182c471SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApply()
741e182c471SBarry Smith              a second time to set the values.
742d4bb536fSBarry Smith 
743b9cd556bSLois Curfman McInnes     Notes:
7440298fd71SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
745d4bb536fSBarry Smith 
7469a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
747413f72f0SBarry 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.
748413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
7490f5bd95cSBarry Smith 
750a997ad1aSLois Curfman McInnes     Level: advanced
751a997ad1aSLois Curfman McInnes 
75232fd6b96SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
75332fd6b96SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
75432fd6b96SBarry Smith 
755273d9f13SBarry Smith     Concepts: mapping^global to local
756d4bb536fSBarry Smith 
7579d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),
758413f72f0SBarry Smith           ISLocalToGlobalMappingDestroy()
759d4bb536fSBarry Smith @*/
760413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
761d4bb536fSBarry Smith {
7629d90f715SBarry Smith   PetscErrorCode ierr;
7639d90f715SBarry Smith 
7649d90f715SBarry Smith   PetscFunctionBegin;
7659d90f715SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
766413f72f0SBarry Smith   if (!mapping->data) {
767413f72f0SBarry Smith     ierr = ISGlobalToLocalMappingSetUp(mapping);CHKERRQ(ierr);
7689d90f715SBarry Smith   }
769413f72f0SBarry Smith   ierr = (*mapping->ops->globaltolocalmappingapply)(mapping,type,n,idx,nout,idxout);CHKERRQ(ierr);
7709d90f715SBarry Smith   PetscFunctionReturn(0);
7719d90f715SBarry Smith }
7729d90f715SBarry Smith 
773d4fe737eSStefano Zampini /*@
774d4fe737eSStefano Zampini     ISGlobalToLocalMappingApplyIS - Creates from an IS in the global numbering
775d4fe737eSStefano Zampini     a new index set using the local numbering defined in an ISLocalToGlobalMapping
776d4fe737eSStefano Zampini     context.
777d4fe737eSStefano Zampini 
778d4fe737eSStefano Zampini     Not collective
779d4fe737eSStefano Zampini 
780d4fe737eSStefano Zampini     Input Parameters:
781d4fe737eSStefano Zampini +   mapping - mapping between local and global numbering
7822785b321SStefano Zampini .   type - IS_GTOLM_MASK - replaces global indices with no local value with -1
7832785b321SStefano Zampini            IS_GTOLM_DROP - drops the indices with no local value from the output list
784d4fe737eSStefano Zampini -   is - index set in global numbering
785d4fe737eSStefano Zampini 
786d4fe737eSStefano Zampini     Output Parameters:
787d4fe737eSStefano Zampini .   newis - index set in local numbering
788d4fe737eSStefano Zampini 
789d4fe737eSStefano Zampini     Level: advanced
790d4fe737eSStefano Zampini 
791d4fe737eSStefano Zampini     Concepts: mapping^local to global
792d4fe737eSStefano Zampini 
793d4fe737eSStefano Zampini .seealso: ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(),
794d4fe737eSStefano Zampini           ISLocalToGlobalMappingDestroy()
795d4fe737eSStefano Zampini @*/
796413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,IS is,IS *newis)
797d4fe737eSStefano Zampini {
798d4fe737eSStefano Zampini   PetscErrorCode ierr;
799d4fe737eSStefano Zampini   PetscInt       n,nout,*idxout;
800d4fe737eSStefano Zampini   const PetscInt *idxin;
801d4fe737eSStefano Zampini 
802d4fe737eSStefano Zampini   PetscFunctionBegin;
803d4fe737eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
804d4fe737eSStefano Zampini   PetscValidHeaderSpecific(is,IS_CLASSID,3);
805d4fe737eSStefano Zampini   PetscValidPointer(newis,4);
806d4fe737eSStefano Zampini 
807d4fe737eSStefano Zampini   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
808d4fe737eSStefano Zampini   ierr = ISGetIndices(is,&idxin);CHKERRQ(ierr);
809d4fe737eSStefano Zampini   if (type == IS_GTOLM_MASK) {
810d4fe737eSStefano Zampini     ierr = PetscMalloc1(n,&idxout);CHKERRQ(ierr);
811d4fe737eSStefano Zampini   } else {
812d4fe737eSStefano Zampini     ierr = ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,NULL);CHKERRQ(ierr);
813d4fe737eSStefano Zampini     ierr = PetscMalloc1(nout,&idxout);CHKERRQ(ierr);
814d4fe737eSStefano Zampini   }
815d4fe737eSStefano Zampini   ierr = ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,idxout);CHKERRQ(ierr);
816d4fe737eSStefano Zampini   ierr = ISRestoreIndices(is,&idxin);CHKERRQ(ierr);
817d4fe737eSStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),nout,idxout,PETSC_OWN_POINTER,newis);CHKERRQ(ierr);
818d4fe737eSStefano Zampini   PetscFunctionReturn(0);
819d4fe737eSStefano Zampini }
820d4fe737eSStefano Zampini 
8219d90f715SBarry Smith /*@
8229d90f715SBarry Smith     ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers
8239d90f715SBarry Smith     specified with a block global numbering.
8249d90f715SBarry Smith 
8259d90f715SBarry Smith     Not collective
8269d90f715SBarry Smith 
8279d90f715SBarry Smith     Input Parameters:
8289d90f715SBarry Smith +   mapping - mapping between local and global numbering
8299d90f715SBarry Smith .   type - IS_GTOLM_MASK - replaces global indices with no local value with -1
8309d90f715SBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
8319d90f715SBarry Smith .   n - number of global indices to map
8329d90f715SBarry Smith -   idx - global indices to map
8339d90f715SBarry Smith 
8349d90f715SBarry Smith     Output Parameters:
8359d90f715SBarry Smith +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
8369d90f715SBarry Smith -   idxout - local index of each global index, one must pass in an array long enough
8379d90f715SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApplyBlock() with
8389d90f715SBarry Smith              idxout == NULL to determine the required length (returned in nout)
8399d90f715SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApplyBlock()
8409d90f715SBarry Smith              a second time to set the values.
8419d90f715SBarry Smith 
8429d90f715SBarry Smith     Notes:
8439d90f715SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
8449d90f715SBarry Smith 
8459a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
846413f72f0SBarry 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.
847413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
8489d90f715SBarry Smith 
8499d90f715SBarry Smith     Level: advanced
8509d90f715SBarry Smith 
8519d90f715SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
8529d90f715SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
8539d90f715SBarry Smith 
8549d90f715SBarry Smith     Concepts: mapping^global to local
8559d90f715SBarry Smith 
8569d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(),
857413f72f0SBarry Smith           ISLocalToGlobalMappingDestroy()
8589d90f715SBarry Smith @*/
859413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,
8609d90f715SBarry Smith                                                  PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
8619d90f715SBarry Smith {
8626849ba73SBarry Smith   PetscErrorCode ierr;
863d4bb536fSBarry Smith 
8643a40ed3dSBarry Smith   PetscFunctionBegin;
8650700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
866413f72f0SBarry Smith   if (!mapping->data) {
867413f72f0SBarry Smith     ierr = ISGlobalToLocalMappingSetUp(mapping);CHKERRQ(ierr);
868d4bb536fSBarry Smith   }
869413f72f0SBarry Smith   ierr = (*mapping->ops->globaltolocalmappingapplyblock)(mapping,type,n,idx,nout,idxout);CHKERRQ(ierr);
8703a40ed3dSBarry Smith   PetscFunctionReturn(0);
871d4bb536fSBarry Smith }
87290f02eecSBarry Smith 
873413f72f0SBarry Smith 
87489d82c54SBarry Smith /*@C
8756a818285SBarry Smith     ISLocalToGlobalMappingGetBlockInfo - Gets the neighbor information for each processor and
87689d82c54SBarry Smith      each index shared by more than one processor
87789d82c54SBarry Smith 
87889d82c54SBarry Smith     Collective on ISLocalToGlobalMapping
87989d82c54SBarry Smith 
88089d82c54SBarry Smith     Input Parameters:
88189d82c54SBarry Smith .   mapping - the mapping from local to global indexing
88289d82c54SBarry Smith 
88389d82c54SBarry Smith     Output Parameter:
88489d82c54SBarry Smith +   nproc - number of processors that are connected to this one
88589d82c54SBarry Smith .   proc - neighboring processors
88607b52d57SBarry Smith .   numproc - number of indices for each subdomain (processor)
8873463a7baSJed Brown -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
88889d82c54SBarry Smith 
88989d82c54SBarry Smith     Level: advanced
89089d82c54SBarry Smith 
891273d9f13SBarry Smith     Concepts: mapping^local to global
89289d82c54SBarry Smith 
8932cfcea29SBarry Smith     Fortran Usage:
8942cfcea29SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
8952cfcea29SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
8962cfcea29SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
8972cfcea29SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
8982cfcea29SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
8992cfcea29SBarry Smith 
9002cfcea29SBarry Smith 
90107b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
90207b52d57SBarry Smith           ISLocalToGlobalMappingRestoreInfo()
90389d82c54SBarry Smith @*/
9046a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
90589d82c54SBarry Smith {
9066849ba73SBarry Smith   PetscErrorCode ierr;
907268a049cSStefano Zampini 
908268a049cSStefano Zampini   PetscFunctionBegin;
909268a049cSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
910268a049cSStefano Zampini   if (mapping->info_cached) {
911268a049cSStefano Zampini     *nproc    = mapping->info_nproc;
912268a049cSStefano Zampini     *procs    = mapping->info_procs;
913268a049cSStefano Zampini     *numprocs = mapping->info_numprocs;
914268a049cSStefano Zampini     *indices  = mapping->info_indices;
915268a049cSStefano Zampini   } else {
916268a049cSStefano Zampini     ierr = ISLocalToGlobalMappingGetBlockInfo_Private(mapping,nproc,procs,numprocs,indices);CHKERRQ(ierr);
917268a049cSStefano Zampini   }
918268a049cSStefano Zampini   PetscFunctionReturn(0);
919268a049cSStefano Zampini }
920268a049cSStefano Zampini 
921268a049cSStefano Zampini static PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
922268a049cSStefano Zampini {
923268a049cSStefano Zampini   PetscErrorCode ierr;
92497f1f81fSBarry Smith   PetscMPIInt    size,rank,tag1,tag2,tag3,*len,*source,imdex;
92532dcc486SBarry Smith   PetscInt       i,n = mapping->n,Ng,ng,max = 0,*lindices = mapping->indices;
92632dcc486SBarry Smith   PetscInt       *nprocs,*owner,nsends,*sends,j,*starts,nmax,nrecvs,*recvs,proc;
92797f1f81fSBarry Smith   PetscInt       cnt,scale,*ownedsenders,*nownedsenders,rstart,nowned;
92832dcc486SBarry Smith   PetscInt       node,nownedm,nt,*sends2,nsends2,*starts2,*lens2,*dest,nrecvs2,*starts3,*recvs2,k,*bprocs,*tmp;
92932dcc486SBarry Smith   PetscInt       first_procs,first_numprocs,*first_indices;
93089d82c54SBarry Smith   MPI_Request    *recv_waits,*send_waits;
93130dcb7c9SBarry Smith   MPI_Status     recv_status,*send_status,*recv_statuses;
932ce94432eSBarry Smith   MPI_Comm       comm;
933ace3abfcSBarry Smith   PetscBool      debug = PETSC_FALSE;
93489d82c54SBarry Smith 
93589d82c54SBarry Smith   PetscFunctionBegin;
936ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)mapping,&comm);CHKERRQ(ierr);
93724cf384cSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
93824cf384cSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
93924cf384cSBarry Smith   if (size == 1) {
94024cf384cSBarry Smith     *nproc         = 0;
9410298fd71SBarry Smith     *procs         = NULL;
94295dccacaSBarry Smith     ierr           = PetscNew(numprocs);CHKERRQ(ierr);
9431e2105dcSBarry Smith     (*numprocs)[0] = 0;
94495dccacaSBarry Smith     ierr           = PetscNew(indices);CHKERRQ(ierr);
9450298fd71SBarry Smith     (*indices)[0]  = NULL;
946268a049cSStefano Zampini     /* save info for reuse */
947268a049cSStefano Zampini     mapping->info_nproc = *nproc;
948268a049cSStefano Zampini     mapping->info_procs = *procs;
949268a049cSStefano Zampini     mapping->info_numprocs = *numprocs;
950268a049cSStefano Zampini     mapping->info_indices = *indices;
951268a049cSStefano Zampini     mapping->info_cached = PETSC_TRUE;
95224cf384cSBarry Smith     PetscFunctionReturn(0);
95324cf384cSBarry Smith   }
95424cf384cSBarry Smith 
955c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)mapping)->options,NULL,"-islocaltoglobalmappinggetinfo_debug",&debug,NULL);CHKERRQ(ierr);
95607b52d57SBarry Smith 
9573677ff5aSBarry Smith   /*
9586a818285SBarry Smith     Notes on ISLocalToGlobalMappingGetBlockInfo
9593677ff5aSBarry Smith 
9603677ff5aSBarry Smith     globally owned node - the nodes that have been assigned to this processor in global
9613677ff5aSBarry Smith            numbering, just for this routine.
9623677ff5aSBarry Smith 
9633677ff5aSBarry Smith     nontrivial globally owned node - node assigned to this processor that is on a subdomain
9643677ff5aSBarry Smith            boundary (i.e. is has more than one local owner)
9653677ff5aSBarry Smith 
9663677ff5aSBarry Smith     locally owned node - node that exists on this processors subdomain
9673677ff5aSBarry Smith 
9683677ff5aSBarry Smith     nontrivial locally owned node - node that is not in the interior (i.e. has more than one
9693677ff5aSBarry Smith            local subdomain
9703677ff5aSBarry Smith   */
97124cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag1);CHKERRQ(ierr);
97224cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag2);CHKERRQ(ierr);
97324cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag3);CHKERRQ(ierr);
97489d82c54SBarry Smith 
97589d82c54SBarry Smith   for (i=0; i<n; i++) {
97689d82c54SBarry Smith     if (lindices[i] > max) max = lindices[i];
97789d82c54SBarry Smith   }
978b2566f29SBarry Smith   ierr   = MPIU_Allreduce(&max,&Ng,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr);
97978058e43SBarry Smith   Ng++;
98089d82c54SBarry Smith   ierr   = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
98189d82c54SBarry Smith   ierr   = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
982bc8ff85bSBarry Smith   scale  = Ng/size + 1;
983a2e34c3dSBarry Smith   ng     = scale; if (rank == size-1) ng = Ng - scale*(size-1); ng = PetscMax(1,ng);
984caba0dd0SBarry Smith   rstart = scale*rank;
98589d82c54SBarry Smith 
98689d82c54SBarry Smith   /* determine ownership ranges of global indices */
987785e854fSJed Brown   ierr = PetscMalloc1(2*size,&nprocs);CHKERRQ(ierr);
98832dcc486SBarry Smith   ierr = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr);
98989d82c54SBarry Smith 
99089d82c54SBarry Smith   /* determine owners of each local node  */
991785e854fSJed Brown   ierr = PetscMalloc1(n,&owner);CHKERRQ(ierr);
99289d82c54SBarry Smith   for (i=0; i<n; i++) {
9933677ff5aSBarry Smith     proc             = lindices[i]/scale; /* processor that globally owns this index */
99427c402fcSBarry Smith     nprocs[2*proc+1] = 1;                 /* processor globally owns at least one of ours */
9953677ff5aSBarry Smith     owner[i]         = proc;
99627c402fcSBarry Smith     nprocs[2*proc]++;                     /* count of how many that processor globally owns of ours */
99789d82c54SBarry Smith   }
99827c402fcSBarry Smith   nsends = 0; for (i=0; i<size; i++) nsends += nprocs[2*i+1];
9997904a332SBarry Smith   ierr = PetscInfo1(mapping,"Number of global owners for my local data %D\n",nsends);CHKERRQ(ierr);
100089d82c54SBarry Smith 
100189d82c54SBarry Smith   /* inform other processors of number of messages and max length*/
100227c402fcSBarry Smith   ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr);
10037904a332SBarry Smith   ierr = PetscInfo1(mapping,"Number of local owners for my global data %D\n",nrecvs);CHKERRQ(ierr);
100489d82c54SBarry Smith 
100589d82c54SBarry Smith   /* post receives for owned rows */
1006785e854fSJed Brown   ierr = PetscMalloc1((2*nrecvs+1)*(nmax+1),&recvs);CHKERRQ(ierr);
1007854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&recv_waits);CHKERRQ(ierr);
100889d82c54SBarry Smith   for (i=0; i<nrecvs; i++) {
100932dcc486SBarry Smith     ierr = MPI_Irecv(recvs+2*nmax*i,2*nmax,MPIU_INT,MPI_ANY_SOURCE,tag1,comm,recv_waits+i);CHKERRQ(ierr);
101089d82c54SBarry Smith   }
101189d82c54SBarry Smith 
101289d82c54SBarry Smith   /* pack messages containing lists of local nodes to owners */
1013854ce69bSBarry Smith   ierr      = PetscMalloc1(2*n+1,&sends);CHKERRQ(ierr);
1014854ce69bSBarry Smith   ierr      = PetscMalloc1(size+1,&starts);CHKERRQ(ierr);
101589d82c54SBarry Smith   starts[0] = 0;
1016f6e5521dSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
101789d82c54SBarry Smith   for (i=0; i<n; i++) {
101889d82c54SBarry Smith     sends[starts[owner[i]]++] = lindices[i];
101930dcb7c9SBarry Smith     sends[starts[owner[i]]++] = i;
102089d82c54SBarry Smith   }
102189d82c54SBarry Smith   ierr = PetscFree(owner);CHKERRQ(ierr);
102289d82c54SBarry Smith   starts[0] = 0;
1023f6e5521dSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
102489d82c54SBarry Smith 
102589d82c54SBarry Smith   /* send the messages */
1026854ce69bSBarry Smith   ierr = PetscMalloc1(nsends+1,&send_waits);CHKERRQ(ierr);
1027854ce69bSBarry Smith   ierr = PetscMalloc1(nsends+1,&dest);CHKERRQ(ierr);
102889d82c54SBarry Smith   cnt = 0;
102989d82c54SBarry Smith   for (i=0; i<size; i++) {
103027c402fcSBarry Smith     if (nprocs[2*i]) {
103132dcc486SBarry Smith       ierr      = MPI_Isend(sends+starts[i],2*nprocs[2*i],MPIU_INT,i,tag1,comm,send_waits+cnt);CHKERRQ(ierr);
103230dcb7c9SBarry Smith       dest[cnt] = i;
103389d82c54SBarry Smith       cnt++;
103489d82c54SBarry Smith     }
103589d82c54SBarry Smith   }
103689d82c54SBarry Smith   ierr = PetscFree(starts);CHKERRQ(ierr);
103789d82c54SBarry Smith 
103889d82c54SBarry Smith   /* wait on receives */
1039854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&source);CHKERRQ(ierr);
1040854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&len);CHKERRQ(ierr);
104189d82c54SBarry Smith   cnt  = nrecvs;
1042854ce69bSBarry Smith   ierr = PetscMalloc1(ng+1,&nownedsenders);CHKERRQ(ierr);
104332dcc486SBarry Smith   ierr = PetscMemzero(nownedsenders,ng*sizeof(PetscInt));CHKERRQ(ierr);
104489d82c54SBarry Smith   while (cnt) {
104589d82c54SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
104689d82c54SBarry Smith     /* unpack receives into our local space */
104732dcc486SBarry Smith     ierr          = MPI_Get_count(&recv_status,MPIU_INT,&len[imdex]);CHKERRQ(ierr);
104889d82c54SBarry Smith     source[imdex] = recv_status.MPI_SOURCE;
104930dcb7c9SBarry Smith     len[imdex]    = len[imdex]/2;
1050caba0dd0SBarry Smith     /* count how many local owners for each of my global owned indices */
105130dcb7c9SBarry Smith     for (i=0; i<len[imdex]; i++) nownedsenders[recvs[2*imdex*nmax+2*i]-rstart]++;
105289d82c54SBarry Smith     cnt--;
105389d82c54SBarry Smith   }
105489d82c54SBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
105589d82c54SBarry Smith 
105630dcb7c9SBarry Smith   /* count how many globally owned indices are on an edge multiplied by how many processors own them. */
1057bc8ff85bSBarry Smith   nowned  = 0;
1058bc8ff85bSBarry Smith   nownedm = 0;
1059bc8ff85bSBarry Smith   for (i=0; i<ng; i++) {
1060bc8ff85bSBarry Smith     if (nownedsenders[i] > 1) {nownedm += nownedsenders[i]; nowned++;}
1061bc8ff85bSBarry Smith   }
1062bc8ff85bSBarry Smith 
1063bc8ff85bSBarry Smith   /* create single array to contain rank of all local owners of each globally owned index */
1064854ce69bSBarry Smith   ierr      = PetscMalloc1(nownedm+1,&ownedsenders);CHKERRQ(ierr);
1065854ce69bSBarry Smith   ierr      = PetscMalloc1(ng+1,&starts);CHKERRQ(ierr);
1066bc8ff85bSBarry Smith   starts[0] = 0;
1067bc8ff85bSBarry Smith   for (i=1; i<ng; i++) {
1068bc8ff85bSBarry Smith     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
1069bc8ff85bSBarry Smith     else starts[i] = starts[i-1];
1070bc8ff85bSBarry Smith   }
1071bc8ff85bSBarry Smith 
107230dcb7c9SBarry Smith   /* for each nontrival globally owned node list all arriving processors */
1073bc8ff85bSBarry Smith   for (i=0; i<nrecvs; i++) {
1074bc8ff85bSBarry Smith     for (j=0; j<len[i]; j++) {
107530dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
1076f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) ownedsenders[starts[node]++] = source[i];
1077bc8ff85bSBarry Smith     }
1078bc8ff85bSBarry Smith   }
1079bc8ff85bSBarry Smith 
108007b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
108130dcb7c9SBarry Smith     starts[0] = 0;
108230dcb7c9SBarry Smith     for (i=1; i<ng; i++) {
108330dcb7c9SBarry Smith       if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
108430dcb7c9SBarry Smith       else starts[i] = starts[i-1];
108530dcb7c9SBarry Smith     }
108630dcb7c9SBarry Smith     for (i=0; i<ng; i++) {
108730dcb7c9SBarry Smith       if (nownedsenders[i] > 1) {
10887904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"[%d] global node %D local owner processors: ",rank,i+rstart);CHKERRQ(ierr);
108930dcb7c9SBarry Smith         for (j=0; j<nownedsenders[i]; j++) {
10907904a332SBarry Smith           ierr = PetscSynchronizedPrintf(comm,"%D ",ownedsenders[starts[i]+j]);CHKERRQ(ierr);
109130dcb7c9SBarry Smith         }
109230dcb7c9SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
109330dcb7c9SBarry Smith       }
109430dcb7c9SBarry Smith     }
10950ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
109607b52d57SBarry Smith   } /* -----------------------------------  */
109730dcb7c9SBarry Smith 
10983677ff5aSBarry Smith   /* wait on original sends */
10993a96401aSBarry Smith   if (nsends) {
1100785e854fSJed Brown     ierr = PetscMalloc1(nsends,&send_status);CHKERRQ(ierr);
11013a96401aSBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
11023a96401aSBarry Smith     ierr = PetscFree(send_status);CHKERRQ(ierr);
11033a96401aSBarry Smith   }
110489d82c54SBarry Smith   ierr = PetscFree(send_waits);CHKERRQ(ierr);
11053a96401aSBarry Smith   ierr = PetscFree(sends);CHKERRQ(ierr);
11063677ff5aSBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
11073677ff5aSBarry Smith 
11083677ff5aSBarry Smith   /* pack messages to send back to local owners */
110930dcb7c9SBarry Smith   starts[0] = 0;
111030dcb7c9SBarry Smith   for (i=1; i<ng; i++) {
111130dcb7c9SBarry Smith     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
111230dcb7c9SBarry Smith     else starts[i] = starts[i-1];
111330dcb7c9SBarry Smith   }
111430dcb7c9SBarry Smith   nsends2 = nrecvs;
1115854ce69bSBarry Smith   ierr    = PetscMalloc1(nsends2+1,&nprocs);CHKERRQ(ierr); /* length of each message */
111630dcb7c9SBarry Smith   for (i=0; i<nrecvs; i++) {
111730dcb7c9SBarry Smith     nprocs[i] = 1;
111830dcb7c9SBarry Smith     for (j=0; j<len[i]; j++) {
111930dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
1120f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) nprocs[i] += 2 + nownedsenders[node];
112130dcb7c9SBarry Smith     }
112230dcb7c9SBarry Smith   }
1123f6e5521dSKarl Rupp   nt = 0;
1124f6e5521dSKarl Rupp   for (i=0; i<nsends2; i++) nt += nprocs[i];
1125f6e5521dSKarl Rupp 
1126854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,&sends2);CHKERRQ(ierr);
1127854ce69bSBarry Smith   ierr = PetscMalloc1(nsends2+1,&starts2);CHKERRQ(ierr);
1128f6e5521dSKarl Rupp 
1129f6e5521dSKarl Rupp   starts2[0] = 0;
1130f6e5521dSKarl Rupp   for (i=1; i<nsends2; i++) starts2[i] = starts2[i-1] + nprocs[i-1];
113130dcb7c9SBarry Smith   /*
113230dcb7c9SBarry Smith      Each message is 1 + nprocs[i] long, and consists of
113330dcb7c9SBarry Smith        (0) the number of nodes being sent back
113430dcb7c9SBarry Smith        (1) the local node number,
113530dcb7c9SBarry Smith        (2) the number of processors sharing it,
113630dcb7c9SBarry Smith        (3) the processors sharing it
113730dcb7c9SBarry Smith   */
113830dcb7c9SBarry Smith   for (i=0; i<nsends2; i++) {
113930dcb7c9SBarry Smith     cnt = 1;
114030dcb7c9SBarry Smith     sends2[starts2[i]] = 0;
114130dcb7c9SBarry Smith     for (j=0; j<len[i]; j++) {
114230dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
114330dcb7c9SBarry Smith       if (nownedsenders[node] > 1) {
114430dcb7c9SBarry Smith         sends2[starts2[i]]++;
114530dcb7c9SBarry Smith         sends2[starts2[i]+cnt++] = recvs[2*i*nmax+2*j+1];
114630dcb7c9SBarry Smith         sends2[starts2[i]+cnt++] = nownedsenders[node];
114732dcc486SBarry Smith         ierr = PetscMemcpy(&sends2[starts2[i]+cnt],&ownedsenders[starts[node]],nownedsenders[node]*sizeof(PetscInt));CHKERRQ(ierr);
114830dcb7c9SBarry Smith         cnt += nownedsenders[node];
114930dcb7c9SBarry Smith       }
115030dcb7c9SBarry Smith     }
115130dcb7c9SBarry Smith   }
115230dcb7c9SBarry Smith 
115330dcb7c9SBarry Smith   /* receive the message lengths */
115430dcb7c9SBarry Smith   nrecvs2 = nsends;
1155854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&lens2);CHKERRQ(ierr);
1156854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&starts3);CHKERRQ(ierr);
1157854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&recv_waits);CHKERRQ(ierr);
115830dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
1159d44834fbSBarry Smith     ierr = MPI_Irecv(&lens2[i],1,MPIU_INT,dest[i],tag2,comm,recv_waits+i);CHKERRQ(ierr);
116030dcb7c9SBarry Smith   }
1161d44834fbSBarry Smith 
11628a8e0b3aSBarry Smith   /* send the message lengths */
11638a8e0b3aSBarry Smith   for (i=0; i<nsends2; i++) {
11648a8e0b3aSBarry Smith     ierr = MPI_Send(&nprocs[i],1,MPIU_INT,source[i],tag2,comm);CHKERRQ(ierr);
11658a8e0b3aSBarry Smith   }
11668a8e0b3aSBarry Smith 
1167d44834fbSBarry Smith   /* wait on receives of lens */
11680c468ba9SBarry Smith   if (nrecvs2) {
1169785e854fSJed Brown     ierr = PetscMalloc1(nrecvs2,&recv_statuses);CHKERRQ(ierr);
1170d44834fbSBarry Smith     ierr = MPI_Waitall(nrecvs2,recv_waits,recv_statuses);CHKERRQ(ierr);
1171d44834fbSBarry Smith     ierr = PetscFree(recv_statuses);CHKERRQ(ierr);
11720c468ba9SBarry Smith   }
1173a2ea699eSBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
1174d44834fbSBarry Smith 
117530dcb7c9SBarry Smith   starts3[0] = 0;
1176d44834fbSBarry Smith   nt         = 0;
117730dcb7c9SBarry Smith   for (i=0; i<nrecvs2-1; i++) {
117830dcb7c9SBarry Smith     starts3[i+1] = starts3[i] + lens2[i];
1179d44834fbSBarry Smith     nt          += lens2[i];
118030dcb7c9SBarry Smith   }
118176466f69SStefano Zampini   if (nrecvs2) nt += lens2[nrecvs2-1];
1182d44834fbSBarry Smith 
1183854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,&recvs2);CHKERRQ(ierr);
1184854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs2+1,&recv_waits);CHKERRQ(ierr);
118552b72c4aSBarry Smith   for (i=0; i<nrecvs2; i++) {
118632dcc486SBarry Smith     ierr = MPI_Irecv(recvs2+starts3[i],lens2[i],MPIU_INT,dest[i],tag3,comm,recv_waits+i);CHKERRQ(ierr);
118730dcb7c9SBarry Smith   }
118830dcb7c9SBarry Smith 
118930dcb7c9SBarry Smith   /* send the messages */
1190854ce69bSBarry Smith   ierr = PetscMalloc1(nsends2+1,&send_waits);CHKERRQ(ierr);
119130dcb7c9SBarry Smith   for (i=0; i<nsends2; i++) {
119232dcc486SBarry Smith     ierr = MPI_Isend(sends2+starts2[i],nprocs[i],MPIU_INT,source[i],tag3,comm,send_waits+i);CHKERRQ(ierr);
119330dcb7c9SBarry Smith   }
119430dcb7c9SBarry Smith 
119530dcb7c9SBarry Smith   /* wait on receives */
11960c468ba9SBarry Smith   if (nrecvs2) {
1197785e854fSJed Brown     ierr = PetscMalloc1(nrecvs2,&recv_statuses);CHKERRQ(ierr);
119830dcb7c9SBarry Smith     ierr = MPI_Waitall(nrecvs2,recv_waits,recv_statuses);CHKERRQ(ierr);
119930dcb7c9SBarry Smith     ierr = PetscFree(recv_statuses);CHKERRQ(ierr);
12000c468ba9SBarry Smith   }
120130dcb7c9SBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
120230dcb7c9SBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
120330dcb7c9SBarry Smith 
120407b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
120530dcb7c9SBarry Smith     cnt = 0;
120630dcb7c9SBarry Smith     for (i=0; i<nrecvs2; i++) {
120730dcb7c9SBarry Smith       nt = recvs2[cnt++];
120830dcb7c9SBarry Smith       for (j=0; j<nt; j++) {
12097904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"[%d] local node %D number of subdomains %D: ",rank,recvs2[cnt],recvs2[cnt+1]);CHKERRQ(ierr);
121030dcb7c9SBarry Smith         for (k=0; k<recvs2[cnt+1]; k++) {
12117904a332SBarry Smith           ierr = PetscSynchronizedPrintf(comm,"%D ",recvs2[cnt+2+k]);CHKERRQ(ierr);
121230dcb7c9SBarry Smith         }
121330dcb7c9SBarry Smith         cnt += 2 + recvs2[cnt+1];
121430dcb7c9SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
121530dcb7c9SBarry Smith       }
121630dcb7c9SBarry Smith     }
12170ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
121807b52d57SBarry Smith   } /* -----------------------------------  */
121930dcb7c9SBarry Smith 
122030dcb7c9SBarry Smith   /* count number subdomains for each local node */
1221785e854fSJed Brown   ierr = PetscMalloc1(size,&nprocs);CHKERRQ(ierr);
122232dcc486SBarry Smith   ierr = PetscMemzero(nprocs,size*sizeof(PetscInt));CHKERRQ(ierr);
122330dcb7c9SBarry Smith   cnt  = 0;
122430dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
122530dcb7c9SBarry Smith     nt = recvs2[cnt++];
122630dcb7c9SBarry Smith     for (j=0; j<nt; j++) {
1227f6e5521dSKarl Rupp       for (k=0; k<recvs2[cnt+1]; k++) nprocs[recvs2[cnt+2+k]]++;
122830dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt+1];
122930dcb7c9SBarry Smith     }
123030dcb7c9SBarry Smith   }
123130dcb7c9SBarry Smith   nt = 0; for (i=0; i<size; i++) nt += (nprocs[i] > 0);
123230dcb7c9SBarry Smith   *nproc    = nt;
1233854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,procs);CHKERRQ(ierr);
1234854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,numprocs);CHKERRQ(ierr);
1235854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,indices);CHKERRQ(ierr);
12360298fd71SBarry Smith   for (i=0;i<nt+1;i++) (*indices)[i]=NULL;
1237785e854fSJed Brown   ierr = PetscMalloc1(size,&bprocs);CHKERRQ(ierr);
123830dcb7c9SBarry Smith   cnt  = 0;
123930dcb7c9SBarry Smith   for (i=0; i<size; i++) {
124030dcb7c9SBarry Smith     if (nprocs[i] > 0) {
124130dcb7c9SBarry Smith       bprocs[i]        = cnt;
124230dcb7c9SBarry Smith       (*procs)[cnt]    = i;
124330dcb7c9SBarry Smith       (*numprocs)[cnt] = nprocs[i];
1244785e854fSJed Brown       ierr             = PetscMalloc1(nprocs[i],&(*indices)[cnt]);CHKERRQ(ierr);
124530dcb7c9SBarry Smith       cnt++;
124630dcb7c9SBarry Smith     }
124730dcb7c9SBarry Smith   }
124830dcb7c9SBarry Smith 
124930dcb7c9SBarry Smith   /* make the list of subdomains for each nontrivial local node */
125032dcc486SBarry Smith   ierr = PetscMemzero(*numprocs,nt*sizeof(PetscInt));CHKERRQ(ierr);
125130dcb7c9SBarry Smith   cnt  = 0;
125230dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
125330dcb7c9SBarry Smith     nt = recvs2[cnt++];
125430dcb7c9SBarry Smith     for (j=0; j<nt; j++) {
1255f6e5521dSKarl Rupp       for (k=0; k<recvs2[cnt+1]; k++) (*indices)[bprocs[recvs2[cnt+2+k]]][(*numprocs)[bprocs[recvs2[cnt+2+k]]]++] = recvs2[cnt];
125630dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt+1];
125730dcb7c9SBarry Smith     }
125830dcb7c9SBarry Smith   }
125930dcb7c9SBarry Smith   ierr = PetscFree(bprocs);CHKERRQ(ierr);
126007b52d57SBarry Smith   ierr = PetscFree(recvs2);CHKERRQ(ierr);
126130dcb7c9SBarry Smith 
126207b52d57SBarry Smith   /* sort the node indexing by their global numbers */
126307b52d57SBarry Smith   nt = *nproc;
126407b52d57SBarry Smith   for (i=0; i<nt; i++) {
1265854ce69bSBarry Smith     ierr = PetscMalloc1((*numprocs)[i],&tmp);CHKERRQ(ierr);
1266f6e5521dSKarl Rupp     for (j=0; j<(*numprocs)[i]; j++) tmp[j] = lindices[(*indices)[i][j]];
126707b52d57SBarry Smith     ierr = PetscSortIntWithArray((*numprocs)[i],tmp,(*indices)[i]);CHKERRQ(ierr);
126807b52d57SBarry Smith     ierr = PetscFree(tmp);CHKERRQ(ierr);
126907b52d57SBarry Smith   }
127007b52d57SBarry Smith 
127107b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
127230dcb7c9SBarry Smith     nt = *nproc;
127330dcb7c9SBarry Smith     for (i=0; i<nt; i++) {
12747904a332SBarry Smith       ierr = PetscSynchronizedPrintf(comm,"[%d] subdomain %D number of indices %D: ",rank,(*procs)[i],(*numprocs)[i]);CHKERRQ(ierr);
127530dcb7c9SBarry Smith       for (j=0; j<(*numprocs)[i]; j++) {
12767904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"%D ",(*indices)[i][j]);CHKERRQ(ierr);
127730dcb7c9SBarry Smith       }
127830dcb7c9SBarry Smith       ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
127930dcb7c9SBarry Smith     }
12800ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
128107b52d57SBarry Smith   } /* -----------------------------------  */
128230dcb7c9SBarry Smith 
128330dcb7c9SBarry Smith   /* wait on sends */
128430dcb7c9SBarry Smith   if (nsends2) {
1285785e854fSJed Brown     ierr = PetscMalloc1(nsends2,&send_status);CHKERRQ(ierr);
128630dcb7c9SBarry Smith     ierr = MPI_Waitall(nsends2,send_waits,send_status);CHKERRQ(ierr);
128730dcb7c9SBarry Smith     ierr = PetscFree(send_status);CHKERRQ(ierr);
128830dcb7c9SBarry Smith   }
128930dcb7c9SBarry Smith 
129030dcb7c9SBarry Smith   ierr = PetscFree(starts3);CHKERRQ(ierr);
129130dcb7c9SBarry Smith   ierr = PetscFree(dest);CHKERRQ(ierr);
129230dcb7c9SBarry Smith   ierr = PetscFree(send_waits);CHKERRQ(ierr);
12933677ff5aSBarry Smith 
1294bc8ff85bSBarry Smith   ierr = PetscFree(nownedsenders);CHKERRQ(ierr);
1295bc8ff85bSBarry Smith   ierr = PetscFree(ownedsenders);CHKERRQ(ierr);
1296bc8ff85bSBarry Smith   ierr = PetscFree(starts);CHKERRQ(ierr);
129730dcb7c9SBarry Smith   ierr = PetscFree(starts2);CHKERRQ(ierr);
129830dcb7c9SBarry Smith   ierr = PetscFree(lens2);CHKERRQ(ierr);
129989d82c54SBarry Smith 
130089d82c54SBarry Smith   ierr = PetscFree(source);CHKERRQ(ierr);
130197f1f81fSBarry Smith   ierr = PetscFree(len);CHKERRQ(ierr);
130289d82c54SBarry Smith   ierr = PetscFree(recvs);CHKERRQ(ierr);
13033a96401aSBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
130430dcb7c9SBarry Smith   ierr = PetscFree(sends2);CHKERRQ(ierr);
130524cf384cSBarry Smith 
130624cf384cSBarry Smith   /* put the information about myself as the first entry in the list */
130724cf384cSBarry Smith   first_procs    = (*procs)[0];
130824cf384cSBarry Smith   first_numprocs = (*numprocs)[0];
130924cf384cSBarry Smith   first_indices  = (*indices)[0];
131024cf384cSBarry Smith   for (i=0; i<*nproc; i++) {
131124cf384cSBarry Smith     if ((*procs)[i] == rank) {
131224cf384cSBarry Smith       (*procs)[0]    = (*procs)[i];
131324cf384cSBarry Smith       (*numprocs)[0] = (*numprocs)[i];
131424cf384cSBarry Smith       (*indices)[0]  = (*indices)[i];
131524cf384cSBarry Smith       (*procs)[i]    = first_procs;
131624cf384cSBarry Smith       (*numprocs)[i] = first_numprocs;
131724cf384cSBarry Smith       (*indices)[i]  = first_indices;
131824cf384cSBarry Smith       break;
131924cf384cSBarry Smith     }
132024cf384cSBarry Smith   }
1321268a049cSStefano Zampini 
1322268a049cSStefano Zampini   /* save info for reuse */
1323268a049cSStefano Zampini   mapping->info_nproc = *nproc;
1324268a049cSStefano Zampini   mapping->info_procs = *procs;
1325268a049cSStefano Zampini   mapping->info_numprocs = *numprocs;
1326268a049cSStefano Zampini   mapping->info_indices = *indices;
1327268a049cSStefano Zampini   mapping->info_cached = PETSC_TRUE;
132889d82c54SBarry Smith   PetscFunctionReturn(0);
132989d82c54SBarry Smith }
133089d82c54SBarry Smith 
13316a818285SBarry Smith /*@C
13326a818285SBarry Smith     ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by ISLocalToGlobalMappingGetBlockInfo()
13336a818285SBarry Smith 
13346a818285SBarry Smith     Collective on ISLocalToGlobalMapping
13356a818285SBarry Smith 
13366a818285SBarry Smith     Input Parameters:
13376a818285SBarry Smith .   mapping - the mapping from local to global indexing
13386a818285SBarry Smith 
13396a818285SBarry Smith     Output Parameter:
13406a818285SBarry Smith +   nproc - number of processors that are connected to this one
13416a818285SBarry Smith .   proc - neighboring processors
13426a818285SBarry Smith .   numproc - number of indices for each processor
13436a818285SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
13446a818285SBarry Smith 
13456a818285SBarry Smith     Level: advanced
13466a818285SBarry Smith 
13476a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
13486a818285SBarry Smith           ISLocalToGlobalMappingGetInfo()
13496a818285SBarry Smith @*/
13506a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
13516a818285SBarry Smith {
13526a818285SBarry Smith   PetscErrorCode ierr;
13536a818285SBarry Smith 
13546a818285SBarry Smith   PetscFunctionBegin;
1355cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1356268a049cSStefano Zampini   if (mapping->info_free) {
13576a818285SBarry Smith     ierr = PetscFree(*numprocs);CHKERRQ(ierr);
13586a818285SBarry Smith     if (*indices) {
1359268a049cSStefano Zampini       PetscInt i;
1360268a049cSStefano Zampini 
13616a818285SBarry Smith       ierr = PetscFree((*indices)[0]);CHKERRQ(ierr);
13626a818285SBarry Smith       for (i=1; i<*nproc; i++) {
13636a818285SBarry Smith         ierr = PetscFree((*indices)[i]);CHKERRQ(ierr);
13646a818285SBarry Smith       }
13656a818285SBarry Smith       ierr = PetscFree(*indices);CHKERRQ(ierr);
13666a818285SBarry Smith     }
1367268a049cSStefano Zampini   }
1368268a049cSStefano Zampini   *nproc    = 0;
1369268a049cSStefano Zampini   *procs    = NULL;
1370268a049cSStefano Zampini   *numprocs = NULL;
1371268a049cSStefano Zampini   *indices  = NULL;
13726a818285SBarry Smith   PetscFunctionReturn(0);
13736a818285SBarry Smith }
13746a818285SBarry Smith 
13756a818285SBarry Smith /*@C
13766a818285SBarry Smith     ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each processor and
13776a818285SBarry Smith      each index shared by more than one processor
13786a818285SBarry Smith 
13796a818285SBarry Smith     Collective on ISLocalToGlobalMapping
13806a818285SBarry Smith 
13816a818285SBarry Smith     Input Parameters:
13826a818285SBarry Smith .   mapping - the mapping from local to global indexing
13836a818285SBarry Smith 
13846a818285SBarry Smith     Output Parameter:
13856a818285SBarry Smith +   nproc - number of processors that are connected to this one
13866a818285SBarry Smith .   proc - neighboring processors
13876a818285SBarry Smith .   numproc - number of indices for each subdomain (processor)
13886a818285SBarry Smith -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
13896a818285SBarry Smith 
13906a818285SBarry Smith     Level: advanced
13916a818285SBarry Smith 
13926a818285SBarry Smith     Concepts: mapping^local to global
13936a818285SBarry Smith 
13946a818285SBarry Smith     Fortran Usage:
13956a818285SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
13966a818285SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
13976a818285SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
13986a818285SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
13996a818285SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
14006a818285SBarry Smith 
14016a818285SBarry Smith 
14026a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
14036a818285SBarry Smith           ISLocalToGlobalMappingRestoreInfo()
14046a818285SBarry Smith @*/
14056a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
14066a818285SBarry Smith {
14076a818285SBarry Smith   PetscErrorCode ierr;
1408268a049cSStefano Zampini   PetscInt       **bindices = NULL,*bnumprocs = NULL,bs = mapping->bs,i,j,k;
14096a818285SBarry Smith 
14106a818285SBarry Smith   PetscFunctionBegin;
14116a818285SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1412268a049cSStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockInfo(mapping,nproc,procs,&bnumprocs,&bindices);CHKERRQ(ierr);
1413268a049cSStefano Zampini   if (bs > 1) { /* we need to expand the cached info */
1414732f65e3SBarry Smith     ierr = PetscCalloc1(*nproc,&*indices);CHKERRQ(ierr);
1415268a049cSStefano Zampini     ierr = PetscCalloc1(*nproc,&*numprocs);CHKERRQ(ierr);
14166a818285SBarry Smith     for (i=0; i<*nproc; i++) {
1417268a049cSStefano Zampini       ierr = PetscMalloc1(bs*bnumprocs[i],&(*indices)[i]);CHKERRQ(ierr);
1418268a049cSStefano Zampini       for (j=0; j<bnumprocs[i]; j++) {
14196a818285SBarry Smith         for (k=0; k<bs; k++) {
14206a818285SBarry Smith           (*indices)[i][j*bs+k] = bs*bindices[i][j] + k;
14216a818285SBarry Smith         }
14226a818285SBarry Smith       }
1423268a049cSStefano Zampini       (*numprocs)[i] = bnumprocs[i]*bs;
14246a818285SBarry Smith     }
1425268a049cSStefano Zampini     mapping->info_free = PETSC_TRUE;
1426268a049cSStefano Zampini   } else {
1427268a049cSStefano Zampini     *numprocs = bnumprocs;
1428268a049cSStefano Zampini     *indices  = bindices;
14296a818285SBarry Smith   }
14306a818285SBarry Smith   PetscFunctionReturn(0);
14316a818285SBarry Smith }
14326a818285SBarry Smith 
143307b52d57SBarry Smith /*@C
143407b52d57SBarry Smith     ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by ISLocalToGlobalMappingGetInfo()
143589d82c54SBarry Smith 
143607b52d57SBarry Smith     Collective on ISLocalToGlobalMapping
143707b52d57SBarry Smith 
143807b52d57SBarry Smith     Input Parameters:
143907b52d57SBarry Smith .   mapping - the mapping from local to global indexing
144007b52d57SBarry Smith 
144107b52d57SBarry Smith     Output Parameter:
144207b52d57SBarry Smith +   nproc - number of processors that are connected to this one
144307b52d57SBarry Smith .   proc - neighboring processors
144407b52d57SBarry Smith .   numproc - number of indices for each processor
144507b52d57SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
144607b52d57SBarry Smith 
144707b52d57SBarry Smith     Level: advanced
144807b52d57SBarry Smith 
144907b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
145007b52d57SBarry Smith           ISLocalToGlobalMappingGetInfo()
145107b52d57SBarry Smith @*/
14527087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
145307b52d57SBarry Smith {
14546849ba73SBarry Smith   PetscErrorCode ierr;
145507b52d57SBarry Smith 
145607b52d57SBarry Smith   PetscFunctionBegin;
14576a818285SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockInfo(mapping,nproc,procs,numprocs,indices);CHKERRQ(ierr);
145807b52d57SBarry Smith   PetscFunctionReturn(0);
145907b52d57SBarry Smith }
146086994e45SJed Brown 
146186994e45SJed Brown /*@C
1462107e9a97SBarry Smith    ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped
146386994e45SJed Brown 
146486994e45SJed Brown    Not Collective
146586994e45SJed Brown 
146686994e45SJed Brown    Input Arguments:
146786994e45SJed Brown . ltog - local to global mapping
146886994e45SJed Brown 
146986994e45SJed Brown    Output Arguments:
1470565245c5SBarry Smith . array - array of indices, the length of this array may be obtained with ISLocalToGlobalMappingGetSize()
147186994e45SJed Brown 
147286994e45SJed Brown    Level: advanced
147386994e45SJed Brown 
1474*95452b02SPatrick Sanan    Notes:
1475*95452b02SPatrick Sanan     ISLocalToGlobalMappingGetSize() returns the length the this array
1476107e9a97SBarry Smith 
1477107e9a97SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreIndices(), ISLocalToGlobalMappingGetBlockIndices(), ISLocalToGlobalMappingRestoreBlockIndices()
147886994e45SJed Brown @*/
14797087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
148086994e45SJed Brown {
148186994e45SJed Brown   PetscFunctionBegin;
148286994e45SJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
148386994e45SJed Brown   PetscValidPointer(array,2);
148445b6f7e9SBarry Smith   if (ltog->bs == 1) {
148586994e45SJed Brown     *array = ltog->indices;
148645b6f7e9SBarry Smith   } else {
148745b6f7e9SBarry Smith     PetscInt       *jj,k,i,j,n = ltog->n, bs = ltog->bs;
148845b6f7e9SBarry Smith     const PetscInt *ii;
148945b6f7e9SBarry Smith     PetscErrorCode ierr;
149045b6f7e9SBarry Smith 
149145b6f7e9SBarry Smith     ierr = PetscMalloc1(bs*n,&jj);CHKERRQ(ierr);
149245b6f7e9SBarry Smith     *array = jj;
149345b6f7e9SBarry Smith     k    = 0;
149445b6f7e9SBarry Smith     ii   = ltog->indices;
149545b6f7e9SBarry Smith     for (i=0; i<n; i++)
149645b6f7e9SBarry Smith       for (j=0; j<bs; j++)
149745b6f7e9SBarry Smith         jj[k++] = bs*ii[i] + j;
149845b6f7e9SBarry Smith   }
149986994e45SJed Brown   PetscFunctionReturn(0);
150086994e45SJed Brown }
150186994e45SJed Brown 
150286994e45SJed Brown /*@C
1503193a2b41SJulian Andrej    ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with ISLocalToGlobalMappingGetIndices()
150486994e45SJed Brown 
150586994e45SJed Brown    Not Collective
150686994e45SJed Brown 
150786994e45SJed Brown    Input Arguments:
150886994e45SJed Brown + ltog - local to global mapping
150986994e45SJed Brown - array - array of indices
151086994e45SJed Brown 
151186994e45SJed Brown    Level: advanced
151286994e45SJed Brown 
151386994e45SJed Brown .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices()
151486994e45SJed Brown @*/
15157087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
151686994e45SJed Brown {
151786994e45SJed Brown   PetscFunctionBegin;
151886994e45SJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
151986994e45SJed Brown   PetscValidPointer(array,2);
152045b6f7e9SBarry Smith   if (ltog->bs == 1 && *array != ltog->indices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
152145b6f7e9SBarry Smith 
152245b6f7e9SBarry Smith   if (ltog->bs > 1) {
152345b6f7e9SBarry Smith     PetscErrorCode ierr;
152445b6f7e9SBarry Smith     ierr = PetscFree(*(void**)array);CHKERRQ(ierr);
152545b6f7e9SBarry Smith   }
152645b6f7e9SBarry Smith   PetscFunctionReturn(0);
152745b6f7e9SBarry Smith }
152845b6f7e9SBarry Smith 
152945b6f7e9SBarry Smith /*@C
153045b6f7e9SBarry Smith    ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block
153145b6f7e9SBarry Smith 
153245b6f7e9SBarry Smith    Not Collective
153345b6f7e9SBarry Smith 
153445b6f7e9SBarry Smith    Input Arguments:
153545b6f7e9SBarry Smith . ltog - local to global mapping
153645b6f7e9SBarry Smith 
153745b6f7e9SBarry Smith    Output Arguments:
153845b6f7e9SBarry Smith . array - array of indices
153945b6f7e9SBarry Smith 
154045b6f7e9SBarry Smith    Level: advanced
154145b6f7e9SBarry Smith 
154245b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreBlockIndices()
154345b6f7e9SBarry Smith @*/
154445b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
154545b6f7e9SBarry Smith {
154645b6f7e9SBarry Smith   PetscFunctionBegin;
154745b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
154845b6f7e9SBarry Smith   PetscValidPointer(array,2);
154945b6f7e9SBarry Smith   *array = ltog->indices;
155045b6f7e9SBarry Smith   PetscFunctionReturn(0);
155145b6f7e9SBarry Smith }
155245b6f7e9SBarry Smith 
155345b6f7e9SBarry Smith /*@C
155445b6f7e9SBarry Smith    ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with ISLocalToGlobalMappingGetBlockIndices()
155545b6f7e9SBarry Smith 
155645b6f7e9SBarry Smith    Not Collective
155745b6f7e9SBarry Smith 
155845b6f7e9SBarry Smith    Input Arguments:
155945b6f7e9SBarry Smith + ltog - local to global mapping
156045b6f7e9SBarry Smith - array - array of indices
156145b6f7e9SBarry Smith 
156245b6f7e9SBarry Smith    Level: advanced
156345b6f7e9SBarry Smith 
156445b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices()
156545b6f7e9SBarry Smith @*/
156645b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
156745b6f7e9SBarry Smith {
156845b6f7e9SBarry Smith   PetscFunctionBegin;
156945b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
157045b6f7e9SBarry Smith   PetscValidPointer(array,2);
157186994e45SJed Brown   if (*array != ltog->indices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
15720298fd71SBarry Smith   *array = NULL;
157386994e45SJed Brown   PetscFunctionReturn(0);
157486994e45SJed Brown }
1575f7efa3c7SJed Brown 
1576f7efa3c7SJed Brown /*@C
1577f7efa3c7SJed Brown    ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings
1578f7efa3c7SJed Brown 
1579f7efa3c7SJed Brown    Not Collective
1580f7efa3c7SJed Brown 
1581f7efa3c7SJed Brown    Input Arguments:
1582f7efa3c7SJed Brown + comm - communicator for the new mapping, must contain the communicator of every mapping to concatenate
1583f7efa3c7SJed Brown . n - number of mappings to concatenate
1584f7efa3c7SJed Brown - ltogs - local to global mappings
1585f7efa3c7SJed Brown 
1586f7efa3c7SJed Brown    Output Arguments:
1587f7efa3c7SJed Brown . ltogcat - new mapping
1588f7efa3c7SJed Brown 
15899d90f715SBarry Smith    Note: this currently always returns a mapping with block size of 1
15909d90f715SBarry Smith 
15919d90f715SBarry Smith    Developer Note: If all the input mapping have the same block size we could easily handle that as a special case
15929d90f715SBarry Smith 
1593f7efa3c7SJed Brown    Level: advanced
1594f7efa3c7SJed Brown 
1595f7efa3c7SJed Brown .seealso: ISLocalToGlobalMappingCreate()
1596f7efa3c7SJed Brown @*/
1597f7efa3c7SJed Brown PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm,PetscInt n,const ISLocalToGlobalMapping ltogs[],ISLocalToGlobalMapping *ltogcat)
1598f7efa3c7SJed Brown {
1599f7efa3c7SJed Brown   PetscInt       i,cnt,m,*idx;
1600f7efa3c7SJed Brown   PetscErrorCode ierr;
1601f7efa3c7SJed Brown 
1602f7efa3c7SJed Brown   PetscFunctionBegin;
1603f7efa3c7SJed Brown   if (n < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have a non-negative number of mappings, given %D",n);
1604f7efa3c7SJed Brown   if (n > 0) PetscValidPointer(ltogs,3);
1605f7efa3c7SJed Brown   for (i=0; i<n; i++) PetscValidHeaderSpecific(ltogs[i],IS_LTOGM_CLASSID,3);
1606f7efa3c7SJed Brown   PetscValidPointer(ltogcat,4);
1607f7efa3c7SJed Brown   for (cnt=0,i=0; i<n; i++) {
1608f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
1609f7efa3c7SJed Brown     cnt += m;
1610f7efa3c7SJed Brown   }
1611785e854fSJed Brown   ierr = PetscMalloc1(cnt,&idx);CHKERRQ(ierr);
1612f7efa3c7SJed Brown   for (cnt=0,i=0; i<n; i++) {
1613f7efa3c7SJed Brown     const PetscInt *subidx;
1614f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
1615f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetIndices(ltogs[i],&subidx);CHKERRQ(ierr);
1616f7efa3c7SJed Brown     ierr = PetscMemcpy(&idx[cnt],subidx,m*sizeof(PetscInt));CHKERRQ(ierr);
1617f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingRestoreIndices(ltogs[i],&subidx);CHKERRQ(ierr);
1618f7efa3c7SJed Brown     cnt += m;
1619f7efa3c7SJed Brown   }
1620f0413b6fSBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,1,cnt,idx,PETSC_OWN_POINTER,ltogcat);CHKERRQ(ierr);
1621f7efa3c7SJed Brown   PetscFunctionReturn(0);
1622f7efa3c7SJed Brown }
162304a59952SBarry Smith 
1624413f72f0SBarry Smith /*MC
1625413f72f0SBarry Smith       ISLOCALTOGLOBALMAPPINGBASIC - basic implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is
1626413f72f0SBarry Smith                                     used this is good for only small and moderate size problems.
1627413f72f0SBarry Smith 
1628413f72f0SBarry Smith    Options Database Keys:
1629413f72f0SBarry Smith +   -islocaltoglobalmapping_type basic - select this method
1630413f72f0SBarry Smith 
1631413f72f0SBarry Smith    Level: beginner
1632413f72f0SBarry Smith 
1633413f72f0SBarry Smith .seealso:  ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType(), ISLOCALTOGLOBALMAPPINGHASH
1634413f72f0SBarry Smith M*/
1635413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Basic(ISLocalToGlobalMapping ltog)
1636413f72f0SBarry Smith {
1637413f72f0SBarry Smith   PetscFunctionBegin;
1638413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Basic;
1639413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Basic;
1640413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Basic;
1641413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Basic;
1642413f72f0SBarry Smith   PetscFunctionReturn(0);
1643413f72f0SBarry Smith }
1644413f72f0SBarry Smith 
1645413f72f0SBarry Smith /*MC
1646413f72f0SBarry Smith       ISLOCALTOGLOBALMAPPINGHASH - hash implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is
1647ed56e8eaSBarry Smith                                     used this is good for large memory problems.
1648413f72f0SBarry Smith 
1649413f72f0SBarry Smith    Options Database Keys:
1650413f72f0SBarry Smith +   -islocaltoglobalmapping_type hash - select this method
1651413f72f0SBarry Smith 
1652ed56e8eaSBarry Smith 
1653*95452b02SPatrick Sanan    Notes:
1654*95452b02SPatrick Sanan     This is selected automatically for large problems if the user does not set the type.
1655ed56e8eaSBarry Smith 
1656413f72f0SBarry Smith    Level: beginner
1657413f72f0SBarry Smith 
1658413f72f0SBarry Smith .seealso:  ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType(), ISLOCALTOGLOBALMAPPINGHASH
1659413f72f0SBarry Smith M*/
1660413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Hash(ISLocalToGlobalMapping ltog)
1661413f72f0SBarry Smith {
1662413f72f0SBarry Smith   PetscFunctionBegin;
1663413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Hash;
1664413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Hash;
1665413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Hash;
1666413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Hash;
1667413f72f0SBarry Smith   PetscFunctionReturn(0);
1668413f72f0SBarry Smith }
1669413f72f0SBarry Smith 
1670413f72f0SBarry Smith 
1671413f72f0SBarry Smith /*@C
1672413f72f0SBarry Smith     ISLocalToGlobalMappingRegister -  Adds a method for applying a global to local mapping with an ISLocalToGlobalMapping
1673413f72f0SBarry Smith 
1674413f72f0SBarry Smith    Not Collective
1675413f72f0SBarry Smith 
1676413f72f0SBarry Smith    Input Parameters:
1677413f72f0SBarry Smith +  sname - name of a new method
1678413f72f0SBarry Smith -  routine_create - routine to create method context
1679413f72f0SBarry Smith 
1680413f72f0SBarry Smith    Notes:
1681ed56e8eaSBarry Smith    ISLocalToGlobalMappingRegister() may be called multiple times to add several user-defined mappings.
1682413f72f0SBarry Smith 
1683413f72f0SBarry Smith    Sample usage:
1684413f72f0SBarry Smith .vb
1685413f72f0SBarry Smith    ISLocalToGlobalMappingRegister("my_mapper",MyCreate);
1686413f72f0SBarry Smith .ve
1687413f72f0SBarry Smith 
1688ed56e8eaSBarry Smith    Then, your mapping can be chosen with the procedural interface via
1689413f72f0SBarry Smith $     ISLocalToGlobalMappingSetType(ltog,"my_mapper")
1690413f72f0SBarry Smith    or at runtime via the option
1691ed56e8eaSBarry Smith $     -islocaltoglobalmapping_type my_mapper
1692413f72f0SBarry Smith 
1693413f72f0SBarry Smith    Level: advanced
1694413f72f0SBarry Smith 
1695413f72f0SBarry Smith .keywords: ISLocalToGlobalMappingType, register
1696413f72f0SBarry Smith 
1697413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingRegisterAll(), ISLocalToGlobalMappingRegisterDestroy(), ISLOCALTOGLOBALMAPPINGBASIC, ISLOCALTOGLOBALMAPPINGHASH
1698413f72f0SBarry Smith 
1699413f72f0SBarry Smith @*/
1700413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRegister(const char sname[],PetscErrorCode (*function)(ISLocalToGlobalMapping))
1701413f72f0SBarry Smith {
1702413f72f0SBarry Smith   PetscErrorCode ierr;
1703413f72f0SBarry Smith 
1704413f72f0SBarry Smith   PetscFunctionBegin;
1705413f72f0SBarry Smith   ierr = PetscFunctionListAdd(&ISLocalToGlobalMappingList,sname,function);CHKERRQ(ierr);
1706413f72f0SBarry Smith   PetscFunctionReturn(0);
1707413f72f0SBarry Smith }
1708413f72f0SBarry Smith 
1709413f72f0SBarry Smith /*@C
1710ed56e8eaSBarry Smith    ISLocalToGlobalMappingSetType - Builds ISLocalToGlobalMapping for a particular global to local mapping approach.
1711413f72f0SBarry Smith 
1712413f72f0SBarry Smith    Logically Collective on ISLocalToGlobalMapping
1713413f72f0SBarry Smith 
1714413f72f0SBarry Smith    Input Parameters:
1715413f72f0SBarry Smith +  ltog - the ISLocalToGlobalMapping object
1716413f72f0SBarry Smith -  type - a known method
1717413f72f0SBarry Smith 
1718413f72f0SBarry Smith    Options Database Key:
1719ed56e8eaSBarry Smith .  -islocaltoglobalmapping_type  <method> - Sets the method; use -help for a list
1720413f72f0SBarry Smith     of available methods (for instance, basic or hash)
1721413f72f0SBarry Smith 
1722413f72f0SBarry Smith    Notes:
1723413f72f0SBarry Smith    See "petsc/include/petscis.h" for available methods
1724413f72f0SBarry Smith 
1725413f72f0SBarry Smith   Normally, it is best to use the ISLocalToGlobalMappingSetFromOptions() command and
1726413f72f0SBarry Smith   then set the ISLocalToGlobalMapping type from the options database rather than by using
1727413f72f0SBarry Smith   this routine.
1728413f72f0SBarry Smith 
1729413f72f0SBarry Smith   Level: intermediate
1730413f72f0SBarry Smith 
1731413f72f0SBarry Smith   Developer Note: ISLocalToGlobalMappingRegister() is used to add new types to ISLocalToGlobalMappingList from which they
1732413f72f0SBarry Smith   are accessed by ISLocalToGlobalMappingSetType().
1733413f72f0SBarry Smith 
1734413f72f0SBarry Smith .keywords: ISLocalToGlobalMapping, set, method
1735413f72f0SBarry Smith 
1736413f72f0SBarry Smith .seealso: PCSetType(), ISLocalToGlobalMappingType, ISLocalToGlobalMappingRegister(), ISLocalToGlobalMappingCreate()
1737413f72f0SBarry Smith 
1738413f72f0SBarry Smith @*/
1739413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType type)
1740413f72f0SBarry Smith {
1741413f72f0SBarry Smith   PetscErrorCode ierr,(*r)(ISLocalToGlobalMapping);
1742413f72f0SBarry Smith   PetscBool      match;
1743413f72f0SBarry Smith 
1744413f72f0SBarry Smith   PetscFunctionBegin;
1745413f72f0SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
1746413f72f0SBarry Smith   PetscValidCharPointer(type,2);
1747413f72f0SBarry Smith 
1748413f72f0SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)ltog,type,&match);CHKERRQ(ierr);
1749413f72f0SBarry Smith   if (match) PetscFunctionReturn(0);
1750413f72f0SBarry Smith 
1751413f72f0SBarry Smith   ierr =  PetscFunctionListFind(ISLocalToGlobalMappingList,type,&r);CHKERRQ(ierr);
1752413f72f0SBarry Smith   if (!r) SETERRQ1(PetscObjectComm((PetscObject)ltog),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested ISLocalToGlobalMapping type %s",type);
1753413f72f0SBarry Smith   /* Destroy the previous private LTOG context */
1754413f72f0SBarry Smith   if (ltog->ops->destroy) {
1755413f72f0SBarry Smith     ierr              = (*ltog->ops->destroy)(ltog);CHKERRQ(ierr);
1756413f72f0SBarry Smith     ltog->ops->destroy = NULL;
1757413f72f0SBarry Smith   }
1758413f72f0SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)ltog,type);CHKERRQ(ierr);
1759413f72f0SBarry Smith   ierr = (*r)(ltog);CHKERRQ(ierr);
1760413f72f0SBarry Smith   PetscFunctionReturn(0);
1761413f72f0SBarry Smith }
1762413f72f0SBarry Smith 
1763413f72f0SBarry Smith PetscBool ISLocalToGlobalMappingRegisterAllCalled = PETSC_FALSE;
1764413f72f0SBarry Smith 
1765413f72f0SBarry Smith /*@C
1766413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAll - Registers all of the local to global mapping components in the IS package.
1767413f72f0SBarry Smith 
1768413f72f0SBarry Smith   Not Collective
1769413f72f0SBarry Smith 
1770413f72f0SBarry Smith   Level: advanced
1771413f72f0SBarry Smith 
1772413f72f0SBarry Smith .keywords: IS, register, all
1773413f72f0SBarry Smith .seealso:  ISRegister(),  ISLocalToGlobalRegister()
1774413f72f0SBarry Smith @*/
1775413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRegisterAll(void)
1776413f72f0SBarry Smith {
1777413f72f0SBarry Smith   PetscErrorCode ierr;
1778413f72f0SBarry Smith 
1779413f72f0SBarry Smith   PetscFunctionBegin;
1780413f72f0SBarry Smith   if (ISLocalToGlobalMappingRegisterAllCalled) PetscFunctionReturn(0);
1781413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAllCalled = PETSC_TRUE;
1782413f72f0SBarry Smith 
1783413f72f0SBarry Smith   ierr = ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGBASIC, ISLocalToGlobalMappingCreate_Basic);CHKERRQ(ierr);
1784413f72f0SBarry Smith   ierr = ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGHASH, ISLocalToGlobalMappingCreate_Hash);CHKERRQ(ierr);
1785413f72f0SBarry Smith   PetscFunctionReturn(0);
1786413f72f0SBarry Smith }
178704a59952SBarry Smith 
1788