xref: /petsc/src/vec/is/utils/isltog.c (revision 1bd0b88e85093c6e165f8be5c7376b198831a222)
12362add9SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/isimpl.h>    /*I "petscis.h"  I*/
3e8f14785SLisandro Dalcin #include <petsc/private/hashmapi.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 {
15e8f14785SLisandro Dalcin   PetscHMapI 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);
85e8f14785SLisandro Dalcin   ierr = PetscHMapICreate(&map->globalht);CHKERRQ(ierr);
86413f72f0SBarry Smith   for (i=0; i<n; i++ ) {
87413f72f0SBarry Smith     if (idx[i] < 0) continue;
88e8f14785SLisandro Dalcin     ierr = PetscHMapISet(map->globalht,idx[i],i);CHKERRQ(ierr);
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);
114e8f14785SLisandro Dalcin   ierr = PetscHMapIDestroy(&map->globalht);CHKERRQ(ierr);
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 {                         \
141e8f14785SLisandro Dalcin     (void)PetscHMapIGet(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 {                         \
150e8f14785SLisandro Dalcin     (void)PetscHMapIGet(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 
26295452b02SPatrick Sanan     Notes:
26395452b02SPatrick 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;
401*1bd0b88eSStefano Zampini 
402*1bd0b88eSStefano Zampini   /* reset the cached information */
403*1bd0b88eSStefano Zampini   ierr = PetscFree(mapping->info_procs);CHKERRQ(ierr);
404*1bd0b88eSStefano Zampini   ierr = PetscFree(mapping->info_numprocs);CHKERRQ(ierr);
405*1bd0b88eSStefano Zampini   if (mapping->info_indices) {
406*1bd0b88eSStefano Zampini     PetscInt i;
407*1bd0b88eSStefano Zampini 
408*1bd0b88eSStefano Zampini     ierr = PetscFree((mapping->info_indices)[0]);CHKERRQ(ierr);
409*1bd0b88eSStefano Zampini     for (i=1; i<mapping->info_nproc; i++) {
410*1bd0b88eSStefano Zampini       ierr = PetscFree(mapping->info_indices[i]);CHKERRQ(ierr);
411*1bd0b88eSStefano Zampini     }
412*1bd0b88eSStefano Zampini     ierr = PetscFree(mapping->info_indices);CHKERRQ(ierr);
413*1bd0b88eSStefano Zampini   }
414*1bd0b88eSStefano Zampini   mapping->info_cached = PETSC_FALSE;
415*1bd0b88eSStefano Zampini 
416413f72f0SBarry Smith   if (mapping->ops->destroy) {
417413f72f0SBarry Smith     ierr = (*mapping->ops->destroy)(mapping);CHKERRQ(ierr);
418413f72f0SBarry Smith   }
41963fa5c83Sstefano_zampini   PetscFunctionReturn(0);
42063fa5c83Sstefano_zampini }
42163fa5c83Sstefano_zampini 
42245b6f7e9SBarry Smith /*@
42345b6f7e9SBarry Smith     ISLocalToGlobalMappingGetBlockSize - Gets the blocksize of the mapping
42445b6f7e9SBarry Smith     ordering and a global parallel ordering.
42545b6f7e9SBarry Smith 
42645b6f7e9SBarry Smith     Not Collective
42745b6f7e9SBarry Smith 
42845b6f7e9SBarry Smith     Input Parameters:
42945b6f7e9SBarry Smith .   mapping - mapping data structure
43045b6f7e9SBarry Smith 
43145b6f7e9SBarry Smith     Output Parameter:
43245b6f7e9SBarry Smith .   bs - the blocksize
43345b6f7e9SBarry Smith 
43445b6f7e9SBarry Smith     Level: advanced
43545b6f7e9SBarry Smith 
43645b6f7e9SBarry Smith     Concepts: mapping^local to global
43745b6f7e9SBarry Smith 
43845b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS()
43945b6f7e9SBarry Smith @*/
44045b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping mapping,PetscInt *bs)
44145b6f7e9SBarry Smith {
44245b6f7e9SBarry Smith   PetscFunctionBegin;
443cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
44445b6f7e9SBarry Smith   *bs = mapping->bs;
44545b6f7e9SBarry Smith   PetscFunctionReturn(0);
44645b6f7e9SBarry Smith }
44745b6f7e9SBarry Smith 
448ba5bb76aSSatish Balay /*@
44990f02eecSBarry Smith     ISLocalToGlobalMappingCreate - Creates a mapping between a local (0 to n)
45090f02eecSBarry Smith     ordering and a global parallel ordering.
4512362add9SBarry Smith 
45289d82c54SBarry Smith     Not Collective, but communicator may have more than one process
453b9cd556bSLois Curfman McInnes 
4542362add9SBarry Smith     Input Parameters:
45589d82c54SBarry Smith +   comm - MPI communicator
456f0413b6fSBarry Smith .   bs - the block size
45728bc9809SBarry Smith .   n - the number of local elements divided by the block size, or equivalently the number of block indices
45828bc9809SBarry 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
459d5ad8652SBarry Smith -   mode - see PetscCopyMode
4602362add9SBarry Smith 
461a997ad1aSLois Curfman McInnes     Output Parameter:
46290f02eecSBarry Smith .   mapping - new mapping data structure
4632362add9SBarry Smith 
46495452b02SPatrick Sanan     Notes:
46595452b02SPatrick Sanan     There is one integer value in indices per block and it represents the actual indices bs*idx + j, where j=0,..,bs-1
466413f72f0SBarry Smith 
4679a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
468413f72f0SBarry 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.
469413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
470413f72f0SBarry Smith 
471a997ad1aSLois Curfman McInnes     Level: advanced
472a997ad1aSLois Curfman McInnes 
473273d9f13SBarry Smith     Concepts: mapping^local to global
4742362add9SBarry Smith 
475413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingSetFromOptions(), ISLOCALTOGLOBALMAPPINGBASIC, ISLOCALTOGLOBALMAPPINGHASH
476413f72f0SBarry Smith           ISLocalToGlobalMappingSetType(), ISLocalToGlobalMappingType
4772362add9SBarry Smith @*/
47860c7cefcSBarry Smith PetscErrorCode  ISLocalToGlobalMappingCreate(MPI_Comm comm,PetscInt bs,PetscInt n,const PetscInt indices[],PetscCopyMode mode,ISLocalToGlobalMapping *mapping)
4792362add9SBarry Smith {
4806849ba73SBarry Smith   PetscErrorCode ierr;
48132dcc486SBarry Smith   PetscInt       *in;
482b46b645bSBarry Smith 
483b46b645bSBarry Smith   PetscFunctionBegin;
48473911063SBarry Smith   if (n) PetscValidIntPointer(indices,3);
4854482741eSBarry Smith   PetscValidPointer(mapping,4);
486b46b645bSBarry Smith 
4870298fd71SBarry Smith   *mapping = NULL;
488607a6623SBarry Smith   ierr = ISInitializePackage();CHKERRQ(ierr);
4892362add9SBarry Smith 
49073107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(*mapping,IS_LTOGM_CLASSID,"ISLocalToGlobalMapping","Local to global mapping","IS",
49160c7cefcSBarry Smith                            comm,ISLocalToGlobalMappingDestroy,ISLocalToGlobalMappingView);CHKERRQ(ierr);
492d4bb536fSBarry Smith   (*mapping)->n             = n;
493f0413b6fSBarry Smith   (*mapping)->bs            = bs;
494268a049cSStefano Zampini   (*mapping)->info_cached   = PETSC_FALSE;
495268a049cSStefano Zampini   (*mapping)->info_free     = PETSC_FALSE;
496268a049cSStefano Zampini   (*mapping)->info_procs    = NULL;
497268a049cSStefano Zampini   (*mapping)->info_numprocs = NULL;
498268a049cSStefano Zampini   (*mapping)->info_indices  = NULL;
499*1bd0b88eSStefano Zampini   (*mapping)->info_nodec    = NULL;
500*1bd0b88eSStefano Zampini   (*mapping)->info_nodei    = NULL;
501413f72f0SBarry Smith 
502413f72f0SBarry Smith   (*mapping)->ops->globaltolocalmappingapply      = NULL;
503413f72f0SBarry Smith   (*mapping)->ops->globaltolocalmappingapplyblock = NULL;
504413f72f0SBarry Smith   (*mapping)->ops->destroy                        = NULL;
505d5ad8652SBarry Smith   if (mode == PETSC_COPY_VALUES) {
506785e854fSJed Brown     ierr = PetscMalloc1(n,&in);CHKERRQ(ierr);
507d5ad8652SBarry Smith     ierr = PetscMemcpy(in,indices,n*sizeof(PetscInt));CHKERRQ(ierr);
508d5ad8652SBarry Smith     (*mapping)->indices = in;
5096389a1a1SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt));CHKERRQ(ierr);
5106389a1a1SBarry Smith   } else if (mode == PETSC_OWN_POINTER) {
5116389a1a1SBarry Smith     (*mapping)->indices = (PetscInt*)indices;
5126389a1a1SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt));CHKERRQ(ierr);
5136389a1a1SBarry Smith   }
51460c7cefcSBarry Smith   else SETERRQ(comm,PETSC_ERR_SUP,"Cannot currently use PETSC_USE_POINTER");
5153a40ed3dSBarry Smith   PetscFunctionReturn(0);
5162362add9SBarry Smith }
5172362add9SBarry Smith 
518413f72f0SBarry Smith PetscFunctionList ISLocalToGlobalMappingList = NULL;
519413f72f0SBarry Smith 
520413f72f0SBarry Smith 
52190f02eecSBarry Smith /*@
5227e99dc12SLawrence Mitchell    ISLocalToGlobalMappingSetFromOptions - Set mapping options from the options database.
5237e99dc12SLawrence Mitchell 
5247e99dc12SLawrence Mitchell    Not collective
5257e99dc12SLawrence Mitchell 
5267e99dc12SLawrence Mitchell    Input Parameters:
5277e99dc12SLawrence Mitchell .  mapping - mapping data structure
5287e99dc12SLawrence Mitchell 
5297e99dc12SLawrence Mitchell    Level: advanced
5307e99dc12SLawrence Mitchell 
5317e99dc12SLawrence Mitchell @*/
5327e99dc12SLawrence Mitchell PetscErrorCode ISLocalToGlobalMappingSetFromOptions(ISLocalToGlobalMapping mapping)
5337e99dc12SLawrence Mitchell {
5347e99dc12SLawrence Mitchell   PetscErrorCode             ierr;
535413f72f0SBarry Smith   char                       type[256];
536413f72f0SBarry Smith   ISLocalToGlobalMappingType defaulttype = "Not set";
5377e99dc12SLawrence Mitchell   PetscBool                  flg;
5387e99dc12SLawrence Mitchell 
5397e99dc12SLawrence Mitchell   PetscFunctionBegin;
5407e99dc12SLawrence Mitchell   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
541413f72f0SBarry Smith   ierr = ISLocalToGlobalMappingRegisterAll();CHKERRQ(ierr);
5427e99dc12SLawrence Mitchell   ierr = PetscObjectOptionsBegin((PetscObject)mapping);CHKERRQ(ierr);
543413f72f0SBarry Smith   ierr = PetscOptionsFList("-islocaltoglobalmapping_type","ISLocalToGlobalMapping method","ISLocalToGlobalMappingSetType",ISLocalToGlobalMappingList,(char*)(((PetscObject)mapping)->type_name) ? ((PetscObject)mapping)->type_name : defaulttype,type,256,&flg);CHKERRQ(ierr);
544413f72f0SBarry Smith   if (flg) {
545413f72f0SBarry Smith     ierr = ISLocalToGlobalMappingSetType(mapping,type);CHKERRQ(ierr);
546413f72f0SBarry Smith   }
5477e99dc12SLawrence Mitchell   ierr = PetscOptionsEnd();CHKERRQ(ierr);
5487e99dc12SLawrence Mitchell   PetscFunctionReturn(0);
5497e99dc12SLawrence Mitchell }
5507e99dc12SLawrence Mitchell 
5517e99dc12SLawrence Mitchell /*@
55290f02eecSBarry Smith    ISLocalToGlobalMappingDestroy - Destroys a mapping between a local (0 to n)
55390f02eecSBarry Smith    ordering and a global parallel ordering.
55490f02eecSBarry Smith 
5550f5bd95cSBarry Smith    Note Collective
556b9cd556bSLois Curfman McInnes 
55790f02eecSBarry Smith    Input Parameters:
55890f02eecSBarry Smith .  mapping - mapping data structure
55990f02eecSBarry Smith 
560a997ad1aSLois Curfman McInnes    Level: advanced
561a997ad1aSLois Curfman McInnes 
5623acfe500SLois Curfman McInnes .seealso: ISLocalToGlobalMappingCreate()
56390f02eecSBarry Smith @*/
5646bf464f9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping *mapping)
56590f02eecSBarry Smith {
566dfbe8321SBarry Smith   PetscErrorCode ierr;
5675fd66863SKarl Rupp 
5683a40ed3dSBarry Smith   PetscFunctionBegin;
5696bf464f9SBarry Smith   if (!*mapping) PetscFunctionReturn(0);
5706bf464f9SBarry Smith   PetscValidHeaderSpecific((*mapping),IS_LTOGM_CLASSID,1);
571997056adSBarry Smith   if (--((PetscObject)(*mapping))->refct > 0) {*mapping = 0;PetscFunctionReturn(0);}
5726bf464f9SBarry Smith   ierr = PetscFree((*mapping)->indices);CHKERRQ(ierr);
573268a049cSStefano Zampini   ierr = PetscFree((*mapping)->info_procs);CHKERRQ(ierr);
574268a049cSStefano Zampini   ierr = PetscFree((*mapping)->info_numprocs);CHKERRQ(ierr);
575268a049cSStefano Zampini   if ((*mapping)->info_indices) {
576268a049cSStefano Zampini     PetscInt i;
577268a049cSStefano Zampini 
578268a049cSStefano Zampini     ierr = PetscFree(((*mapping)->info_indices)[0]);CHKERRQ(ierr);
579268a049cSStefano Zampini     for (i=1; i<(*mapping)->info_nproc; i++) {
580268a049cSStefano Zampini       ierr = PetscFree(((*mapping)->info_indices)[i]);CHKERRQ(ierr);
581268a049cSStefano Zampini     }
582268a049cSStefano Zampini     ierr = PetscFree((*mapping)->info_indices);CHKERRQ(ierr);
583268a049cSStefano Zampini   }
584*1bd0b88eSStefano Zampini   if ((*mapping)->info_nodei) {
585*1bd0b88eSStefano Zampini     ierr = PetscFree(((*mapping)->info_nodei)[0]);CHKERRQ(ierr);
586*1bd0b88eSStefano Zampini   }
587*1bd0b88eSStefano Zampini   ierr = PetscFree((*mapping)->info_nodei);CHKERRQ(ierr);
588*1bd0b88eSStefano Zampini   ierr = PetscFree((*mapping)->info_nodec);CHKERRQ(ierr);
589413f72f0SBarry Smith   if ((*mapping)->ops->destroy) {
590413f72f0SBarry Smith     ierr = (*(*mapping)->ops->destroy)(*mapping);CHKERRQ(ierr);
591413f72f0SBarry Smith   }
592d38fa0fbSBarry Smith   ierr     = PetscHeaderDestroy(mapping);CHKERRQ(ierr);
593992144d0SBarry Smith   *mapping = 0;
5943a40ed3dSBarry Smith   PetscFunctionReturn(0);
59590f02eecSBarry Smith }
59690f02eecSBarry Smith 
59790f02eecSBarry Smith /*@
5983acfe500SLois Curfman McInnes     ISLocalToGlobalMappingApplyIS - Creates from an IS in the local numbering
5993acfe500SLois Curfman McInnes     a new index set using the global numbering defined in an ISLocalToGlobalMapping
6003acfe500SLois Curfman McInnes     context.
60190f02eecSBarry Smith 
602b9cd556bSLois Curfman McInnes     Not collective
603b9cd556bSLois Curfman McInnes 
60490f02eecSBarry Smith     Input Parameters:
605b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
606b9cd556bSLois Curfman McInnes -   is - index set in local numbering
60790f02eecSBarry Smith 
60890f02eecSBarry Smith     Output Parameters:
60990f02eecSBarry Smith .   newis - index set in global numbering
61090f02eecSBarry Smith 
611a997ad1aSLois Curfman McInnes     Level: advanced
612a997ad1aSLois Curfman McInnes 
613273d9f13SBarry Smith     Concepts: mapping^local to global
6143acfe500SLois Curfman McInnes 
61590f02eecSBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),
616d4bb536fSBarry Smith           ISLocalToGlobalMappingDestroy(), ISGlobalToLocalMappingApply()
61790f02eecSBarry Smith @*/
6187087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping,IS is,IS *newis)
61990f02eecSBarry Smith {
6206849ba73SBarry Smith   PetscErrorCode ierr;
621e24637baSBarry Smith   PetscInt       n,*idxout;
6225d0c19d7SBarry Smith   const PetscInt *idxin;
6233a40ed3dSBarry Smith 
6243a40ed3dSBarry Smith   PetscFunctionBegin;
6250700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
6260700a824SBarry Smith   PetscValidHeaderSpecific(is,IS_CLASSID,2);
6274482741eSBarry Smith   PetscValidPointer(newis,3);
62890f02eecSBarry Smith 
6293b9aefa3SBarry Smith   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
63090f02eecSBarry Smith   ierr = ISGetIndices(is,&idxin);CHKERRQ(ierr);
631785e854fSJed Brown   ierr = PetscMalloc1(n,&idxout);CHKERRQ(ierr);
632e24637baSBarry Smith   ierr = ISLocalToGlobalMappingApply(mapping,n,idxin,idxout);CHKERRQ(ierr);
6333b9aefa3SBarry Smith   ierr = ISRestoreIndices(is,&idxin);CHKERRQ(ierr);
634543f3098SMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),n,idxout,PETSC_OWN_POINTER,newis);CHKERRQ(ierr);
6353a40ed3dSBarry Smith   PetscFunctionReturn(0);
63690f02eecSBarry Smith }
63790f02eecSBarry Smith 
638b89cb25eSSatish Balay /*@
6393acfe500SLois Curfman McInnes    ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering
6403acfe500SLois Curfman McInnes    and converts them to the global numbering.
64190f02eecSBarry Smith 
642b9cd556bSLois Curfman McInnes    Not collective
643b9cd556bSLois Curfman McInnes 
644bb25748dSBarry Smith    Input Parameters:
645b9cd556bSLois Curfman McInnes +  mapping - the local to global mapping context
646bb25748dSBarry Smith .  N - number of integers
647b9cd556bSLois Curfman McInnes -  in - input indices in local numbering
648bb25748dSBarry Smith 
649bb25748dSBarry Smith    Output Parameter:
650bb25748dSBarry Smith .  out - indices in global numbering
651bb25748dSBarry Smith 
652b9cd556bSLois Curfman McInnes    Notes:
653b9cd556bSLois Curfman McInnes    The in and out array parameters may be identical.
654d4bb536fSBarry Smith 
655a997ad1aSLois Curfman McInnes    Level: advanced
656a997ad1aSLois Curfman McInnes 
65745b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
6580752156aSBarry Smith           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
659d4bb536fSBarry Smith           AOPetscToApplication(), ISGlobalToLocalMappingApply()
660bb25748dSBarry Smith 
661273d9f13SBarry Smith     Concepts: mapping^local to global
662afcb2eb5SJed Brown @*/
663afcb2eb5SJed Brown PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
664afcb2eb5SJed Brown {
665cbc1caf0SMatthew G. Knepley   PetscInt i,bs,Nmax;
66645b6f7e9SBarry Smith 
66745b6f7e9SBarry Smith   PetscFunctionBegin;
668cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
669cbc1caf0SMatthew G. Knepley   bs   = mapping->bs;
670cbc1caf0SMatthew G. Knepley   Nmax = bs*mapping->n;
67145b6f7e9SBarry Smith   if (bs == 1) {
672cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
67345b6f7e9SBarry Smith     for (i=0; i<N; i++) {
67445b6f7e9SBarry Smith       if (in[i] < 0) {
67545b6f7e9SBarry Smith         out[i] = in[i];
67645b6f7e9SBarry Smith         continue;
67745b6f7e9SBarry Smith       }
678e24637baSBarry 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);
67945b6f7e9SBarry Smith       out[i] = idx[in[i]];
68045b6f7e9SBarry Smith     }
68145b6f7e9SBarry Smith   } else {
682cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
68345b6f7e9SBarry Smith     for (i=0; i<N; i++) {
68445b6f7e9SBarry Smith       if (in[i] < 0) {
68545b6f7e9SBarry Smith         out[i] = in[i];
68645b6f7e9SBarry Smith         continue;
68745b6f7e9SBarry Smith       }
688e24637baSBarry 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);
68945b6f7e9SBarry Smith       out[i] = idx[in[i]/bs]*bs + (in[i] % bs);
69045b6f7e9SBarry Smith     }
69145b6f7e9SBarry Smith   }
69245b6f7e9SBarry Smith   PetscFunctionReturn(0);
69345b6f7e9SBarry Smith }
69445b6f7e9SBarry Smith 
69545b6f7e9SBarry Smith /*@
6966006e8d2SBarry Smith    ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local block numbering and converts them to the global block numbering
69745b6f7e9SBarry Smith 
69845b6f7e9SBarry Smith    Not collective
69945b6f7e9SBarry Smith 
70045b6f7e9SBarry Smith    Input Parameters:
70145b6f7e9SBarry Smith +  mapping - the local to global mapping context
70245b6f7e9SBarry Smith .  N - number of integers
7036006e8d2SBarry Smith -  in - input indices in local block numbering
70445b6f7e9SBarry Smith 
70545b6f7e9SBarry Smith    Output Parameter:
7066006e8d2SBarry Smith .  out - indices in global block numbering
70745b6f7e9SBarry Smith 
70845b6f7e9SBarry Smith    Notes:
70945b6f7e9SBarry Smith    The in and out array parameters may be identical.
71045b6f7e9SBarry Smith 
7116006e8d2SBarry Smith    Example:
7126006e8d2SBarry 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
7136006e8d2SBarry Smith      (the first block) would produce 0 and the mapping applied to 1 (the second block) would produce 3.
7146006e8d2SBarry Smith 
71545b6f7e9SBarry Smith    Level: advanced
71645b6f7e9SBarry Smith 
71745b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
71845b6f7e9SBarry Smith           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
71945b6f7e9SBarry Smith           AOPetscToApplication(), ISGlobalToLocalMappingApply()
72045b6f7e9SBarry Smith 
72145b6f7e9SBarry Smith     Concepts: mapping^local to global
72245b6f7e9SBarry Smith @*/
72345b6f7e9SBarry Smith PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
72445b6f7e9SBarry Smith {
725cbc1caf0SMatthew G. Knepley 
726cbc1caf0SMatthew G. Knepley   PetscFunctionBegin;
727cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
728cbc1caf0SMatthew G. Knepley   {
729afcb2eb5SJed Brown     PetscInt i,Nmax = mapping->n;
730afcb2eb5SJed Brown     const PetscInt *idx = mapping->indices;
731d4bb536fSBarry Smith 
732afcb2eb5SJed Brown     for (i=0; i<N; i++) {
733afcb2eb5SJed Brown       if (in[i] < 0) {
734afcb2eb5SJed Brown         out[i] = in[i];
735afcb2eb5SJed Brown         continue;
736afcb2eb5SJed Brown       }
737e24637baSBarry 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);
738afcb2eb5SJed Brown       out[i] = idx[in[i]];
739afcb2eb5SJed Brown     }
740cbc1caf0SMatthew G. Knepley   }
741afcb2eb5SJed Brown   PetscFunctionReturn(0);
742afcb2eb5SJed Brown }
743d4bb536fSBarry Smith 
7447e99dc12SLawrence Mitchell /*@
745a997ad1aSLois Curfman McInnes     ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers
746a997ad1aSLois Curfman McInnes     specified with a global numbering.
747d4bb536fSBarry Smith 
748b9cd556bSLois Curfman McInnes     Not collective
749b9cd556bSLois Curfman McInnes 
750d4bb536fSBarry Smith     Input Parameters:
751b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
752d4bb536fSBarry Smith .   type - IS_GTOLM_MASK - replaces global indices with no local value with -1
753d4bb536fSBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
754d4bb536fSBarry Smith .   n - number of global indices to map
755b9cd556bSLois Curfman McInnes -   idx - global indices to map
756d4bb536fSBarry Smith 
757d4bb536fSBarry Smith     Output Parameters:
758b9cd556bSLois Curfman McInnes +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
759b9cd556bSLois Curfman McInnes -   idxout - local index of each global index, one must pass in an array long enough
760e182c471SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApply() with
7610298fd71SBarry Smith              idxout == NULL to determine the required length (returned in nout)
762e182c471SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApply()
763e182c471SBarry Smith              a second time to set the values.
764d4bb536fSBarry Smith 
765b9cd556bSLois Curfman McInnes     Notes:
7660298fd71SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
767d4bb536fSBarry Smith 
7689a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
769413f72f0SBarry 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.
770413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
7710f5bd95cSBarry Smith 
772a997ad1aSLois Curfman McInnes     Level: advanced
773a997ad1aSLois Curfman McInnes 
77432fd6b96SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
77532fd6b96SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
77632fd6b96SBarry Smith 
777273d9f13SBarry Smith     Concepts: mapping^global to local
778d4bb536fSBarry Smith 
7799d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),
780413f72f0SBarry Smith           ISLocalToGlobalMappingDestroy()
781d4bb536fSBarry Smith @*/
782413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
783d4bb536fSBarry Smith {
7849d90f715SBarry Smith   PetscErrorCode ierr;
7859d90f715SBarry Smith 
7869d90f715SBarry Smith   PetscFunctionBegin;
7879d90f715SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
788413f72f0SBarry Smith   if (!mapping->data) {
789413f72f0SBarry Smith     ierr = ISGlobalToLocalMappingSetUp(mapping);CHKERRQ(ierr);
7909d90f715SBarry Smith   }
791413f72f0SBarry Smith   ierr = (*mapping->ops->globaltolocalmappingapply)(mapping,type,n,idx,nout,idxout);CHKERRQ(ierr);
7929d90f715SBarry Smith   PetscFunctionReturn(0);
7939d90f715SBarry Smith }
7949d90f715SBarry Smith 
795d4fe737eSStefano Zampini /*@
796d4fe737eSStefano Zampini     ISGlobalToLocalMappingApplyIS - Creates from an IS in the global numbering
797d4fe737eSStefano Zampini     a new index set using the local numbering defined in an ISLocalToGlobalMapping
798d4fe737eSStefano Zampini     context.
799d4fe737eSStefano Zampini 
800d4fe737eSStefano Zampini     Not collective
801d4fe737eSStefano Zampini 
802d4fe737eSStefano Zampini     Input Parameters:
803d4fe737eSStefano Zampini +   mapping - mapping between local and global numbering
8042785b321SStefano Zampini .   type - IS_GTOLM_MASK - replaces global indices with no local value with -1
8052785b321SStefano Zampini            IS_GTOLM_DROP - drops the indices with no local value from the output list
806d4fe737eSStefano Zampini -   is - index set in global numbering
807d4fe737eSStefano Zampini 
808d4fe737eSStefano Zampini     Output Parameters:
809d4fe737eSStefano Zampini .   newis - index set in local numbering
810d4fe737eSStefano Zampini 
811d4fe737eSStefano Zampini     Level: advanced
812d4fe737eSStefano Zampini 
813d4fe737eSStefano Zampini     Concepts: mapping^local to global
814d4fe737eSStefano Zampini 
815d4fe737eSStefano Zampini .seealso: ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(),
816d4fe737eSStefano Zampini           ISLocalToGlobalMappingDestroy()
817d4fe737eSStefano Zampini @*/
818413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,IS is,IS *newis)
819d4fe737eSStefano Zampini {
820d4fe737eSStefano Zampini   PetscErrorCode ierr;
821d4fe737eSStefano Zampini   PetscInt       n,nout,*idxout;
822d4fe737eSStefano Zampini   const PetscInt *idxin;
823d4fe737eSStefano Zampini 
824d4fe737eSStefano Zampini   PetscFunctionBegin;
825d4fe737eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
826d4fe737eSStefano Zampini   PetscValidHeaderSpecific(is,IS_CLASSID,3);
827d4fe737eSStefano Zampini   PetscValidPointer(newis,4);
828d4fe737eSStefano Zampini 
829d4fe737eSStefano Zampini   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
830d4fe737eSStefano Zampini   ierr = ISGetIndices(is,&idxin);CHKERRQ(ierr);
831d4fe737eSStefano Zampini   if (type == IS_GTOLM_MASK) {
832d4fe737eSStefano Zampini     ierr = PetscMalloc1(n,&idxout);CHKERRQ(ierr);
833d4fe737eSStefano Zampini   } else {
834d4fe737eSStefano Zampini     ierr = ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,NULL);CHKERRQ(ierr);
835d4fe737eSStefano Zampini     ierr = PetscMalloc1(nout,&idxout);CHKERRQ(ierr);
836d4fe737eSStefano Zampini   }
837d4fe737eSStefano Zampini   ierr = ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,idxout);CHKERRQ(ierr);
838d4fe737eSStefano Zampini   ierr = ISRestoreIndices(is,&idxin);CHKERRQ(ierr);
839d4fe737eSStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),nout,idxout,PETSC_OWN_POINTER,newis);CHKERRQ(ierr);
840d4fe737eSStefano Zampini   PetscFunctionReturn(0);
841d4fe737eSStefano Zampini }
842d4fe737eSStefano Zampini 
8439d90f715SBarry Smith /*@
8449d90f715SBarry Smith     ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers
8459d90f715SBarry Smith     specified with a block global numbering.
8469d90f715SBarry Smith 
8479d90f715SBarry Smith     Not collective
8489d90f715SBarry Smith 
8499d90f715SBarry Smith     Input Parameters:
8509d90f715SBarry Smith +   mapping - mapping between local and global numbering
8519d90f715SBarry Smith .   type - IS_GTOLM_MASK - replaces global indices with no local value with -1
8529d90f715SBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
8539d90f715SBarry Smith .   n - number of global indices to map
8549d90f715SBarry Smith -   idx - global indices to map
8559d90f715SBarry Smith 
8569d90f715SBarry Smith     Output Parameters:
8579d90f715SBarry Smith +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
8589d90f715SBarry Smith -   idxout - local index of each global index, one must pass in an array long enough
8599d90f715SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApplyBlock() with
8609d90f715SBarry Smith              idxout == NULL to determine the required length (returned in nout)
8619d90f715SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApplyBlock()
8629d90f715SBarry Smith              a second time to set the values.
8639d90f715SBarry Smith 
8649d90f715SBarry Smith     Notes:
8659d90f715SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
8669d90f715SBarry Smith 
8679a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
868413f72f0SBarry 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.
869413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
8709d90f715SBarry Smith 
8719d90f715SBarry Smith     Level: advanced
8729d90f715SBarry Smith 
8739d90f715SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
8749d90f715SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
8759d90f715SBarry Smith 
8769d90f715SBarry Smith     Concepts: mapping^global to local
8779d90f715SBarry Smith 
8789d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(),
879413f72f0SBarry Smith           ISLocalToGlobalMappingDestroy()
8809d90f715SBarry Smith @*/
881413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,
8829d90f715SBarry Smith                                                  PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
8839d90f715SBarry Smith {
8846849ba73SBarry Smith   PetscErrorCode ierr;
885d4bb536fSBarry Smith 
8863a40ed3dSBarry Smith   PetscFunctionBegin;
8870700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
888413f72f0SBarry Smith   if (!mapping->data) {
889413f72f0SBarry Smith     ierr = ISGlobalToLocalMappingSetUp(mapping);CHKERRQ(ierr);
890d4bb536fSBarry Smith   }
891413f72f0SBarry Smith   ierr = (*mapping->ops->globaltolocalmappingapplyblock)(mapping,type,n,idx,nout,idxout);CHKERRQ(ierr);
8923a40ed3dSBarry Smith   PetscFunctionReturn(0);
893d4bb536fSBarry Smith }
89490f02eecSBarry Smith 
895413f72f0SBarry Smith 
89689d82c54SBarry Smith /*@C
8976a818285SBarry Smith     ISLocalToGlobalMappingGetBlockInfo - Gets the neighbor information for each processor and
89889d82c54SBarry Smith      each index shared by more than one processor
89989d82c54SBarry Smith 
90089d82c54SBarry Smith     Collective on ISLocalToGlobalMapping
90189d82c54SBarry Smith 
90289d82c54SBarry Smith     Input Parameters:
90389d82c54SBarry Smith .   mapping - the mapping from local to global indexing
90489d82c54SBarry Smith 
90589d82c54SBarry Smith     Output Parameter:
90689d82c54SBarry Smith +   nproc - number of processors that are connected to this one
90789d82c54SBarry Smith .   proc - neighboring processors
90807b52d57SBarry Smith .   numproc - number of indices for each subdomain (processor)
9093463a7baSJed Brown -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
91089d82c54SBarry Smith 
91189d82c54SBarry Smith     Level: advanced
91289d82c54SBarry Smith 
913273d9f13SBarry Smith     Concepts: mapping^local to global
91489d82c54SBarry Smith 
9152cfcea29SBarry Smith     Fortran Usage:
9162cfcea29SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
9172cfcea29SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
9182cfcea29SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
9192cfcea29SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
9202cfcea29SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
9212cfcea29SBarry Smith 
9222cfcea29SBarry Smith 
92307b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
92407b52d57SBarry Smith           ISLocalToGlobalMappingRestoreInfo()
92589d82c54SBarry Smith @*/
9266a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
92789d82c54SBarry Smith {
9286849ba73SBarry Smith   PetscErrorCode ierr;
929268a049cSStefano Zampini 
930268a049cSStefano Zampini   PetscFunctionBegin;
931268a049cSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
932268a049cSStefano Zampini   if (mapping->info_cached) {
933268a049cSStefano Zampini     *nproc    = mapping->info_nproc;
934268a049cSStefano Zampini     *procs    = mapping->info_procs;
935268a049cSStefano Zampini     *numprocs = mapping->info_numprocs;
936268a049cSStefano Zampini     *indices  = mapping->info_indices;
937268a049cSStefano Zampini   } else {
938268a049cSStefano Zampini     ierr = ISLocalToGlobalMappingGetBlockInfo_Private(mapping,nproc,procs,numprocs,indices);CHKERRQ(ierr);
939268a049cSStefano Zampini   }
940268a049cSStefano Zampini   PetscFunctionReturn(0);
941268a049cSStefano Zampini }
942268a049cSStefano Zampini 
943268a049cSStefano Zampini static PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
944268a049cSStefano Zampini {
945268a049cSStefano Zampini   PetscErrorCode ierr;
94697f1f81fSBarry Smith   PetscMPIInt    size,rank,tag1,tag2,tag3,*len,*source,imdex;
94732dcc486SBarry Smith   PetscInt       i,n = mapping->n,Ng,ng,max = 0,*lindices = mapping->indices;
94832dcc486SBarry Smith   PetscInt       *nprocs,*owner,nsends,*sends,j,*starts,nmax,nrecvs,*recvs,proc;
94997f1f81fSBarry Smith   PetscInt       cnt,scale,*ownedsenders,*nownedsenders,rstart,nowned;
95032dcc486SBarry Smith   PetscInt       node,nownedm,nt,*sends2,nsends2,*starts2,*lens2,*dest,nrecvs2,*starts3,*recvs2,k,*bprocs,*tmp;
95132dcc486SBarry Smith   PetscInt       first_procs,first_numprocs,*first_indices;
95289d82c54SBarry Smith   MPI_Request    *recv_waits,*send_waits;
95330dcb7c9SBarry Smith   MPI_Status     recv_status,*send_status,*recv_statuses;
954ce94432eSBarry Smith   MPI_Comm       comm;
955ace3abfcSBarry Smith   PetscBool      debug = PETSC_FALSE;
95689d82c54SBarry Smith 
95789d82c54SBarry Smith   PetscFunctionBegin;
958ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)mapping,&comm);CHKERRQ(ierr);
95924cf384cSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
96024cf384cSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
96124cf384cSBarry Smith   if (size == 1) {
96224cf384cSBarry Smith     *nproc         = 0;
9630298fd71SBarry Smith     *procs         = NULL;
96495dccacaSBarry Smith     ierr           = PetscNew(numprocs);CHKERRQ(ierr);
9651e2105dcSBarry Smith     (*numprocs)[0] = 0;
96695dccacaSBarry Smith     ierr           = PetscNew(indices);CHKERRQ(ierr);
9670298fd71SBarry Smith     (*indices)[0]  = NULL;
968268a049cSStefano Zampini     /* save info for reuse */
969268a049cSStefano Zampini     mapping->info_nproc = *nproc;
970268a049cSStefano Zampini     mapping->info_procs = *procs;
971268a049cSStefano Zampini     mapping->info_numprocs = *numprocs;
972268a049cSStefano Zampini     mapping->info_indices = *indices;
973268a049cSStefano Zampini     mapping->info_cached = PETSC_TRUE;
97424cf384cSBarry Smith     PetscFunctionReturn(0);
97524cf384cSBarry Smith   }
97624cf384cSBarry Smith 
977c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)mapping)->options,NULL,"-islocaltoglobalmappinggetinfo_debug",&debug,NULL);CHKERRQ(ierr);
97807b52d57SBarry Smith 
9793677ff5aSBarry Smith   /*
9806a818285SBarry Smith     Notes on ISLocalToGlobalMappingGetBlockInfo
9813677ff5aSBarry Smith 
9823677ff5aSBarry Smith     globally owned node - the nodes that have been assigned to this processor in global
9833677ff5aSBarry Smith            numbering, just for this routine.
9843677ff5aSBarry Smith 
9853677ff5aSBarry Smith     nontrivial globally owned node - node assigned to this processor that is on a subdomain
9863677ff5aSBarry Smith            boundary (i.e. is has more than one local owner)
9873677ff5aSBarry Smith 
9883677ff5aSBarry Smith     locally owned node - node that exists on this processors subdomain
9893677ff5aSBarry Smith 
9903677ff5aSBarry Smith     nontrivial locally owned node - node that is not in the interior (i.e. has more than one
9913677ff5aSBarry Smith            local subdomain
9923677ff5aSBarry Smith   */
99324cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag1);CHKERRQ(ierr);
99424cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag2);CHKERRQ(ierr);
99524cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag3);CHKERRQ(ierr);
99689d82c54SBarry Smith 
99789d82c54SBarry Smith   for (i=0; i<n; i++) {
99889d82c54SBarry Smith     if (lindices[i] > max) max = lindices[i];
99989d82c54SBarry Smith   }
1000b2566f29SBarry Smith   ierr   = MPIU_Allreduce(&max,&Ng,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr);
100178058e43SBarry Smith   Ng++;
100289d82c54SBarry Smith   ierr   = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
100389d82c54SBarry Smith   ierr   = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1004bc8ff85bSBarry Smith   scale  = Ng/size + 1;
1005a2e34c3dSBarry Smith   ng     = scale; if (rank == size-1) ng = Ng - scale*(size-1); ng = PetscMax(1,ng);
1006caba0dd0SBarry Smith   rstart = scale*rank;
100789d82c54SBarry Smith 
100889d82c54SBarry Smith   /* determine ownership ranges of global indices */
1009785e854fSJed Brown   ierr = PetscMalloc1(2*size,&nprocs);CHKERRQ(ierr);
101032dcc486SBarry Smith   ierr = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr);
101189d82c54SBarry Smith 
101289d82c54SBarry Smith   /* determine owners of each local node  */
1013785e854fSJed Brown   ierr = PetscMalloc1(n,&owner);CHKERRQ(ierr);
101489d82c54SBarry Smith   for (i=0; i<n; i++) {
10153677ff5aSBarry Smith     proc             = lindices[i]/scale; /* processor that globally owns this index */
101627c402fcSBarry Smith     nprocs[2*proc+1] = 1;                 /* processor globally owns at least one of ours */
10173677ff5aSBarry Smith     owner[i]         = proc;
101827c402fcSBarry Smith     nprocs[2*proc]++;                     /* count of how many that processor globally owns of ours */
101989d82c54SBarry Smith   }
102027c402fcSBarry Smith   nsends = 0; for (i=0; i<size; i++) nsends += nprocs[2*i+1];
10217904a332SBarry Smith   ierr = PetscInfo1(mapping,"Number of global owners for my local data %D\n",nsends);CHKERRQ(ierr);
102289d82c54SBarry Smith 
102389d82c54SBarry Smith   /* inform other processors of number of messages and max length*/
102427c402fcSBarry Smith   ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr);
10257904a332SBarry Smith   ierr = PetscInfo1(mapping,"Number of local owners for my global data %D\n",nrecvs);CHKERRQ(ierr);
102689d82c54SBarry Smith 
102789d82c54SBarry Smith   /* post receives for owned rows */
1028785e854fSJed Brown   ierr = PetscMalloc1((2*nrecvs+1)*(nmax+1),&recvs);CHKERRQ(ierr);
1029854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&recv_waits);CHKERRQ(ierr);
103089d82c54SBarry Smith   for (i=0; i<nrecvs; i++) {
103132dcc486SBarry Smith     ierr = MPI_Irecv(recvs+2*nmax*i,2*nmax,MPIU_INT,MPI_ANY_SOURCE,tag1,comm,recv_waits+i);CHKERRQ(ierr);
103289d82c54SBarry Smith   }
103389d82c54SBarry Smith 
103489d82c54SBarry Smith   /* pack messages containing lists of local nodes to owners */
1035854ce69bSBarry Smith   ierr      = PetscMalloc1(2*n+1,&sends);CHKERRQ(ierr);
1036854ce69bSBarry Smith   ierr      = PetscMalloc1(size+1,&starts);CHKERRQ(ierr);
103789d82c54SBarry Smith   starts[0] = 0;
1038f6e5521dSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
103989d82c54SBarry Smith   for (i=0; i<n; i++) {
104089d82c54SBarry Smith     sends[starts[owner[i]]++] = lindices[i];
104130dcb7c9SBarry Smith     sends[starts[owner[i]]++] = i;
104289d82c54SBarry Smith   }
104389d82c54SBarry Smith   ierr = PetscFree(owner);CHKERRQ(ierr);
104489d82c54SBarry Smith   starts[0] = 0;
1045f6e5521dSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
104689d82c54SBarry Smith 
104789d82c54SBarry Smith   /* send the messages */
1048854ce69bSBarry Smith   ierr = PetscMalloc1(nsends+1,&send_waits);CHKERRQ(ierr);
1049854ce69bSBarry Smith   ierr = PetscMalloc1(nsends+1,&dest);CHKERRQ(ierr);
105089d82c54SBarry Smith   cnt = 0;
105189d82c54SBarry Smith   for (i=0; i<size; i++) {
105227c402fcSBarry Smith     if (nprocs[2*i]) {
105332dcc486SBarry Smith       ierr      = MPI_Isend(sends+starts[i],2*nprocs[2*i],MPIU_INT,i,tag1,comm,send_waits+cnt);CHKERRQ(ierr);
105430dcb7c9SBarry Smith       dest[cnt] = i;
105589d82c54SBarry Smith       cnt++;
105689d82c54SBarry Smith     }
105789d82c54SBarry Smith   }
105889d82c54SBarry Smith   ierr = PetscFree(starts);CHKERRQ(ierr);
105989d82c54SBarry Smith 
106089d82c54SBarry Smith   /* wait on receives */
1061854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&source);CHKERRQ(ierr);
1062854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&len);CHKERRQ(ierr);
106389d82c54SBarry Smith   cnt  = nrecvs;
1064854ce69bSBarry Smith   ierr = PetscMalloc1(ng+1,&nownedsenders);CHKERRQ(ierr);
106532dcc486SBarry Smith   ierr = PetscMemzero(nownedsenders,ng*sizeof(PetscInt));CHKERRQ(ierr);
106689d82c54SBarry Smith   while (cnt) {
106789d82c54SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
106889d82c54SBarry Smith     /* unpack receives into our local space */
106932dcc486SBarry Smith     ierr          = MPI_Get_count(&recv_status,MPIU_INT,&len[imdex]);CHKERRQ(ierr);
107089d82c54SBarry Smith     source[imdex] = recv_status.MPI_SOURCE;
107130dcb7c9SBarry Smith     len[imdex]    = len[imdex]/2;
1072caba0dd0SBarry Smith     /* count how many local owners for each of my global owned indices */
107330dcb7c9SBarry Smith     for (i=0; i<len[imdex]; i++) nownedsenders[recvs[2*imdex*nmax+2*i]-rstart]++;
107489d82c54SBarry Smith     cnt--;
107589d82c54SBarry Smith   }
107689d82c54SBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
107789d82c54SBarry Smith 
107830dcb7c9SBarry Smith   /* count how many globally owned indices are on an edge multiplied by how many processors own them. */
1079bc8ff85bSBarry Smith   nowned  = 0;
1080bc8ff85bSBarry Smith   nownedm = 0;
1081bc8ff85bSBarry Smith   for (i=0; i<ng; i++) {
1082bc8ff85bSBarry Smith     if (nownedsenders[i] > 1) {nownedm += nownedsenders[i]; nowned++;}
1083bc8ff85bSBarry Smith   }
1084bc8ff85bSBarry Smith 
1085bc8ff85bSBarry Smith   /* create single array to contain rank of all local owners of each globally owned index */
1086854ce69bSBarry Smith   ierr      = PetscMalloc1(nownedm+1,&ownedsenders);CHKERRQ(ierr);
1087854ce69bSBarry Smith   ierr      = PetscMalloc1(ng+1,&starts);CHKERRQ(ierr);
1088bc8ff85bSBarry Smith   starts[0] = 0;
1089bc8ff85bSBarry Smith   for (i=1; i<ng; i++) {
1090bc8ff85bSBarry Smith     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
1091bc8ff85bSBarry Smith     else starts[i] = starts[i-1];
1092bc8ff85bSBarry Smith   }
1093bc8ff85bSBarry Smith 
109430dcb7c9SBarry Smith   /* for each nontrival globally owned node list all arriving processors */
1095bc8ff85bSBarry Smith   for (i=0; i<nrecvs; i++) {
1096bc8ff85bSBarry Smith     for (j=0; j<len[i]; j++) {
109730dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
1098f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) ownedsenders[starts[node]++] = source[i];
1099bc8ff85bSBarry Smith     }
1100bc8ff85bSBarry Smith   }
1101bc8ff85bSBarry Smith 
110207b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
110330dcb7c9SBarry Smith     starts[0] = 0;
110430dcb7c9SBarry Smith     for (i=1; i<ng; i++) {
110530dcb7c9SBarry Smith       if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
110630dcb7c9SBarry Smith       else starts[i] = starts[i-1];
110730dcb7c9SBarry Smith     }
110830dcb7c9SBarry Smith     for (i=0; i<ng; i++) {
110930dcb7c9SBarry Smith       if (nownedsenders[i] > 1) {
11107904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"[%d] global node %D local owner processors: ",rank,i+rstart);CHKERRQ(ierr);
111130dcb7c9SBarry Smith         for (j=0; j<nownedsenders[i]; j++) {
11127904a332SBarry Smith           ierr = PetscSynchronizedPrintf(comm,"%D ",ownedsenders[starts[i]+j]);CHKERRQ(ierr);
111330dcb7c9SBarry Smith         }
111430dcb7c9SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
111530dcb7c9SBarry Smith       }
111630dcb7c9SBarry Smith     }
11170ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
111807b52d57SBarry Smith   } /* -----------------------------------  */
111930dcb7c9SBarry Smith 
11203677ff5aSBarry Smith   /* wait on original sends */
11213a96401aSBarry Smith   if (nsends) {
1122785e854fSJed Brown     ierr = PetscMalloc1(nsends,&send_status);CHKERRQ(ierr);
11233a96401aSBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
11243a96401aSBarry Smith     ierr = PetscFree(send_status);CHKERRQ(ierr);
11253a96401aSBarry Smith   }
112689d82c54SBarry Smith   ierr = PetscFree(send_waits);CHKERRQ(ierr);
11273a96401aSBarry Smith   ierr = PetscFree(sends);CHKERRQ(ierr);
11283677ff5aSBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
11293677ff5aSBarry Smith 
11303677ff5aSBarry Smith   /* pack messages to send back to local owners */
113130dcb7c9SBarry Smith   starts[0] = 0;
113230dcb7c9SBarry Smith   for (i=1; i<ng; i++) {
113330dcb7c9SBarry Smith     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
113430dcb7c9SBarry Smith     else starts[i] = starts[i-1];
113530dcb7c9SBarry Smith   }
113630dcb7c9SBarry Smith   nsends2 = nrecvs;
1137854ce69bSBarry Smith   ierr    = PetscMalloc1(nsends2+1,&nprocs);CHKERRQ(ierr); /* length of each message */
113830dcb7c9SBarry Smith   for (i=0; i<nrecvs; i++) {
113930dcb7c9SBarry Smith     nprocs[i] = 1;
114030dcb7c9SBarry Smith     for (j=0; j<len[i]; j++) {
114130dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
1142f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) nprocs[i] += 2 + nownedsenders[node];
114330dcb7c9SBarry Smith     }
114430dcb7c9SBarry Smith   }
1145f6e5521dSKarl Rupp   nt = 0;
1146f6e5521dSKarl Rupp   for (i=0; i<nsends2; i++) nt += nprocs[i];
1147f6e5521dSKarl Rupp 
1148854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,&sends2);CHKERRQ(ierr);
1149854ce69bSBarry Smith   ierr = PetscMalloc1(nsends2+1,&starts2);CHKERRQ(ierr);
1150f6e5521dSKarl Rupp 
1151f6e5521dSKarl Rupp   starts2[0] = 0;
1152f6e5521dSKarl Rupp   for (i=1; i<nsends2; i++) starts2[i] = starts2[i-1] + nprocs[i-1];
115330dcb7c9SBarry Smith   /*
115430dcb7c9SBarry Smith      Each message is 1 + nprocs[i] long, and consists of
115530dcb7c9SBarry Smith        (0) the number of nodes being sent back
115630dcb7c9SBarry Smith        (1) the local node number,
115730dcb7c9SBarry Smith        (2) the number of processors sharing it,
115830dcb7c9SBarry Smith        (3) the processors sharing it
115930dcb7c9SBarry Smith   */
116030dcb7c9SBarry Smith   for (i=0; i<nsends2; i++) {
116130dcb7c9SBarry Smith     cnt = 1;
116230dcb7c9SBarry Smith     sends2[starts2[i]] = 0;
116330dcb7c9SBarry Smith     for (j=0; j<len[i]; j++) {
116430dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
116530dcb7c9SBarry Smith       if (nownedsenders[node] > 1) {
116630dcb7c9SBarry Smith         sends2[starts2[i]]++;
116730dcb7c9SBarry Smith         sends2[starts2[i]+cnt++] = recvs[2*i*nmax+2*j+1];
116830dcb7c9SBarry Smith         sends2[starts2[i]+cnt++] = nownedsenders[node];
116932dcc486SBarry Smith         ierr = PetscMemcpy(&sends2[starts2[i]+cnt],&ownedsenders[starts[node]],nownedsenders[node]*sizeof(PetscInt));CHKERRQ(ierr);
117030dcb7c9SBarry Smith         cnt += nownedsenders[node];
117130dcb7c9SBarry Smith       }
117230dcb7c9SBarry Smith     }
117330dcb7c9SBarry Smith   }
117430dcb7c9SBarry Smith 
117530dcb7c9SBarry Smith   /* receive the message lengths */
117630dcb7c9SBarry Smith   nrecvs2 = nsends;
1177854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&lens2);CHKERRQ(ierr);
1178854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&starts3);CHKERRQ(ierr);
1179854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&recv_waits);CHKERRQ(ierr);
118030dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
1181d44834fbSBarry Smith     ierr = MPI_Irecv(&lens2[i],1,MPIU_INT,dest[i],tag2,comm,recv_waits+i);CHKERRQ(ierr);
118230dcb7c9SBarry Smith   }
1183d44834fbSBarry Smith 
11848a8e0b3aSBarry Smith   /* send the message lengths */
11858a8e0b3aSBarry Smith   for (i=0; i<nsends2; i++) {
11868a8e0b3aSBarry Smith     ierr = MPI_Send(&nprocs[i],1,MPIU_INT,source[i],tag2,comm);CHKERRQ(ierr);
11878a8e0b3aSBarry Smith   }
11888a8e0b3aSBarry Smith 
1189d44834fbSBarry Smith   /* wait on receives of lens */
11900c468ba9SBarry Smith   if (nrecvs2) {
1191785e854fSJed Brown     ierr = PetscMalloc1(nrecvs2,&recv_statuses);CHKERRQ(ierr);
1192d44834fbSBarry Smith     ierr = MPI_Waitall(nrecvs2,recv_waits,recv_statuses);CHKERRQ(ierr);
1193d44834fbSBarry Smith     ierr = PetscFree(recv_statuses);CHKERRQ(ierr);
11940c468ba9SBarry Smith   }
1195a2ea699eSBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
1196d44834fbSBarry Smith 
119730dcb7c9SBarry Smith   starts3[0] = 0;
1198d44834fbSBarry Smith   nt         = 0;
119930dcb7c9SBarry Smith   for (i=0; i<nrecvs2-1; i++) {
120030dcb7c9SBarry Smith     starts3[i+1] = starts3[i] + lens2[i];
1201d44834fbSBarry Smith     nt          += lens2[i];
120230dcb7c9SBarry Smith   }
120376466f69SStefano Zampini   if (nrecvs2) nt += lens2[nrecvs2-1];
1204d44834fbSBarry Smith 
1205854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,&recvs2);CHKERRQ(ierr);
1206854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs2+1,&recv_waits);CHKERRQ(ierr);
120752b72c4aSBarry Smith   for (i=0; i<nrecvs2; i++) {
120832dcc486SBarry Smith     ierr = MPI_Irecv(recvs2+starts3[i],lens2[i],MPIU_INT,dest[i],tag3,comm,recv_waits+i);CHKERRQ(ierr);
120930dcb7c9SBarry Smith   }
121030dcb7c9SBarry Smith 
121130dcb7c9SBarry Smith   /* send the messages */
1212854ce69bSBarry Smith   ierr = PetscMalloc1(nsends2+1,&send_waits);CHKERRQ(ierr);
121330dcb7c9SBarry Smith   for (i=0; i<nsends2; i++) {
121432dcc486SBarry Smith     ierr = MPI_Isend(sends2+starts2[i],nprocs[i],MPIU_INT,source[i],tag3,comm,send_waits+i);CHKERRQ(ierr);
121530dcb7c9SBarry Smith   }
121630dcb7c9SBarry Smith 
121730dcb7c9SBarry Smith   /* wait on receives */
12180c468ba9SBarry Smith   if (nrecvs2) {
1219785e854fSJed Brown     ierr = PetscMalloc1(nrecvs2,&recv_statuses);CHKERRQ(ierr);
122030dcb7c9SBarry Smith     ierr = MPI_Waitall(nrecvs2,recv_waits,recv_statuses);CHKERRQ(ierr);
122130dcb7c9SBarry Smith     ierr = PetscFree(recv_statuses);CHKERRQ(ierr);
12220c468ba9SBarry Smith   }
122330dcb7c9SBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
122430dcb7c9SBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
122530dcb7c9SBarry Smith 
122607b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
122730dcb7c9SBarry Smith     cnt = 0;
122830dcb7c9SBarry Smith     for (i=0; i<nrecvs2; i++) {
122930dcb7c9SBarry Smith       nt = recvs2[cnt++];
123030dcb7c9SBarry Smith       for (j=0; j<nt; j++) {
12317904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"[%d] local node %D number of subdomains %D: ",rank,recvs2[cnt],recvs2[cnt+1]);CHKERRQ(ierr);
123230dcb7c9SBarry Smith         for (k=0; k<recvs2[cnt+1]; k++) {
12337904a332SBarry Smith           ierr = PetscSynchronizedPrintf(comm,"%D ",recvs2[cnt+2+k]);CHKERRQ(ierr);
123430dcb7c9SBarry Smith         }
123530dcb7c9SBarry Smith         cnt += 2 + recvs2[cnt+1];
123630dcb7c9SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
123730dcb7c9SBarry Smith       }
123830dcb7c9SBarry Smith     }
12390ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
124007b52d57SBarry Smith   } /* -----------------------------------  */
124130dcb7c9SBarry Smith 
124230dcb7c9SBarry Smith   /* count number subdomains for each local node */
1243785e854fSJed Brown   ierr = PetscMalloc1(size,&nprocs);CHKERRQ(ierr);
124432dcc486SBarry Smith   ierr = PetscMemzero(nprocs,size*sizeof(PetscInt));CHKERRQ(ierr);
124530dcb7c9SBarry Smith   cnt  = 0;
124630dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
124730dcb7c9SBarry Smith     nt = recvs2[cnt++];
124830dcb7c9SBarry Smith     for (j=0; j<nt; j++) {
1249f6e5521dSKarl Rupp       for (k=0; k<recvs2[cnt+1]; k++) nprocs[recvs2[cnt+2+k]]++;
125030dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt+1];
125130dcb7c9SBarry Smith     }
125230dcb7c9SBarry Smith   }
125330dcb7c9SBarry Smith   nt = 0; for (i=0; i<size; i++) nt += (nprocs[i] > 0);
125430dcb7c9SBarry Smith   *nproc    = nt;
1255854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,procs);CHKERRQ(ierr);
1256854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,numprocs);CHKERRQ(ierr);
1257854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,indices);CHKERRQ(ierr);
12580298fd71SBarry Smith   for (i=0;i<nt+1;i++) (*indices)[i]=NULL;
1259785e854fSJed Brown   ierr = PetscMalloc1(size,&bprocs);CHKERRQ(ierr);
126030dcb7c9SBarry Smith   cnt  = 0;
126130dcb7c9SBarry Smith   for (i=0; i<size; i++) {
126230dcb7c9SBarry Smith     if (nprocs[i] > 0) {
126330dcb7c9SBarry Smith       bprocs[i]        = cnt;
126430dcb7c9SBarry Smith       (*procs)[cnt]    = i;
126530dcb7c9SBarry Smith       (*numprocs)[cnt] = nprocs[i];
1266785e854fSJed Brown       ierr             = PetscMalloc1(nprocs[i],&(*indices)[cnt]);CHKERRQ(ierr);
126730dcb7c9SBarry Smith       cnt++;
126830dcb7c9SBarry Smith     }
126930dcb7c9SBarry Smith   }
127030dcb7c9SBarry Smith 
127130dcb7c9SBarry Smith   /* make the list of subdomains for each nontrivial local node */
127232dcc486SBarry Smith   ierr = PetscMemzero(*numprocs,nt*sizeof(PetscInt));CHKERRQ(ierr);
127330dcb7c9SBarry Smith   cnt  = 0;
127430dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
127530dcb7c9SBarry Smith     nt = recvs2[cnt++];
127630dcb7c9SBarry Smith     for (j=0; j<nt; j++) {
1277f6e5521dSKarl Rupp       for (k=0; k<recvs2[cnt+1]; k++) (*indices)[bprocs[recvs2[cnt+2+k]]][(*numprocs)[bprocs[recvs2[cnt+2+k]]]++] = recvs2[cnt];
127830dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt+1];
127930dcb7c9SBarry Smith     }
128030dcb7c9SBarry Smith   }
128130dcb7c9SBarry Smith   ierr = PetscFree(bprocs);CHKERRQ(ierr);
128207b52d57SBarry Smith   ierr = PetscFree(recvs2);CHKERRQ(ierr);
128330dcb7c9SBarry Smith 
128407b52d57SBarry Smith   /* sort the node indexing by their global numbers */
128507b52d57SBarry Smith   nt = *nproc;
128607b52d57SBarry Smith   for (i=0; i<nt; i++) {
1287854ce69bSBarry Smith     ierr = PetscMalloc1((*numprocs)[i],&tmp);CHKERRQ(ierr);
1288f6e5521dSKarl Rupp     for (j=0; j<(*numprocs)[i]; j++) tmp[j] = lindices[(*indices)[i][j]];
128907b52d57SBarry Smith     ierr = PetscSortIntWithArray((*numprocs)[i],tmp,(*indices)[i]);CHKERRQ(ierr);
129007b52d57SBarry Smith     ierr = PetscFree(tmp);CHKERRQ(ierr);
129107b52d57SBarry Smith   }
129207b52d57SBarry Smith 
129307b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
129430dcb7c9SBarry Smith     nt = *nproc;
129530dcb7c9SBarry Smith     for (i=0; i<nt; i++) {
12967904a332SBarry Smith       ierr = PetscSynchronizedPrintf(comm,"[%d] subdomain %D number of indices %D: ",rank,(*procs)[i],(*numprocs)[i]);CHKERRQ(ierr);
129730dcb7c9SBarry Smith       for (j=0; j<(*numprocs)[i]; j++) {
12987904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"%D ",(*indices)[i][j]);CHKERRQ(ierr);
129930dcb7c9SBarry Smith       }
130030dcb7c9SBarry Smith       ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
130130dcb7c9SBarry Smith     }
13020ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
130307b52d57SBarry Smith   } /* -----------------------------------  */
130430dcb7c9SBarry Smith 
130530dcb7c9SBarry Smith   /* wait on sends */
130630dcb7c9SBarry Smith   if (nsends2) {
1307785e854fSJed Brown     ierr = PetscMalloc1(nsends2,&send_status);CHKERRQ(ierr);
130830dcb7c9SBarry Smith     ierr = MPI_Waitall(nsends2,send_waits,send_status);CHKERRQ(ierr);
130930dcb7c9SBarry Smith     ierr = PetscFree(send_status);CHKERRQ(ierr);
131030dcb7c9SBarry Smith   }
131130dcb7c9SBarry Smith 
131230dcb7c9SBarry Smith   ierr = PetscFree(starts3);CHKERRQ(ierr);
131330dcb7c9SBarry Smith   ierr = PetscFree(dest);CHKERRQ(ierr);
131430dcb7c9SBarry Smith   ierr = PetscFree(send_waits);CHKERRQ(ierr);
13153677ff5aSBarry Smith 
1316bc8ff85bSBarry Smith   ierr = PetscFree(nownedsenders);CHKERRQ(ierr);
1317bc8ff85bSBarry Smith   ierr = PetscFree(ownedsenders);CHKERRQ(ierr);
1318bc8ff85bSBarry Smith   ierr = PetscFree(starts);CHKERRQ(ierr);
131930dcb7c9SBarry Smith   ierr = PetscFree(starts2);CHKERRQ(ierr);
132030dcb7c9SBarry Smith   ierr = PetscFree(lens2);CHKERRQ(ierr);
132189d82c54SBarry Smith 
132289d82c54SBarry Smith   ierr = PetscFree(source);CHKERRQ(ierr);
132397f1f81fSBarry Smith   ierr = PetscFree(len);CHKERRQ(ierr);
132489d82c54SBarry Smith   ierr = PetscFree(recvs);CHKERRQ(ierr);
13253a96401aSBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
132630dcb7c9SBarry Smith   ierr = PetscFree(sends2);CHKERRQ(ierr);
132724cf384cSBarry Smith 
132824cf384cSBarry Smith   /* put the information about myself as the first entry in the list */
132924cf384cSBarry Smith   first_procs    = (*procs)[0];
133024cf384cSBarry Smith   first_numprocs = (*numprocs)[0];
133124cf384cSBarry Smith   first_indices  = (*indices)[0];
133224cf384cSBarry Smith   for (i=0; i<*nproc; i++) {
133324cf384cSBarry Smith     if ((*procs)[i] == rank) {
133424cf384cSBarry Smith       (*procs)[0]    = (*procs)[i];
133524cf384cSBarry Smith       (*numprocs)[0] = (*numprocs)[i];
133624cf384cSBarry Smith       (*indices)[0]  = (*indices)[i];
133724cf384cSBarry Smith       (*procs)[i]    = first_procs;
133824cf384cSBarry Smith       (*numprocs)[i] = first_numprocs;
133924cf384cSBarry Smith       (*indices)[i]  = first_indices;
134024cf384cSBarry Smith       break;
134124cf384cSBarry Smith     }
134224cf384cSBarry Smith   }
1343268a049cSStefano Zampini 
1344268a049cSStefano Zampini   /* save info for reuse */
1345268a049cSStefano Zampini   mapping->info_nproc = *nproc;
1346268a049cSStefano Zampini   mapping->info_procs = *procs;
1347268a049cSStefano Zampini   mapping->info_numprocs = *numprocs;
1348268a049cSStefano Zampini   mapping->info_indices = *indices;
1349268a049cSStefano Zampini   mapping->info_cached = PETSC_TRUE;
135089d82c54SBarry Smith   PetscFunctionReturn(0);
135189d82c54SBarry Smith }
135289d82c54SBarry Smith 
13536a818285SBarry Smith /*@C
13546a818285SBarry Smith     ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by ISLocalToGlobalMappingGetBlockInfo()
13556a818285SBarry Smith 
13566a818285SBarry Smith     Collective on ISLocalToGlobalMapping
13576a818285SBarry Smith 
13586a818285SBarry Smith     Input Parameters:
13596a818285SBarry Smith .   mapping - the mapping from local to global indexing
13606a818285SBarry Smith 
13616a818285SBarry Smith     Output Parameter:
13626a818285SBarry Smith +   nproc - number of processors that are connected to this one
13636a818285SBarry Smith .   proc - neighboring processors
13646a818285SBarry Smith .   numproc - number of indices for each processor
13656a818285SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
13666a818285SBarry Smith 
13676a818285SBarry Smith     Level: advanced
13686a818285SBarry Smith 
13696a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
13706a818285SBarry Smith           ISLocalToGlobalMappingGetInfo()
13716a818285SBarry Smith @*/
13726a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
13736a818285SBarry Smith {
13746a818285SBarry Smith   PetscErrorCode ierr;
13756a818285SBarry Smith 
13766a818285SBarry Smith   PetscFunctionBegin;
1377cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1378268a049cSStefano Zampini   if (mapping->info_free) {
13796a818285SBarry Smith     ierr = PetscFree(*numprocs);CHKERRQ(ierr);
13806a818285SBarry Smith     if (*indices) {
1381268a049cSStefano Zampini       PetscInt i;
1382268a049cSStefano Zampini 
13836a818285SBarry Smith       ierr = PetscFree((*indices)[0]);CHKERRQ(ierr);
13846a818285SBarry Smith       for (i=1; i<*nproc; i++) {
13856a818285SBarry Smith         ierr = PetscFree((*indices)[i]);CHKERRQ(ierr);
13866a818285SBarry Smith       }
13876a818285SBarry Smith       ierr = PetscFree(*indices);CHKERRQ(ierr);
13886a818285SBarry Smith     }
1389268a049cSStefano Zampini   }
1390268a049cSStefano Zampini   *nproc    = 0;
1391268a049cSStefano Zampini   *procs    = NULL;
1392268a049cSStefano Zampini   *numprocs = NULL;
1393268a049cSStefano Zampini   *indices  = NULL;
13946a818285SBarry Smith   PetscFunctionReturn(0);
13956a818285SBarry Smith }
13966a818285SBarry Smith 
13976a818285SBarry Smith /*@C
13986a818285SBarry Smith     ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each processor and
13996a818285SBarry Smith      each index shared by more than one processor
14006a818285SBarry Smith 
14016a818285SBarry Smith     Collective on ISLocalToGlobalMapping
14026a818285SBarry Smith 
14036a818285SBarry Smith     Input Parameters:
14046a818285SBarry Smith .   mapping - the mapping from local to global indexing
14056a818285SBarry Smith 
14066a818285SBarry Smith     Output Parameter:
14076a818285SBarry Smith +   nproc - number of processors that are connected to this one
14086a818285SBarry Smith .   proc - neighboring processors
14096a818285SBarry Smith .   numproc - number of indices for each subdomain (processor)
14106a818285SBarry Smith -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
14116a818285SBarry Smith 
14126a818285SBarry Smith     Level: advanced
14136a818285SBarry Smith 
14146a818285SBarry Smith     Concepts: mapping^local to global
14156a818285SBarry Smith 
1416*1bd0b88eSStefano Zampini     Notes: The user needs to call ISLocalToGlobalMappingRestoreInfo when the data is no longer needed.
1417*1bd0b88eSStefano Zampini 
14186a818285SBarry Smith     Fortran Usage:
14196a818285SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
14206a818285SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
14216a818285SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
14226a818285SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
14236a818285SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
14246a818285SBarry Smith 
14256a818285SBarry Smith 
14266a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
14276a818285SBarry Smith           ISLocalToGlobalMappingRestoreInfo()
14286a818285SBarry Smith @*/
14296a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
14306a818285SBarry Smith {
14316a818285SBarry Smith   PetscErrorCode ierr;
1432268a049cSStefano Zampini   PetscInt       **bindices = NULL,*bnumprocs = NULL,bs = mapping->bs,i,j,k;
14336a818285SBarry Smith 
14346a818285SBarry Smith   PetscFunctionBegin;
14356a818285SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1436268a049cSStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockInfo(mapping,nproc,procs,&bnumprocs,&bindices);CHKERRQ(ierr);
1437268a049cSStefano Zampini   if (bs > 1) { /* we need to expand the cached info */
1438732f65e3SBarry Smith     ierr = PetscCalloc1(*nproc,&*indices);CHKERRQ(ierr);
1439268a049cSStefano Zampini     ierr = PetscCalloc1(*nproc,&*numprocs);CHKERRQ(ierr);
14406a818285SBarry Smith     for (i=0; i<*nproc; i++) {
1441268a049cSStefano Zampini       ierr = PetscMalloc1(bs*bnumprocs[i],&(*indices)[i]);CHKERRQ(ierr);
1442268a049cSStefano Zampini       for (j=0; j<bnumprocs[i]; j++) {
14436a818285SBarry Smith         for (k=0; k<bs; k++) {
14446a818285SBarry Smith           (*indices)[i][j*bs+k] = bs*bindices[i][j] + k;
14456a818285SBarry Smith         }
14466a818285SBarry Smith       }
1447268a049cSStefano Zampini       (*numprocs)[i] = bnumprocs[i]*bs;
14486a818285SBarry Smith     }
1449268a049cSStefano Zampini     mapping->info_free = PETSC_TRUE;
1450268a049cSStefano Zampini   } else {
1451268a049cSStefano Zampini     *numprocs = bnumprocs;
1452268a049cSStefano Zampini     *indices  = bindices;
14536a818285SBarry Smith   }
14546a818285SBarry Smith   PetscFunctionReturn(0);
14556a818285SBarry Smith }
14566a818285SBarry Smith 
145707b52d57SBarry Smith /*@C
145807b52d57SBarry Smith     ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by ISLocalToGlobalMappingGetInfo()
145989d82c54SBarry Smith 
146007b52d57SBarry Smith     Collective on ISLocalToGlobalMapping
146107b52d57SBarry Smith 
146207b52d57SBarry Smith     Input Parameters:
146307b52d57SBarry Smith .   mapping - the mapping from local to global indexing
146407b52d57SBarry Smith 
146507b52d57SBarry Smith     Output Parameter:
146607b52d57SBarry Smith +   nproc - number of processors that are connected to this one
146707b52d57SBarry Smith .   proc - neighboring processors
146807b52d57SBarry Smith .   numproc - number of indices for each processor
146907b52d57SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
147007b52d57SBarry Smith 
147107b52d57SBarry Smith     Level: advanced
147207b52d57SBarry Smith 
147307b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
147407b52d57SBarry Smith           ISLocalToGlobalMappingGetInfo()
147507b52d57SBarry Smith @*/
14767087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
147707b52d57SBarry Smith {
14786849ba73SBarry Smith   PetscErrorCode ierr;
147907b52d57SBarry Smith 
148007b52d57SBarry Smith   PetscFunctionBegin;
14816a818285SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockInfo(mapping,nproc,procs,numprocs,indices);CHKERRQ(ierr);
148207b52d57SBarry Smith   PetscFunctionReturn(0);
148307b52d57SBarry Smith }
148486994e45SJed Brown 
148586994e45SJed Brown /*@C
1486*1bd0b88eSStefano Zampini     ISLocalToGlobalMappingGetNodeInfo - Gets the neighbor information for each node
1487*1bd0b88eSStefano Zampini 
1488*1bd0b88eSStefano Zampini     Collective on ISLocalToGlobalMapping
1489*1bd0b88eSStefano Zampini 
1490*1bd0b88eSStefano Zampini     Input Parameters:
1491*1bd0b88eSStefano Zampini .   mapping - the mapping from local to global indexing
1492*1bd0b88eSStefano Zampini 
1493*1bd0b88eSStefano Zampini     Output Parameter:
1494*1bd0b88eSStefano Zampini +   nnodes - number of local nodes (same ISLocalToGlobalMappingGetSize())
1495*1bd0b88eSStefano Zampini .   count - number of neighboring processors per node
1496*1bd0b88eSStefano Zampini -   indices - indices of processes sharing the node (sorted)
1497*1bd0b88eSStefano Zampini 
1498*1bd0b88eSStefano Zampini     Level: advanced
1499*1bd0b88eSStefano Zampini 
1500*1bd0b88eSStefano Zampini     Concepts: mapping^local to global
1501*1bd0b88eSStefano Zampini 
1502*1bd0b88eSStefano Zampini     Notes: The user needs to call ISLocalToGlobalMappingRestoreInfo when the data is no longer needed.
1503*1bd0b88eSStefano Zampini 
1504*1bd0b88eSStefano Zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
1505*1bd0b88eSStefano Zampini           ISLocalToGlobalMappingGetInfo(), ISLocalToGlobalMappingRestoreNodeInfo()
1506*1bd0b88eSStefano Zampini @*/
1507*1bd0b88eSStefano Zampini PetscErrorCode  ISLocalToGlobalMappingGetNodeInfo(ISLocalToGlobalMapping mapping,PetscInt *nnodes,PetscInt *count[],PetscInt **indices[])
1508*1bd0b88eSStefano Zampini {
1509*1bd0b88eSStefano Zampini   PetscInt       n;
1510*1bd0b88eSStefano Zampini   PetscErrorCode ierr;
1511*1bd0b88eSStefano Zampini 
1512*1bd0b88eSStefano Zampini   PetscFunctionBegin;
1513*1bd0b88eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1514*1bd0b88eSStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(mapping,&n);CHKERRQ(ierr);
1515*1bd0b88eSStefano Zampini   if (!mapping->info_nodec) {
1516*1bd0b88eSStefano Zampini     PetscInt i,m,n_neigh,*neigh,*n_shared,**shared;
1517*1bd0b88eSStefano Zampini 
1518*1bd0b88eSStefano Zampini     ierr = PetscCalloc1(n+1,&mapping->info_nodec);CHKERRQ(ierr); /* always allocate to flag setup */
1519*1bd0b88eSStefano Zampini     ierr = PetscMalloc1(n,&mapping->info_nodei);CHKERRQ(ierr);
1520*1bd0b88eSStefano Zampini     ierr = ISLocalToGlobalMappingGetInfo(mapping,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
1521*1bd0b88eSStefano Zampini     for (i=0,m=0;i<n;i++) { mapping->info_nodec[i] = 1; m++; }
1522*1bd0b88eSStefano Zampini     for (i=1;i<n_neigh;i++) {
1523*1bd0b88eSStefano Zampini       PetscInt j;
1524*1bd0b88eSStefano Zampini 
1525*1bd0b88eSStefano Zampini       m += n_shared[i];
1526*1bd0b88eSStefano Zampini       for (j=0;j<n_shared[i];j++) mapping->info_nodec[shared[i][j]] += 1;
1527*1bd0b88eSStefano Zampini     }
1528*1bd0b88eSStefano Zampini     if (n) { ierr = PetscMalloc1(m,&mapping->info_nodei[0]);CHKERRQ(ierr); }
1529*1bd0b88eSStefano Zampini     for (i=1;i<n;i++) mapping->info_nodei[i] = mapping->info_nodei[i-1] + mapping->info_nodec[i-1];
1530*1bd0b88eSStefano Zampini     ierr = PetscMemzero(mapping->info_nodec,n*sizeof(PetscInt));CHKERRQ(ierr);
1531*1bd0b88eSStefano Zampini     for (i=0;i<n;i++) { mapping->info_nodec[i] = 1; mapping->info_nodei[i][0] = neigh[0]; }
1532*1bd0b88eSStefano Zampini     for (i=1;i<n_neigh;i++) {
1533*1bd0b88eSStefano Zampini       PetscInt j;
1534*1bd0b88eSStefano Zampini 
1535*1bd0b88eSStefano Zampini       for (j=0;j<n_shared[i];j++) {
1536*1bd0b88eSStefano Zampini         PetscInt k = shared[i][j];
1537*1bd0b88eSStefano Zampini 
1538*1bd0b88eSStefano Zampini         mapping->info_nodei[k][mapping->info_nodec[k]] = neigh[i];
1539*1bd0b88eSStefano Zampini         mapping->info_nodec[k] += 1;
1540*1bd0b88eSStefano Zampini       }
1541*1bd0b88eSStefano Zampini     }
1542*1bd0b88eSStefano Zampini     for (i=0;i<n;i++) { ierr = PetscSortRemoveDupsInt(&mapping->info_nodec[i],mapping->info_nodei[i]);CHKERRQ(ierr); }
1543*1bd0b88eSStefano Zampini     ierr = ISLocalToGlobalMappingRestoreInfo(mapping,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
1544*1bd0b88eSStefano Zampini   }
1545*1bd0b88eSStefano Zampini   if (nnodes)  *nnodes  = n;
1546*1bd0b88eSStefano Zampini   if (count)   *count   = mapping->info_nodec;
1547*1bd0b88eSStefano Zampini   if (indices) *indices = mapping->info_nodei;
1548*1bd0b88eSStefano Zampini   PetscFunctionReturn(0);
1549*1bd0b88eSStefano Zampini }
1550*1bd0b88eSStefano Zampini 
1551*1bd0b88eSStefano Zampini /*@C
1552*1bd0b88eSStefano Zampini     ISLocalToGlobalMappingRestoreNodeInfo - Frees the memory allocated by ISLocalToGlobalMappingGetNodeInfo()
1553*1bd0b88eSStefano Zampini 
1554*1bd0b88eSStefano Zampini     Collective on ISLocalToGlobalMapping
1555*1bd0b88eSStefano Zampini 
1556*1bd0b88eSStefano Zampini     Input Parameters:
1557*1bd0b88eSStefano Zampini .   mapping - the mapping from local to global indexing
1558*1bd0b88eSStefano Zampini 
1559*1bd0b88eSStefano Zampini     Output Parameter:
1560*1bd0b88eSStefano Zampini +   nnodes - number of local nodes
1561*1bd0b88eSStefano Zampini .   count - number of neighboring processors per node
1562*1bd0b88eSStefano Zampini -   indices - indices of processes sharing the node (sorted)
1563*1bd0b88eSStefano Zampini 
1564*1bd0b88eSStefano Zampini     Level: advanced
1565*1bd0b88eSStefano Zampini 
1566*1bd0b88eSStefano Zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
1567*1bd0b88eSStefano Zampini           ISLocalToGlobalMappingGetInfo()
1568*1bd0b88eSStefano Zampini @*/
1569*1bd0b88eSStefano Zampini PetscErrorCode  ISLocalToGlobalMappingRestoreNodeInfo(ISLocalToGlobalMapping mapping,PetscInt *nnodes,PetscInt *count[],PetscInt **indices[])
1570*1bd0b88eSStefano Zampini {
1571*1bd0b88eSStefano Zampini   PetscFunctionBegin;
1572*1bd0b88eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1573*1bd0b88eSStefano Zampini   if (nnodes)  *nnodes  = 0;
1574*1bd0b88eSStefano Zampini   if (count)   *count   = NULL;
1575*1bd0b88eSStefano Zampini   if (indices) *indices = NULL;
1576*1bd0b88eSStefano Zampini   PetscFunctionReturn(0);
1577*1bd0b88eSStefano Zampini }
1578*1bd0b88eSStefano Zampini 
1579*1bd0b88eSStefano Zampini /*@C
1580107e9a97SBarry Smith    ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped
158186994e45SJed Brown 
158286994e45SJed Brown    Not Collective
158386994e45SJed Brown 
158486994e45SJed Brown    Input Arguments:
158586994e45SJed Brown . ltog - local to global mapping
158686994e45SJed Brown 
158786994e45SJed Brown    Output Arguments:
1588565245c5SBarry Smith . array - array of indices, the length of this array may be obtained with ISLocalToGlobalMappingGetSize()
158986994e45SJed Brown 
159086994e45SJed Brown    Level: advanced
159186994e45SJed Brown 
159295452b02SPatrick Sanan    Notes:
159395452b02SPatrick Sanan     ISLocalToGlobalMappingGetSize() returns the length the this array
1594107e9a97SBarry Smith 
1595107e9a97SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreIndices(), ISLocalToGlobalMappingGetBlockIndices(), ISLocalToGlobalMappingRestoreBlockIndices()
159686994e45SJed Brown @*/
15977087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
159886994e45SJed Brown {
159986994e45SJed Brown   PetscFunctionBegin;
160086994e45SJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
160186994e45SJed Brown   PetscValidPointer(array,2);
160245b6f7e9SBarry Smith   if (ltog->bs == 1) {
160386994e45SJed Brown     *array = ltog->indices;
160445b6f7e9SBarry Smith   } else {
160545b6f7e9SBarry Smith     PetscInt       *jj,k,i,j,n = ltog->n, bs = ltog->bs;
160645b6f7e9SBarry Smith     const PetscInt *ii;
160745b6f7e9SBarry Smith     PetscErrorCode ierr;
160845b6f7e9SBarry Smith 
160945b6f7e9SBarry Smith     ierr = PetscMalloc1(bs*n,&jj);CHKERRQ(ierr);
161045b6f7e9SBarry Smith     *array = jj;
161145b6f7e9SBarry Smith     k    = 0;
161245b6f7e9SBarry Smith     ii   = ltog->indices;
161345b6f7e9SBarry Smith     for (i=0; i<n; i++)
161445b6f7e9SBarry Smith       for (j=0; j<bs; j++)
161545b6f7e9SBarry Smith         jj[k++] = bs*ii[i] + j;
161645b6f7e9SBarry Smith   }
161786994e45SJed Brown   PetscFunctionReturn(0);
161886994e45SJed Brown }
161986994e45SJed Brown 
162086994e45SJed Brown /*@C
1621193a2b41SJulian Andrej    ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with ISLocalToGlobalMappingGetIndices()
162286994e45SJed Brown 
162386994e45SJed Brown    Not Collective
162486994e45SJed Brown 
162586994e45SJed Brown    Input Arguments:
162686994e45SJed Brown + ltog - local to global mapping
162786994e45SJed Brown - array - array of indices
162886994e45SJed Brown 
162986994e45SJed Brown    Level: advanced
163086994e45SJed Brown 
163186994e45SJed Brown .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices()
163286994e45SJed Brown @*/
16337087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
163486994e45SJed Brown {
163586994e45SJed Brown   PetscFunctionBegin;
163686994e45SJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
163786994e45SJed Brown   PetscValidPointer(array,2);
163845b6f7e9SBarry Smith   if (ltog->bs == 1 && *array != ltog->indices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
163945b6f7e9SBarry Smith 
164045b6f7e9SBarry Smith   if (ltog->bs > 1) {
164145b6f7e9SBarry Smith     PetscErrorCode ierr;
164245b6f7e9SBarry Smith     ierr = PetscFree(*(void**)array);CHKERRQ(ierr);
164345b6f7e9SBarry Smith   }
164445b6f7e9SBarry Smith   PetscFunctionReturn(0);
164545b6f7e9SBarry Smith }
164645b6f7e9SBarry Smith 
164745b6f7e9SBarry Smith /*@C
164845b6f7e9SBarry Smith    ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block
164945b6f7e9SBarry Smith 
165045b6f7e9SBarry Smith    Not Collective
165145b6f7e9SBarry Smith 
165245b6f7e9SBarry Smith    Input Arguments:
165345b6f7e9SBarry Smith . ltog - local to global mapping
165445b6f7e9SBarry Smith 
165545b6f7e9SBarry Smith    Output Arguments:
165645b6f7e9SBarry Smith . array - array of indices
165745b6f7e9SBarry Smith 
165845b6f7e9SBarry Smith    Level: advanced
165945b6f7e9SBarry Smith 
166045b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreBlockIndices()
166145b6f7e9SBarry Smith @*/
166245b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
166345b6f7e9SBarry Smith {
166445b6f7e9SBarry Smith   PetscFunctionBegin;
166545b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
166645b6f7e9SBarry Smith   PetscValidPointer(array,2);
166745b6f7e9SBarry Smith   *array = ltog->indices;
166845b6f7e9SBarry Smith   PetscFunctionReturn(0);
166945b6f7e9SBarry Smith }
167045b6f7e9SBarry Smith 
167145b6f7e9SBarry Smith /*@C
167245b6f7e9SBarry Smith    ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with ISLocalToGlobalMappingGetBlockIndices()
167345b6f7e9SBarry Smith 
167445b6f7e9SBarry Smith    Not Collective
167545b6f7e9SBarry Smith 
167645b6f7e9SBarry Smith    Input Arguments:
167745b6f7e9SBarry Smith + ltog - local to global mapping
167845b6f7e9SBarry Smith - array - array of indices
167945b6f7e9SBarry Smith 
168045b6f7e9SBarry Smith    Level: advanced
168145b6f7e9SBarry Smith 
168245b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices()
168345b6f7e9SBarry Smith @*/
168445b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
168545b6f7e9SBarry Smith {
168645b6f7e9SBarry Smith   PetscFunctionBegin;
168745b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
168845b6f7e9SBarry Smith   PetscValidPointer(array,2);
168986994e45SJed Brown   if (*array != ltog->indices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
16900298fd71SBarry Smith   *array = NULL;
169186994e45SJed Brown   PetscFunctionReturn(0);
169286994e45SJed Brown }
1693f7efa3c7SJed Brown 
1694f7efa3c7SJed Brown /*@C
1695f7efa3c7SJed Brown    ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings
1696f7efa3c7SJed Brown 
1697f7efa3c7SJed Brown    Not Collective
1698f7efa3c7SJed Brown 
1699f7efa3c7SJed Brown    Input Arguments:
1700f7efa3c7SJed Brown + comm - communicator for the new mapping, must contain the communicator of every mapping to concatenate
1701f7efa3c7SJed Brown . n - number of mappings to concatenate
1702f7efa3c7SJed Brown - ltogs - local to global mappings
1703f7efa3c7SJed Brown 
1704f7efa3c7SJed Brown    Output Arguments:
1705f7efa3c7SJed Brown . ltogcat - new mapping
1706f7efa3c7SJed Brown 
17079d90f715SBarry Smith    Note: this currently always returns a mapping with block size of 1
17089d90f715SBarry Smith 
17099d90f715SBarry Smith    Developer Note: If all the input mapping have the same block size we could easily handle that as a special case
17109d90f715SBarry Smith 
1711f7efa3c7SJed Brown    Level: advanced
1712f7efa3c7SJed Brown 
1713f7efa3c7SJed Brown .seealso: ISLocalToGlobalMappingCreate()
1714f7efa3c7SJed Brown @*/
1715f7efa3c7SJed Brown PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm,PetscInt n,const ISLocalToGlobalMapping ltogs[],ISLocalToGlobalMapping *ltogcat)
1716f7efa3c7SJed Brown {
1717f7efa3c7SJed Brown   PetscInt       i,cnt,m,*idx;
1718f7efa3c7SJed Brown   PetscErrorCode ierr;
1719f7efa3c7SJed Brown 
1720f7efa3c7SJed Brown   PetscFunctionBegin;
1721f7efa3c7SJed Brown   if (n < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have a non-negative number of mappings, given %D",n);
1722f7efa3c7SJed Brown   if (n > 0) PetscValidPointer(ltogs,3);
1723f7efa3c7SJed Brown   for (i=0; i<n; i++) PetscValidHeaderSpecific(ltogs[i],IS_LTOGM_CLASSID,3);
1724f7efa3c7SJed Brown   PetscValidPointer(ltogcat,4);
1725f7efa3c7SJed Brown   for (cnt=0,i=0; i<n; i++) {
1726f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
1727f7efa3c7SJed Brown     cnt += m;
1728f7efa3c7SJed Brown   }
1729785e854fSJed Brown   ierr = PetscMalloc1(cnt,&idx);CHKERRQ(ierr);
1730f7efa3c7SJed Brown   for (cnt=0,i=0; i<n; i++) {
1731f7efa3c7SJed Brown     const PetscInt *subidx;
1732f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
1733f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetIndices(ltogs[i],&subidx);CHKERRQ(ierr);
1734f7efa3c7SJed Brown     ierr = PetscMemcpy(&idx[cnt],subidx,m*sizeof(PetscInt));CHKERRQ(ierr);
1735f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingRestoreIndices(ltogs[i],&subidx);CHKERRQ(ierr);
1736f7efa3c7SJed Brown     cnt += m;
1737f7efa3c7SJed Brown   }
1738f0413b6fSBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,1,cnt,idx,PETSC_OWN_POINTER,ltogcat);CHKERRQ(ierr);
1739f7efa3c7SJed Brown   PetscFunctionReturn(0);
1740f7efa3c7SJed Brown }
174104a59952SBarry Smith 
1742413f72f0SBarry Smith /*MC
1743413f72f0SBarry Smith       ISLOCALTOGLOBALMAPPINGBASIC - basic implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is
1744413f72f0SBarry Smith                                     used this is good for only small and moderate size problems.
1745413f72f0SBarry Smith 
1746413f72f0SBarry Smith    Options Database Keys:
1747413f72f0SBarry Smith +   -islocaltoglobalmapping_type basic - select this method
1748413f72f0SBarry Smith 
1749413f72f0SBarry Smith    Level: beginner
1750413f72f0SBarry Smith 
1751413f72f0SBarry Smith .seealso:  ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType(), ISLOCALTOGLOBALMAPPINGHASH
1752413f72f0SBarry Smith M*/
1753413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Basic(ISLocalToGlobalMapping ltog)
1754413f72f0SBarry Smith {
1755413f72f0SBarry Smith   PetscFunctionBegin;
1756413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Basic;
1757413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Basic;
1758413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Basic;
1759413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Basic;
1760413f72f0SBarry Smith   PetscFunctionReturn(0);
1761413f72f0SBarry Smith }
1762413f72f0SBarry Smith 
1763413f72f0SBarry Smith /*MC
1764413f72f0SBarry Smith       ISLOCALTOGLOBALMAPPINGHASH - hash implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is
1765ed56e8eaSBarry Smith                                     used this is good for large memory problems.
1766413f72f0SBarry Smith 
1767413f72f0SBarry Smith    Options Database Keys:
1768413f72f0SBarry Smith +   -islocaltoglobalmapping_type hash - select this method
1769413f72f0SBarry Smith 
1770ed56e8eaSBarry Smith 
177195452b02SPatrick Sanan    Notes:
177295452b02SPatrick Sanan     This is selected automatically for large problems if the user does not set the type.
1773ed56e8eaSBarry Smith 
1774413f72f0SBarry Smith    Level: beginner
1775413f72f0SBarry Smith 
1776413f72f0SBarry Smith .seealso:  ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType(), ISLOCALTOGLOBALMAPPINGHASH
1777413f72f0SBarry Smith M*/
1778413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Hash(ISLocalToGlobalMapping ltog)
1779413f72f0SBarry Smith {
1780413f72f0SBarry Smith   PetscFunctionBegin;
1781413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Hash;
1782413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Hash;
1783413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Hash;
1784413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Hash;
1785413f72f0SBarry Smith   PetscFunctionReturn(0);
1786413f72f0SBarry Smith }
1787413f72f0SBarry Smith 
1788413f72f0SBarry Smith 
1789413f72f0SBarry Smith /*@C
1790413f72f0SBarry Smith     ISLocalToGlobalMappingRegister -  Adds a method for applying a global to local mapping with an ISLocalToGlobalMapping
1791413f72f0SBarry Smith 
1792413f72f0SBarry Smith    Not Collective
1793413f72f0SBarry Smith 
1794413f72f0SBarry Smith    Input Parameters:
1795413f72f0SBarry Smith +  sname - name of a new method
1796413f72f0SBarry Smith -  routine_create - routine to create method context
1797413f72f0SBarry Smith 
1798413f72f0SBarry Smith    Notes:
1799ed56e8eaSBarry Smith    ISLocalToGlobalMappingRegister() may be called multiple times to add several user-defined mappings.
1800413f72f0SBarry Smith 
1801413f72f0SBarry Smith    Sample usage:
1802413f72f0SBarry Smith .vb
1803413f72f0SBarry Smith    ISLocalToGlobalMappingRegister("my_mapper",MyCreate);
1804413f72f0SBarry Smith .ve
1805413f72f0SBarry Smith 
1806ed56e8eaSBarry Smith    Then, your mapping can be chosen with the procedural interface via
1807413f72f0SBarry Smith $     ISLocalToGlobalMappingSetType(ltog,"my_mapper")
1808413f72f0SBarry Smith    or at runtime via the option
1809ed56e8eaSBarry Smith $     -islocaltoglobalmapping_type my_mapper
1810413f72f0SBarry Smith 
1811413f72f0SBarry Smith    Level: advanced
1812413f72f0SBarry Smith 
1813413f72f0SBarry Smith .keywords: ISLocalToGlobalMappingType, register
1814413f72f0SBarry Smith 
1815413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingRegisterAll(), ISLocalToGlobalMappingRegisterDestroy(), ISLOCALTOGLOBALMAPPINGBASIC, ISLOCALTOGLOBALMAPPINGHASH
1816413f72f0SBarry Smith 
1817413f72f0SBarry Smith @*/
1818413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRegister(const char sname[],PetscErrorCode (*function)(ISLocalToGlobalMapping))
1819413f72f0SBarry Smith {
1820413f72f0SBarry Smith   PetscErrorCode ierr;
1821413f72f0SBarry Smith 
1822413f72f0SBarry Smith   PetscFunctionBegin;
1823413f72f0SBarry Smith   ierr = PetscFunctionListAdd(&ISLocalToGlobalMappingList,sname,function);CHKERRQ(ierr);
1824413f72f0SBarry Smith   PetscFunctionReturn(0);
1825413f72f0SBarry Smith }
1826413f72f0SBarry Smith 
1827413f72f0SBarry Smith /*@C
1828ed56e8eaSBarry Smith    ISLocalToGlobalMappingSetType - Builds ISLocalToGlobalMapping for a particular global to local mapping approach.
1829413f72f0SBarry Smith 
1830413f72f0SBarry Smith    Logically Collective on ISLocalToGlobalMapping
1831413f72f0SBarry Smith 
1832413f72f0SBarry Smith    Input Parameters:
1833413f72f0SBarry Smith +  ltog - the ISLocalToGlobalMapping object
1834413f72f0SBarry Smith -  type - a known method
1835413f72f0SBarry Smith 
1836413f72f0SBarry Smith    Options Database Key:
1837ed56e8eaSBarry Smith .  -islocaltoglobalmapping_type  <method> - Sets the method; use -help for a list
1838413f72f0SBarry Smith     of available methods (for instance, basic or hash)
1839413f72f0SBarry Smith 
1840413f72f0SBarry Smith    Notes:
1841413f72f0SBarry Smith    See "petsc/include/petscis.h" for available methods
1842413f72f0SBarry Smith 
1843413f72f0SBarry Smith   Normally, it is best to use the ISLocalToGlobalMappingSetFromOptions() command and
1844413f72f0SBarry Smith   then set the ISLocalToGlobalMapping type from the options database rather than by using
1845413f72f0SBarry Smith   this routine.
1846413f72f0SBarry Smith 
1847413f72f0SBarry Smith   Level: intermediate
1848413f72f0SBarry Smith 
1849413f72f0SBarry Smith   Developer Note: ISLocalToGlobalMappingRegister() is used to add new types to ISLocalToGlobalMappingList from which they
1850413f72f0SBarry Smith   are accessed by ISLocalToGlobalMappingSetType().
1851413f72f0SBarry Smith 
1852413f72f0SBarry Smith .keywords: ISLocalToGlobalMapping, set, method
1853413f72f0SBarry Smith 
1854413f72f0SBarry Smith .seealso: PCSetType(), ISLocalToGlobalMappingType, ISLocalToGlobalMappingRegister(), ISLocalToGlobalMappingCreate()
1855413f72f0SBarry Smith 
1856413f72f0SBarry Smith @*/
1857413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType type)
1858413f72f0SBarry Smith {
1859413f72f0SBarry Smith   PetscErrorCode ierr,(*r)(ISLocalToGlobalMapping);
1860413f72f0SBarry Smith   PetscBool      match;
1861413f72f0SBarry Smith 
1862413f72f0SBarry Smith   PetscFunctionBegin;
1863413f72f0SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
1864413f72f0SBarry Smith   PetscValidCharPointer(type,2);
1865413f72f0SBarry Smith 
1866413f72f0SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)ltog,type,&match);CHKERRQ(ierr);
1867413f72f0SBarry Smith   if (match) PetscFunctionReturn(0);
1868413f72f0SBarry Smith 
1869413f72f0SBarry Smith   ierr =  PetscFunctionListFind(ISLocalToGlobalMappingList,type,&r);CHKERRQ(ierr);
1870413f72f0SBarry Smith   if (!r) SETERRQ1(PetscObjectComm((PetscObject)ltog),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested ISLocalToGlobalMapping type %s",type);
1871413f72f0SBarry Smith   /* Destroy the previous private LTOG context */
1872413f72f0SBarry Smith   if (ltog->ops->destroy) {
1873413f72f0SBarry Smith     ierr              = (*ltog->ops->destroy)(ltog);CHKERRQ(ierr);
1874413f72f0SBarry Smith     ltog->ops->destroy = NULL;
1875413f72f0SBarry Smith   }
1876413f72f0SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)ltog,type);CHKERRQ(ierr);
1877413f72f0SBarry Smith   ierr = (*r)(ltog);CHKERRQ(ierr);
1878413f72f0SBarry Smith   PetscFunctionReturn(0);
1879413f72f0SBarry Smith }
1880413f72f0SBarry Smith 
1881413f72f0SBarry Smith PetscBool ISLocalToGlobalMappingRegisterAllCalled = PETSC_FALSE;
1882413f72f0SBarry Smith 
1883413f72f0SBarry Smith /*@C
1884413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAll - Registers all of the local to global mapping components in the IS package.
1885413f72f0SBarry Smith 
1886413f72f0SBarry Smith   Not Collective
1887413f72f0SBarry Smith 
1888413f72f0SBarry Smith   Level: advanced
1889413f72f0SBarry Smith 
1890413f72f0SBarry Smith .keywords: IS, register, all
1891413f72f0SBarry Smith .seealso:  ISRegister(),  ISLocalToGlobalRegister()
1892413f72f0SBarry Smith @*/
1893413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRegisterAll(void)
1894413f72f0SBarry Smith {
1895413f72f0SBarry Smith   PetscErrorCode ierr;
1896413f72f0SBarry Smith 
1897413f72f0SBarry Smith   PetscFunctionBegin;
1898413f72f0SBarry Smith   if (ISLocalToGlobalMappingRegisterAllCalled) PetscFunctionReturn(0);
1899413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAllCalled = PETSC_TRUE;
1900413f72f0SBarry Smith 
1901413f72f0SBarry Smith   ierr = ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGBASIC, ISLocalToGlobalMappingCreate_Basic);CHKERRQ(ierr);
1902413f72f0SBarry Smith   ierr = ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGHASH, ISLocalToGlobalMappingCreate_Hash);CHKERRQ(ierr);
1903413f72f0SBarry Smith   PetscFunctionReturn(0);
1904413f72f0SBarry Smith }
190504a59952SBarry Smith 
1906