xref: /petsc/src/vec/is/utils/isltog.c (revision 071fcb05f2a6aea8aef2c2090580530b9e9df77e)
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 
199305a4c7SMatthew G. Knepley PetscErrorCode ISGetPointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points)
209305a4c7SMatthew G. Knepley {
219305a4c7SMatthew G. Knepley   PetscInt       numCells, step = 1;
229305a4c7SMatthew G. Knepley   PetscBool      isStride;
239305a4c7SMatthew G. Knepley   PetscErrorCode ierr;
249305a4c7SMatthew G. Knepley 
259305a4c7SMatthew G. Knepley   PetscFunctionBeginHot;
269305a4c7SMatthew G. Knepley   *pStart = 0;
279305a4c7SMatthew G. Knepley   *points = NULL;
289305a4c7SMatthew G. Knepley   ierr = ISGetLocalSize(pointIS, &numCells);CHKERRQ(ierr);
299305a4c7SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) pointIS, ISSTRIDE, &isStride);CHKERRQ(ierr);
309305a4c7SMatthew G. Knepley   if (isStride) {ierr = ISStrideGetInfo(pointIS, pStart, &step);CHKERRQ(ierr);}
319305a4c7SMatthew G. Knepley   *pEnd   = *pStart + numCells;
329305a4c7SMatthew G. Knepley   if (!isStride || step != 1) {ierr = ISGetIndices(pointIS, points);CHKERRQ(ierr);}
339305a4c7SMatthew G. Knepley   PetscFunctionReturn(0);
349305a4c7SMatthew G. Knepley }
359305a4c7SMatthew G. Knepley 
369305a4c7SMatthew G. Knepley PetscErrorCode ISRestorePointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points)
379305a4c7SMatthew G. Knepley {
389305a4c7SMatthew G. Knepley   PetscInt       step = 1;
399305a4c7SMatthew G. Knepley   PetscBool      isStride;
409305a4c7SMatthew G. Knepley   PetscErrorCode ierr;
419305a4c7SMatthew G. Knepley 
429305a4c7SMatthew G. Knepley   PetscFunctionBeginHot;
439305a4c7SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) pointIS, ISSTRIDE, &isStride);CHKERRQ(ierr);
449305a4c7SMatthew G. Knepley   if (isStride) {ierr = ISStrideGetInfo(pointIS, pStart, &step);CHKERRQ(ierr);}
459305a4c7SMatthew G. Knepley   if (!isStride || step != 1) {ierr = ISGetIndices(pointIS, points);CHKERRQ(ierr);}
469305a4c7SMatthew G. Knepley   PetscFunctionReturn(0);
479305a4c7SMatthew G. Knepley }
489305a4c7SMatthew G. Knepley 
499305a4c7SMatthew G. Knepley PetscErrorCode ISGetPointSubrange(IS subpointIS, PetscInt pStart, PetscInt pEnd, const PetscInt *points)
509305a4c7SMatthew G. Knepley {
519305a4c7SMatthew G. Knepley   PetscErrorCode ierr;
529305a4c7SMatthew G. Knepley 
539305a4c7SMatthew G. Knepley   PetscFunctionBeginHot;
549305a4c7SMatthew G. Knepley   if (points) {
559305a4c7SMatthew G. Knepley     ierr = ISSetType(subpointIS, ISGENERAL);CHKERRQ(ierr);
569305a4c7SMatthew G. Knepley     ierr = ISGeneralSetIndices(subpointIS, pEnd-pStart, &points[pStart], PETSC_USE_POINTER);CHKERRQ(ierr);
579305a4c7SMatthew G. Knepley   } else {
589305a4c7SMatthew G. Knepley     ierr = ISSetType(subpointIS, ISSTRIDE);CHKERRQ(ierr);
599305a4c7SMatthew G. Knepley     ierr = ISStrideSetStride(subpointIS, pEnd-pStart, pStart, 1);CHKERRQ(ierr);
609305a4c7SMatthew G. Knepley   }
619305a4c7SMatthew G. Knepley   PetscFunctionReturn(0);
629305a4c7SMatthew G. Knepley }
639305a4c7SMatthew G. Knepley 
64413f72f0SBarry Smith /* -----------------------------------------------------------------------------------------*/
65413f72f0SBarry Smith 
66413f72f0SBarry Smith /*
67413f72f0SBarry Smith     Creates the global mapping information in the ISLocalToGlobalMapping structure
68413f72f0SBarry Smith 
69413f72f0SBarry Smith     If the user has not selected how to handle the global to local mapping then use HASH for "large" problems
70413f72f0SBarry Smith */
71413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp(ISLocalToGlobalMapping mapping)
72413f72f0SBarry Smith {
73413f72f0SBarry Smith   PetscInt       i,*idx = mapping->indices,n = mapping->n,end,start;
74413f72f0SBarry Smith   PetscErrorCode ierr;
75413f72f0SBarry Smith 
76413f72f0SBarry Smith   PetscFunctionBegin;
77413f72f0SBarry Smith   if (mapping->data) PetscFunctionReturn(0);
78413f72f0SBarry Smith   end   = 0;
79413f72f0SBarry Smith   start = PETSC_MAX_INT;
80413f72f0SBarry Smith 
81413f72f0SBarry Smith   for (i=0; i<n; i++) {
82413f72f0SBarry Smith     if (idx[i] < 0) continue;
83413f72f0SBarry Smith     if (idx[i] < start) start = idx[i];
84413f72f0SBarry Smith     if (idx[i] > end)   end   = idx[i];
85413f72f0SBarry Smith   }
86413f72f0SBarry Smith   if (start > end) {start = 0; end = -1;}
87413f72f0SBarry Smith   mapping->globalstart = start;
88413f72f0SBarry Smith   mapping->globalend   = end;
89413f72f0SBarry Smith   if (!((PetscObject)mapping)->type_name) {
90413f72f0SBarry Smith     if ((end - start) > PetscMax(4*n,1000000)) {
917f79407eSBarry Smith       ierr = ISLocalToGlobalMappingSetType(mapping,ISLOCALTOGLOBALMAPPINGHASH);CHKERRQ(ierr);
92413f72f0SBarry Smith     } else {
937f79407eSBarry Smith       ierr = ISLocalToGlobalMappingSetType(mapping,ISLOCALTOGLOBALMAPPINGBASIC);CHKERRQ(ierr);
94413f72f0SBarry Smith     }
95413f72f0SBarry Smith   }
96413f72f0SBarry Smith   ierr = (*mapping->ops->globaltolocalmappingsetup)(mapping);CHKERRQ(ierr);
97413f72f0SBarry Smith   PetscFunctionReturn(0);
98413f72f0SBarry Smith }
99413f72f0SBarry Smith 
100413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp_Basic(ISLocalToGlobalMapping mapping)
101413f72f0SBarry Smith {
102413f72f0SBarry Smith   PetscErrorCode              ierr;
103413f72f0SBarry Smith   PetscInt                    i,*idx = mapping->indices,n = mapping->n,end,start,*globals;
104413f72f0SBarry Smith   ISLocalToGlobalMapping_Basic *map;
105413f72f0SBarry Smith 
106413f72f0SBarry Smith   PetscFunctionBegin;
107413f72f0SBarry Smith   start            = mapping->globalstart;
108413f72f0SBarry Smith   end              = mapping->globalend;
109413f72f0SBarry Smith   ierr             = PetscNew(&map);CHKERRQ(ierr);
110413f72f0SBarry Smith   ierr             = PetscMalloc1(end-start+2,&globals);CHKERRQ(ierr);
111413f72f0SBarry Smith   map->globals     = globals;
112413f72f0SBarry Smith   for (i=0; i<end-start+1; i++) globals[i] = -1;
113413f72f0SBarry Smith   for (i=0; i<n; i++) {
114413f72f0SBarry Smith     if (idx[i] < 0) continue;
115413f72f0SBarry Smith     globals[idx[i] - start] = i;
116413f72f0SBarry Smith   }
117413f72f0SBarry Smith   mapping->data = (void*)map;
118413f72f0SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)mapping,(end-start+1)*sizeof(PetscInt));CHKERRQ(ierr);
119413f72f0SBarry Smith   PetscFunctionReturn(0);
120413f72f0SBarry Smith }
121413f72f0SBarry Smith 
122413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp_Hash(ISLocalToGlobalMapping mapping)
123413f72f0SBarry Smith {
124413f72f0SBarry Smith   PetscErrorCode              ierr;
125413f72f0SBarry Smith   PetscInt                    i,*idx = mapping->indices,n = mapping->n;
126413f72f0SBarry Smith   ISLocalToGlobalMapping_Hash *map;
127413f72f0SBarry Smith 
128413f72f0SBarry Smith   PetscFunctionBegin;
129413f72f0SBarry Smith   ierr = PetscNew(&map);CHKERRQ(ierr);
130e8f14785SLisandro Dalcin   ierr = PetscHMapICreate(&map->globalht);CHKERRQ(ierr);
131413f72f0SBarry Smith   for (i=0; i<n; i++ ) {
132413f72f0SBarry Smith     if (idx[i] < 0) continue;
133e8f14785SLisandro Dalcin     ierr = PetscHMapISet(map->globalht,idx[i],i);CHKERRQ(ierr);
134413f72f0SBarry Smith   }
135413f72f0SBarry Smith   mapping->data = (void*)map;
136413f72f0SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)mapping,2*n*sizeof(PetscInt));CHKERRQ(ierr);
137413f72f0SBarry Smith   PetscFunctionReturn(0);
138413f72f0SBarry Smith }
139413f72f0SBarry Smith 
140413f72f0SBarry Smith static PetscErrorCode ISLocalToGlobalMappingDestroy_Basic(ISLocalToGlobalMapping mapping)
141413f72f0SBarry Smith {
142413f72f0SBarry Smith   PetscErrorCode              ierr;
143413f72f0SBarry Smith   ISLocalToGlobalMapping_Basic *map  = (ISLocalToGlobalMapping_Basic *)mapping->data;
144413f72f0SBarry Smith 
145413f72f0SBarry Smith   PetscFunctionBegin;
146413f72f0SBarry Smith   if (!map) PetscFunctionReturn(0);
147413f72f0SBarry Smith   ierr = PetscFree(map->globals);CHKERRQ(ierr);
148413f72f0SBarry Smith   ierr = PetscFree(mapping->data);CHKERRQ(ierr);
149413f72f0SBarry Smith   PetscFunctionReturn(0);
150413f72f0SBarry Smith }
151413f72f0SBarry Smith 
152413f72f0SBarry Smith static PetscErrorCode ISLocalToGlobalMappingDestroy_Hash(ISLocalToGlobalMapping mapping)
153413f72f0SBarry Smith {
154413f72f0SBarry Smith   PetscErrorCode              ierr;
155413f72f0SBarry Smith   ISLocalToGlobalMapping_Hash *map  = (ISLocalToGlobalMapping_Hash*)mapping->data;
156413f72f0SBarry Smith 
157413f72f0SBarry Smith   PetscFunctionBegin;
158413f72f0SBarry Smith   if (!map) PetscFunctionReturn(0);
159e8f14785SLisandro Dalcin   ierr = PetscHMapIDestroy(&map->globalht);CHKERRQ(ierr);
160413f72f0SBarry Smith   ierr = PetscFree(mapping->data);CHKERRQ(ierr);
161413f72f0SBarry Smith   PetscFunctionReturn(0);
162413f72f0SBarry Smith }
163413f72f0SBarry Smith 
164413f72f0SBarry Smith #define GTOLTYPE _Basic
165413f72f0SBarry Smith #define GTOLNAME _Basic
166541bf97eSAdrian Croucher #define GTOLBS mapping->bs
167413f72f0SBarry Smith #define GTOL(g, local) do {                  \
168413f72f0SBarry Smith     local = map->globals[g/bs - start];      \
169413f72f0SBarry Smith     local = bs*local + (g % bs);             \
170413f72f0SBarry Smith   } while (0)
171541bf97eSAdrian Croucher 
172413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
173413f72f0SBarry Smith 
174413f72f0SBarry Smith #define GTOLTYPE _Basic
175413f72f0SBarry Smith #define GTOLNAME Block_Basic
176541bf97eSAdrian Croucher #define GTOLBS 1
177413f72f0SBarry Smith #define GTOL(g, local) do {                  \
178413f72f0SBarry Smith     local = map->globals[g - start];         \
179413f72f0SBarry Smith   } while (0)
180413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
181413f72f0SBarry Smith 
182413f72f0SBarry Smith #define GTOLTYPE _Hash
183413f72f0SBarry Smith #define GTOLNAME _Hash
184541bf97eSAdrian Croucher #define GTOLBS mapping->bs
185413f72f0SBarry Smith #define GTOL(g, local) do {                         \
186e8f14785SLisandro Dalcin     (void)PetscHMapIGet(map->globalht,g/bs,&local); \
187413f72f0SBarry Smith     local = bs*local + (g % bs);                    \
188413f72f0SBarry Smith    } while (0)
189413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
190413f72f0SBarry Smith 
191413f72f0SBarry Smith #define GTOLTYPE _Hash
192413f72f0SBarry Smith #define GTOLNAME Block_Hash
193541bf97eSAdrian Croucher #define GTOLBS 1
194413f72f0SBarry Smith #define GTOL(g, local) do {                         \
195e8f14785SLisandro Dalcin     (void)PetscHMapIGet(map->globalht,g,&local);    \
196413f72f0SBarry Smith   } while (0)
197413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
198413f72f0SBarry Smith 
1996658fb44Sstefano_zampini /*@
2006658fb44Sstefano_zampini     ISLocalToGlobalMappingDuplicate - Duplicates the local to global mapping object
2016658fb44Sstefano_zampini 
2026658fb44Sstefano_zampini     Not Collective
2036658fb44Sstefano_zampini 
2046658fb44Sstefano_zampini     Input Parameter:
2056658fb44Sstefano_zampini .   ltog - local to global mapping
2066658fb44Sstefano_zampini 
2076658fb44Sstefano_zampini     Output Parameter:
2086658fb44Sstefano_zampini .   nltog - the duplicated local to global mapping
2096658fb44Sstefano_zampini 
2106658fb44Sstefano_zampini     Level: advanced
2116658fb44Sstefano_zampini 
2126658fb44Sstefano_zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
2136658fb44Sstefano_zampini @*/
2146658fb44Sstefano_zampini PetscErrorCode  ISLocalToGlobalMappingDuplicate(ISLocalToGlobalMapping ltog,ISLocalToGlobalMapping* nltog)
2156658fb44Sstefano_zampini {
2166658fb44Sstefano_zampini   PetscErrorCode ierr;
2176658fb44Sstefano_zampini 
2186658fb44Sstefano_zampini   PetscFunctionBegin;
2196658fb44Sstefano_zampini   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
2206658fb44Sstefano_zampini   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)ltog),ltog->bs,ltog->n,ltog->indices,PETSC_COPY_VALUES,nltog);CHKERRQ(ierr);
2216658fb44Sstefano_zampini   PetscFunctionReturn(0);
2226658fb44Sstefano_zampini }
2236658fb44Sstefano_zampini 
224565245c5SBarry Smith /*@
225107e9a97SBarry Smith     ISLocalToGlobalMappingGetSize - Gets the local size of a local to global mapping
2263b9aefa3SBarry Smith 
2273b9aefa3SBarry Smith     Not Collective
2283b9aefa3SBarry Smith 
2293b9aefa3SBarry Smith     Input Parameter:
2303b9aefa3SBarry Smith .   ltog - local to global mapping
2313b9aefa3SBarry Smith 
2323b9aefa3SBarry Smith     Output Parameter:
233107e9a97SBarry Smith .   n - the number of entries in the local mapping, ISLocalToGlobalMappingGetIndices() returns an array of this length
2343b9aefa3SBarry Smith 
2353b9aefa3SBarry Smith     Level: advanced
2363b9aefa3SBarry Smith 
2373b9aefa3SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
2383b9aefa3SBarry Smith @*/
2397087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping mapping,PetscInt *n)
2403b9aefa3SBarry Smith {
2413b9aefa3SBarry Smith   PetscFunctionBegin;
2420700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
2434482741eSBarry Smith   PetscValidIntPointer(n,2);
244107e9a97SBarry Smith   *n = mapping->bs*mapping->n;
2453b9aefa3SBarry Smith   PetscFunctionReturn(0);
2463b9aefa3SBarry Smith }
2473b9aefa3SBarry Smith 
2485a5d4f66SBarry Smith /*@C
2495a5d4f66SBarry Smith     ISLocalToGlobalMappingView - View a local to global mapping
2505a5d4f66SBarry Smith 
251b9cd556bSLois Curfman McInnes     Not Collective
252b9cd556bSLois Curfman McInnes 
2535a5d4f66SBarry Smith     Input Parameters:
2543b9aefa3SBarry Smith +   ltog - local to global mapping
2553b9aefa3SBarry Smith -   viewer - viewer
2565a5d4f66SBarry Smith 
257a997ad1aSLois Curfman McInnes     Level: advanced
258a997ad1aSLois Curfman McInnes 
2595a5d4f66SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
2605a5d4f66SBarry Smith @*/
2617087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingView(ISLocalToGlobalMapping mapping,PetscViewer viewer)
2625a5d4f66SBarry Smith {
26332dcc486SBarry Smith   PetscInt       i;
26432dcc486SBarry Smith   PetscMPIInt    rank;
265ace3abfcSBarry Smith   PetscBool      iascii;
2666849ba73SBarry Smith   PetscErrorCode ierr;
2675a5d4f66SBarry Smith 
2685a5d4f66SBarry Smith   PetscFunctionBegin;
2690700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
2703050cee2SBarry Smith   if (!viewer) {
271ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mapping),&viewer);CHKERRQ(ierr);
2723050cee2SBarry Smith   }
2730700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
2745a5d4f66SBarry Smith 
275ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mapping),&rank);CHKERRQ(ierr);
276251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
27732077d6dSBarry Smith   if (iascii) {
27898c3331eSBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)mapping,viewer);CHKERRQ(ierr);
2791575c14dSBarry Smith     ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
2805a5d4f66SBarry Smith     for (i=0; i<mapping->n; i++) {
2817904a332SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D %D\n",rank,i,mapping->indices[i]);CHKERRQ(ierr);
2826831982aSBarry Smith     }
283b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
2841575c14dSBarry Smith     ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
2851575c14dSBarry Smith   }
2865a5d4f66SBarry Smith   PetscFunctionReturn(0);
2875a5d4f66SBarry Smith }
2885a5d4f66SBarry Smith 
2891f428162SBarry Smith /*@
2902bdab257SBarry Smith     ISLocalToGlobalMappingCreateIS - Creates a mapping between a local (0 to n)
2912bdab257SBarry Smith     ordering and a global parallel ordering.
2922bdab257SBarry Smith 
2930f5bd95cSBarry Smith     Not collective
294b9cd556bSLois Curfman McInnes 
295a997ad1aSLois Curfman McInnes     Input Parameter:
2968c03b21aSDmitry Karpeev .   is - index set containing the global numbers for each local number
2972bdab257SBarry Smith 
298a997ad1aSLois Curfman McInnes     Output Parameter:
2992bdab257SBarry Smith .   mapping - new mapping data structure
3002bdab257SBarry Smith 
30195452b02SPatrick Sanan     Notes:
30295452b02SPatrick Sanan     the block size of the IS determines the block size of the mapping
303a997ad1aSLois Curfman McInnes     Level: advanced
304a997ad1aSLois Curfman McInnes 
3057e99dc12SLawrence Mitchell .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetFromOptions()
3062bdab257SBarry Smith @*/
3077087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingCreateIS(IS is,ISLocalToGlobalMapping *mapping)
3082bdab257SBarry Smith {
3096849ba73SBarry Smith   PetscErrorCode ierr;
3103bbf0e92SBarry Smith   PetscInt       n,bs;
3115d0c19d7SBarry Smith   const PetscInt *indices;
3122bdab257SBarry Smith   MPI_Comm       comm;
3133bbf0e92SBarry Smith   PetscBool      isblock;
3143a40ed3dSBarry Smith 
3153a40ed3dSBarry Smith   PetscFunctionBegin;
3160700a824SBarry Smith   PetscValidHeaderSpecific(is,IS_CLASSID,1);
3174482741eSBarry Smith   PetscValidPointer(mapping,2);
3182bdab257SBarry Smith 
3192bdab257SBarry Smith   ierr = PetscObjectGetComm((PetscObject)is,&comm);CHKERRQ(ierr);
3203b9aefa3SBarry Smith   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
3213bbf0e92SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)is,ISBLOCK,&isblock);CHKERRQ(ierr);
3226006e8d2SBarry Smith   if (!isblock) {
323f0413b6fSBarry Smith     ierr = ISGetIndices(is,&indices);CHKERRQ(ierr);
324f0413b6fSBarry Smith     ierr = ISLocalToGlobalMappingCreate(comm,1,n,indices,PETSC_COPY_VALUES,mapping);CHKERRQ(ierr);
3252bdab257SBarry Smith     ierr = ISRestoreIndices(is,&indices);CHKERRQ(ierr);
3266006e8d2SBarry Smith   } else {
3276006e8d2SBarry Smith     ierr = ISGetBlockSize(is,&bs);CHKERRQ(ierr);
328f0413b6fSBarry Smith     ierr = ISBlockGetIndices(is,&indices);CHKERRQ(ierr);
32928bc9809SBarry Smith     ierr = ISLocalToGlobalMappingCreate(comm,bs,n/bs,indices,PETSC_COPY_VALUES,mapping);CHKERRQ(ierr);
330f0413b6fSBarry Smith     ierr = ISBlockRestoreIndices(is,&indices);CHKERRQ(ierr);
3316006e8d2SBarry Smith   }
3323a40ed3dSBarry Smith   PetscFunctionReturn(0);
3332bdab257SBarry Smith }
3345a5d4f66SBarry Smith 
335a4d96a55SJed Brown /*@C
336a4d96a55SJed Brown     ISLocalToGlobalMappingCreateSF - Creates a mapping between a local (0 to n)
337a4d96a55SJed Brown     ordering and a global parallel ordering.
338a4d96a55SJed Brown 
339a4d96a55SJed Brown     Collective
340a4d96a55SJed Brown 
341a4d96a55SJed Brown     Input Parameter:
342a4d96a55SJed Brown +   sf - star forest mapping contiguous local indices to (rank, offset)
343a4d96a55SJed Brown -   start - first global index on this process
344a4d96a55SJed Brown 
345a4d96a55SJed Brown     Output Parameter:
346a4d96a55SJed Brown .   mapping - new mapping data structure
347a4d96a55SJed Brown 
348a4d96a55SJed Brown     Level: advanced
349a4d96a55SJed Brown 
3507e99dc12SLawrence Mitchell .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingSetFromOptions()
351a4d96a55SJed Brown @*/
352a4d96a55SJed Brown PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF sf,PetscInt start,ISLocalToGlobalMapping *mapping)
353a4d96a55SJed Brown {
354a4d96a55SJed Brown   PetscErrorCode ierr;
355a4d96a55SJed Brown   PetscInt       i,maxlocal,nroots,nleaves,*globals,*ltog;
356a4d96a55SJed Brown   const PetscInt *ilocal;
357a4d96a55SJed Brown   MPI_Comm       comm;
358a4d96a55SJed Brown 
359a4d96a55SJed Brown   PetscFunctionBegin;
360a4d96a55SJed Brown   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
361a4d96a55SJed Brown   PetscValidPointer(mapping,3);
362a4d96a55SJed Brown 
363a4d96a55SJed Brown   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
3640298fd71SBarry Smith   ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr);
365f6e5521dSKarl Rupp   if (ilocal) {
366f6e5521dSKarl Rupp     for (i=0,maxlocal=0; i<nleaves; i++) maxlocal = PetscMax(maxlocal,ilocal[i]+1);
367f6e5521dSKarl Rupp   }
368a4d96a55SJed Brown   else maxlocal = nleaves;
369785e854fSJed Brown   ierr = PetscMalloc1(nroots,&globals);CHKERRQ(ierr);
370785e854fSJed Brown   ierr = PetscMalloc1(maxlocal,&ltog);CHKERRQ(ierr);
371a4d96a55SJed Brown   for (i=0; i<nroots; i++) globals[i] = start + i;
372a4d96a55SJed Brown   for (i=0; i<maxlocal; i++) ltog[i] = -1;
373a4d96a55SJed Brown   ierr = PetscSFBcastBegin(sf,MPIU_INT,globals,ltog);CHKERRQ(ierr);
374a4d96a55SJed Brown   ierr = PetscSFBcastEnd(sf,MPIU_INT,globals,ltog);CHKERRQ(ierr);
375f0413b6fSBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,1,maxlocal,ltog,PETSC_OWN_POINTER,mapping);CHKERRQ(ierr);
376a4d96a55SJed Brown   ierr = PetscFree(globals);CHKERRQ(ierr);
377a4d96a55SJed Brown   PetscFunctionReturn(0);
378a4d96a55SJed Brown }
379b46b645bSBarry Smith 
38063fa5c83Sstefano_zampini /*@
38163fa5c83Sstefano_zampini     ISLocalToGlobalMappingSetBlockSize - Sets the blocksize of the mapping
38263fa5c83Sstefano_zampini 
38363fa5c83Sstefano_zampini     Not collective
38463fa5c83Sstefano_zampini 
38563fa5c83Sstefano_zampini     Input Parameters:
38663fa5c83Sstefano_zampini .   mapping - mapping data structure
38763fa5c83Sstefano_zampini .   bs - the blocksize
38863fa5c83Sstefano_zampini 
38963fa5c83Sstefano_zampini     Level: advanced
39063fa5c83Sstefano_zampini 
39163fa5c83Sstefano_zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS()
39263fa5c83Sstefano_zampini @*/
39363fa5c83Sstefano_zampini PetscErrorCode  ISLocalToGlobalMappingSetBlockSize(ISLocalToGlobalMapping mapping,PetscInt bs)
39463fa5c83Sstefano_zampini {
395a59f3c4dSstefano_zampini   PetscInt       *nid;
396a59f3c4dSstefano_zampini   const PetscInt *oid;
397a59f3c4dSstefano_zampini   PetscInt       i,cn,on,obs,nn;
39863fa5c83Sstefano_zampini   PetscErrorCode ierr;
39963fa5c83Sstefano_zampini 
40063fa5c83Sstefano_zampini   PetscFunctionBegin;
40163fa5c83Sstefano_zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
40263fa5c83Sstefano_zampini   if (bs < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid block size %D",bs);
40363fa5c83Sstefano_zampini   if (bs == mapping->bs) PetscFunctionReturn(0);
40463fa5c83Sstefano_zampini   on  = mapping->n;
40563fa5c83Sstefano_zampini   obs = mapping->bs;
40663fa5c83Sstefano_zampini   oid = mapping->indices;
40763fa5c83Sstefano_zampini   nn  = (on*obs)/bs;
40863fa5c83Sstefano_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);
409a59f3c4dSstefano_zampini 
41063fa5c83Sstefano_zampini   ierr = PetscMalloc1(nn,&nid);CHKERRQ(ierr);
411a59f3c4dSstefano_zampini   ierr = ISLocalToGlobalMappingGetIndices(mapping,&oid);CHKERRQ(ierr);
412a59f3c4dSstefano_zampini   for (i=0;i<nn;i++) {
413a59f3c4dSstefano_zampini     PetscInt j;
414a59f3c4dSstefano_zampini     for (j=0,cn=0;j<bs-1;j++) {
415a59f3c4dSstefano_zampini       if (oid[i*bs+j] < 0) { cn++; continue; }
416a59f3c4dSstefano_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]);
417a59f3c4dSstefano_zampini     }
418a59f3c4dSstefano_zampini     if (oid[i*bs+j] < 0) cn++;
4198b7cb0e6Sstefano_zampini     if (cn) {
420a59f3c4dSstefano_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);
421a59f3c4dSstefano_zampini       nid[i] = -1;
4228b7cb0e6Sstefano_zampini     } else {
423a59f3c4dSstefano_zampini       nid[i] = oid[i*bs]/bs;
42463fa5c83Sstefano_zampini     }
42563fa5c83Sstefano_zampini   }
426a59f3c4dSstefano_zampini   ierr = ISLocalToGlobalMappingRestoreIndices(mapping,&oid);CHKERRQ(ierr);
427a59f3c4dSstefano_zampini 
42863fa5c83Sstefano_zampini   mapping->n           = nn;
42963fa5c83Sstefano_zampini   mapping->bs          = bs;
43063fa5c83Sstefano_zampini   ierr                 = PetscFree(mapping->indices);CHKERRQ(ierr);
43163fa5c83Sstefano_zampini   mapping->indices     = nid;
432c9345713Sstefano_zampini   mapping->globalstart = 0;
433c9345713Sstefano_zampini   mapping->globalend   = 0;
4341bd0b88eSStefano Zampini 
4351bd0b88eSStefano Zampini   /* reset the cached information */
4361bd0b88eSStefano Zampini   ierr = PetscFree(mapping->info_procs);CHKERRQ(ierr);
4371bd0b88eSStefano Zampini   ierr = PetscFree(mapping->info_numprocs);CHKERRQ(ierr);
4381bd0b88eSStefano Zampini   if (mapping->info_indices) {
4391bd0b88eSStefano Zampini     PetscInt i;
4401bd0b88eSStefano Zampini 
4411bd0b88eSStefano Zampini     ierr = PetscFree((mapping->info_indices)[0]);CHKERRQ(ierr);
4421bd0b88eSStefano Zampini     for (i=1; i<mapping->info_nproc; i++) {
4431bd0b88eSStefano Zampini       ierr = PetscFree(mapping->info_indices[i]);CHKERRQ(ierr);
4441bd0b88eSStefano Zampini     }
4451bd0b88eSStefano Zampini     ierr = PetscFree(mapping->info_indices);CHKERRQ(ierr);
4461bd0b88eSStefano Zampini   }
4471bd0b88eSStefano Zampini   mapping->info_cached = PETSC_FALSE;
4481bd0b88eSStefano Zampini 
449413f72f0SBarry Smith   if (mapping->ops->destroy) {
450413f72f0SBarry Smith     ierr = (*mapping->ops->destroy)(mapping);CHKERRQ(ierr);
451413f72f0SBarry Smith   }
45263fa5c83Sstefano_zampini   PetscFunctionReturn(0);
45363fa5c83Sstefano_zampini }
45463fa5c83Sstefano_zampini 
45545b6f7e9SBarry Smith /*@
45645b6f7e9SBarry Smith     ISLocalToGlobalMappingGetBlockSize - Gets the blocksize of the mapping
45745b6f7e9SBarry Smith     ordering and a global parallel ordering.
45845b6f7e9SBarry Smith 
45945b6f7e9SBarry Smith     Not Collective
46045b6f7e9SBarry Smith 
46145b6f7e9SBarry Smith     Input Parameters:
46245b6f7e9SBarry Smith .   mapping - mapping data structure
46345b6f7e9SBarry Smith 
46445b6f7e9SBarry Smith     Output Parameter:
46545b6f7e9SBarry Smith .   bs - the blocksize
46645b6f7e9SBarry Smith 
46745b6f7e9SBarry Smith     Level: advanced
46845b6f7e9SBarry Smith 
46945b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS()
47045b6f7e9SBarry Smith @*/
47145b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping mapping,PetscInt *bs)
47245b6f7e9SBarry Smith {
47345b6f7e9SBarry Smith   PetscFunctionBegin;
474cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
47545b6f7e9SBarry Smith   *bs = mapping->bs;
47645b6f7e9SBarry Smith   PetscFunctionReturn(0);
47745b6f7e9SBarry Smith }
47845b6f7e9SBarry Smith 
479ba5bb76aSSatish Balay /*@
48090f02eecSBarry Smith     ISLocalToGlobalMappingCreate - Creates a mapping between a local (0 to n)
48190f02eecSBarry Smith     ordering and a global parallel ordering.
4822362add9SBarry Smith 
48389d82c54SBarry Smith     Not Collective, but communicator may have more than one process
484b9cd556bSLois Curfman McInnes 
4852362add9SBarry Smith     Input Parameters:
48689d82c54SBarry Smith +   comm - MPI communicator
487f0413b6fSBarry Smith .   bs - the block size
48828bc9809SBarry Smith .   n - the number of local elements divided by the block size, or equivalently the number of block indices
48928bc9809SBarry 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
490d5ad8652SBarry Smith -   mode - see PetscCopyMode
4912362add9SBarry Smith 
492a997ad1aSLois Curfman McInnes     Output Parameter:
49390f02eecSBarry Smith .   mapping - new mapping data structure
4942362add9SBarry Smith 
49595452b02SPatrick Sanan     Notes:
49695452b02SPatrick Sanan     There is one integer value in indices per block and it represents the actual indices bs*idx + j, where j=0,..,bs-1
497413f72f0SBarry Smith 
4989a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
499413f72f0SBarry 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.
500413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
501413f72f0SBarry Smith 
502a997ad1aSLois Curfman McInnes     Level: advanced
503a997ad1aSLois Curfman McInnes 
504413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingSetFromOptions(), ISLOCALTOGLOBALMAPPINGBASIC, ISLOCALTOGLOBALMAPPINGHASH
505413f72f0SBarry Smith           ISLocalToGlobalMappingSetType(), ISLocalToGlobalMappingType
5062362add9SBarry Smith @*/
50760c7cefcSBarry Smith PetscErrorCode  ISLocalToGlobalMappingCreate(MPI_Comm comm,PetscInt bs,PetscInt n,const PetscInt indices[],PetscCopyMode mode,ISLocalToGlobalMapping *mapping)
5082362add9SBarry Smith {
5096849ba73SBarry Smith   PetscErrorCode ierr;
51032dcc486SBarry Smith   PetscInt       *in;
511b46b645bSBarry Smith 
512b46b645bSBarry Smith   PetscFunctionBegin;
51373911063SBarry Smith   if (n) PetscValidIntPointer(indices,3);
5144482741eSBarry Smith   PetscValidPointer(mapping,4);
515b46b645bSBarry Smith 
5160298fd71SBarry Smith   *mapping = NULL;
517607a6623SBarry Smith   ierr = ISInitializePackage();CHKERRQ(ierr);
5182362add9SBarry Smith 
51973107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(*mapping,IS_LTOGM_CLASSID,"ISLocalToGlobalMapping","Local to global mapping","IS",
52060c7cefcSBarry Smith                            comm,ISLocalToGlobalMappingDestroy,ISLocalToGlobalMappingView);CHKERRQ(ierr);
521d4bb536fSBarry Smith   (*mapping)->n             = n;
522f0413b6fSBarry Smith   (*mapping)->bs            = bs;
523268a049cSStefano Zampini   (*mapping)->info_cached   = PETSC_FALSE;
524268a049cSStefano Zampini   (*mapping)->info_free     = PETSC_FALSE;
525268a049cSStefano Zampini   (*mapping)->info_procs    = NULL;
526268a049cSStefano Zampini   (*mapping)->info_numprocs = NULL;
527268a049cSStefano Zampini   (*mapping)->info_indices  = NULL;
5281bd0b88eSStefano Zampini   (*mapping)->info_nodec    = NULL;
5291bd0b88eSStefano Zampini   (*mapping)->info_nodei    = NULL;
530413f72f0SBarry Smith 
531413f72f0SBarry Smith   (*mapping)->ops->globaltolocalmappingapply      = NULL;
532413f72f0SBarry Smith   (*mapping)->ops->globaltolocalmappingapplyblock = NULL;
533413f72f0SBarry Smith   (*mapping)->ops->destroy                        = NULL;
534d5ad8652SBarry Smith   if (mode == PETSC_COPY_VALUES) {
535785e854fSJed Brown     ierr = PetscMalloc1(n,&in);CHKERRQ(ierr);
536580bdb30SBarry Smith     ierr = PetscArraycpy(in,indices,n);CHKERRQ(ierr);
537d5ad8652SBarry Smith     (*mapping)->indices = in;
5386389a1a1SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt));CHKERRQ(ierr);
5396389a1a1SBarry Smith   } else if (mode == PETSC_OWN_POINTER) {
5406389a1a1SBarry Smith     (*mapping)->indices = (PetscInt*)indices;
5416389a1a1SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt));CHKERRQ(ierr);
5426389a1a1SBarry Smith   }
54360c7cefcSBarry Smith   else SETERRQ(comm,PETSC_ERR_SUP,"Cannot currently use PETSC_USE_POINTER");
5443a40ed3dSBarry Smith   PetscFunctionReturn(0);
5452362add9SBarry Smith }
5462362add9SBarry Smith 
547413f72f0SBarry Smith PetscFunctionList ISLocalToGlobalMappingList = NULL;
548413f72f0SBarry Smith 
549413f72f0SBarry Smith 
55090f02eecSBarry Smith /*@
5517e99dc12SLawrence Mitchell    ISLocalToGlobalMappingSetFromOptions - Set mapping options from the options database.
5527e99dc12SLawrence Mitchell 
5537e99dc12SLawrence Mitchell    Not collective
5547e99dc12SLawrence Mitchell 
5557e99dc12SLawrence Mitchell    Input Parameters:
5567e99dc12SLawrence Mitchell .  mapping - mapping data structure
5577e99dc12SLawrence Mitchell 
5587e99dc12SLawrence Mitchell    Level: advanced
5597e99dc12SLawrence Mitchell 
5607e99dc12SLawrence Mitchell @*/
5617e99dc12SLawrence Mitchell PetscErrorCode ISLocalToGlobalMappingSetFromOptions(ISLocalToGlobalMapping mapping)
5627e99dc12SLawrence Mitchell {
5637e99dc12SLawrence Mitchell   PetscErrorCode             ierr;
564413f72f0SBarry Smith   char                       type[256];
565413f72f0SBarry Smith   ISLocalToGlobalMappingType defaulttype = "Not set";
5667e99dc12SLawrence Mitchell   PetscBool                  flg;
5677e99dc12SLawrence Mitchell 
5687e99dc12SLawrence Mitchell   PetscFunctionBegin;
5697e99dc12SLawrence Mitchell   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
570413f72f0SBarry Smith   ierr = ISLocalToGlobalMappingRegisterAll();CHKERRQ(ierr);
5717e99dc12SLawrence Mitchell   ierr = PetscObjectOptionsBegin((PetscObject)mapping);CHKERRQ(ierr);
572413f72f0SBarry Smith   ierr = PetscOptionsFList("-islocaltoglobalmapping_type","ISLocalToGlobalMapping method","ISLocalToGlobalMappingSetType",ISLocalToGlobalMappingList,(char*)(((PetscObject)mapping)->type_name) ? ((PetscObject)mapping)->type_name : defaulttype,type,256,&flg);CHKERRQ(ierr);
573413f72f0SBarry Smith   if (flg) {
574413f72f0SBarry Smith     ierr = ISLocalToGlobalMappingSetType(mapping,type);CHKERRQ(ierr);
575413f72f0SBarry Smith   }
5767e99dc12SLawrence Mitchell   ierr = PetscOptionsEnd();CHKERRQ(ierr);
5777e99dc12SLawrence Mitchell   PetscFunctionReturn(0);
5787e99dc12SLawrence Mitchell }
5797e99dc12SLawrence Mitchell 
5807e99dc12SLawrence Mitchell /*@
58190f02eecSBarry Smith    ISLocalToGlobalMappingDestroy - Destroys a mapping between a local (0 to n)
58290f02eecSBarry Smith    ordering and a global parallel ordering.
58390f02eecSBarry Smith 
5840f5bd95cSBarry Smith    Note Collective
585b9cd556bSLois Curfman McInnes 
58690f02eecSBarry Smith    Input Parameters:
58790f02eecSBarry Smith .  mapping - mapping data structure
58890f02eecSBarry Smith 
589a997ad1aSLois Curfman McInnes    Level: advanced
590a997ad1aSLois Curfman McInnes 
5913acfe500SLois Curfman McInnes .seealso: ISLocalToGlobalMappingCreate()
59290f02eecSBarry Smith @*/
5936bf464f9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping *mapping)
59490f02eecSBarry Smith {
595dfbe8321SBarry Smith   PetscErrorCode ierr;
5965fd66863SKarl Rupp 
5973a40ed3dSBarry Smith   PetscFunctionBegin;
5986bf464f9SBarry Smith   if (!*mapping) PetscFunctionReturn(0);
5996bf464f9SBarry Smith   PetscValidHeaderSpecific((*mapping),IS_LTOGM_CLASSID,1);
600997056adSBarry Smith   if (--((PetscObject)(*mapping))->refct > 0) {*mapping = 0;PetscFunctionReturn(0);}
6016bf464f9SBarry Smith   ierr = PetscFree((*mapping)->indices);CHKERRQ(ierr);
602268a049cSStefano Zampini   ierr = PetscFree((*mapping)->info_procs);CHKERRQ(ierr);
603268a049cSStefano Zampini   ierr = PetscFree((*mapping)->info_numprocs);CHKERRQ(ierr);
604268a049cSStefano Zampini   if ((*mapping)->info_indices) {
605268a049cSStefano Zampini     PetscInt i;
606268a049cSStefano Zampini 
607268a049cSStefano Zampini     ierr = PetscFree(((*mapping)->info_indices)[0]);CHKERRQ(ierr);
608268a049cSStefano Zampini     for (i=1; i<(*mapping)->info_nproc; i++) {
609268a049cSStefano Zampini       ierr = PetscFree(((*mapping)->info_indices)[i]);CHKERRQ(ierr);
610268a049cSStefano Zampini     }
611268a049cSStefano Zampini     ierr = PetscFree((*mapping)->info_indices);CHKERRQ(ierr);
612268a049cSStefano Zampini   }
6131bd0b88eSStefano Zampini   if ((*mapping)->info_nodei) {
6141bd0b88eSStefano Zampini     ierr = PetscFree(((*mapping)->info_nodei)[0]);CHKERRQ(ierr);
6151bd0b88eSStefano Zampini   }
616*071fcb05SBarry Smith   ierr = PetscFree2((*mapping)->info_nodec,(*mapping)->info_nodei);CHKERRQ(ierr);
617413f72f0SBarry Smith   if ((*mapping)->ops->destroy) {
618413f72f0SBarry Smith     ierr = (*(*mapping)->ops->destroy)(*mapping);CHKERRQ(ierr);
619413f72f0SBarry Smith   }
620d38fa0fbSBarry Smith   ierr     = PetscHeaderDestroy(mapping);CHKERRQ(ierr);
621992144d0SBarry Smith   *mapping = 0;
6223a40ed3dSBarry Smith   PetscFunctionReturn(0);
62390f02eecSBarry Smith }
62490f02eecSBarry Smith 
62590f02eecSBarry Smith /*@
6263acfe500SLois Curfman McInnes     ISLocalToGlobalMappingApplyIS - Creates from an IS in the local numbering
6273acfe500SLois Curfman McInnes     a new index set using the global numbering defined in an ISLocalToGlobalMapping
6283acfe500SLois Curfman McInnes     context.
62990f02eecSBarry Smith 
6304cb36875SStefano Zampini     Collective on is
631b9cd556bSLois Curfman McInnes 
63290f02eecSBarry Smith     Input Parameters:
633b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
634b9cd556bSLois Curfman McInnes -   is - index set in local numbering
63590f02eecSBarry Smith 
63690f02eecSBarry Smith     Output Parameters:
63790f02eecSBarry Smith .   newis - index set in global numbering
63890f02eecSBarry Smith 
6394cb36875SStefano Zampini     Notes:
6404cb36875SStefano Zampini     The output IS will have the same communicator of the input IS.
6414cb36875SStefano Zampini 
642a997ad1aSLois Curfman McInnes     Level: advanced
643a997ad1aSLois Curfman McInnes 
64490f02eecSBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),
645d4bb536fSBarry Smith           ISLocalToGlobalMappingDestroy(), ISGlobalToLocalMappingApply()
64690f02eecSBarry Smith @*/
6477087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping,IS is,IS *newis)
64890f02eecSBarry Smith {
6496849ba73SBarry Smith   PetscErrorCode ierr;
650e24637baSBarry Smith   PetscInt       n,*idxout;
6515d0c19d7SBarry Smith   const PetscInt *idxin;
6523a40ed3dSBarry Smith 
6533a40ed3dSBarry Smith   PetscFunctionBegin;
6540700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
6550700a824SBarry Smith   PetscValidHeaderSpecific(is,IS_CLASSID,2);
6564482741eSBarry Smith   PetscValidPointer(newis,3);
65790f02eecSBarry Smith 
6583b9aefa3SBarry Smith   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
65990f02eecSBarry Smith   ierr = ISGetIndices(is,&idxin);CHKERRQ(ierr);
660785e854fSJed Brown   ierr = PetscMalloc1(n,&idxout);CHKERRQ(ierr);
661e24637baSBarry Smith   ierr = ISLocalToGlobalMappingApply(mapping,n,idxin,idxout);CHKERRQ(ierr);
6623b9aefa3SBarry Smith   ierr = ISRestoreIndices(is,&idxin);CHKERRQ(ierr);
663543f3098SMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),n,idxout,PETSC_OWN_POINTER,newis);CHKERRQ(ierr);
6643a40ed3dSBarry Smith   PetscFunctionReturn(0);
66590f02eecSBarry Smith }
66690f02eecSBarry Smith 
667b89cb25eSSatish Balay /*@
6683acfe500SLois Curfman McInnes    ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering
6693acfe500SLois Curfman McInnes    and converts them to the global numbering.
67090f02eecSBarry Smith 
671b9cd556bSLois Curfman McInnes    Not collective
672b9cd556bSLois Curfman McInnes 
673bb25748dSBarry Smith    Input Parameters:
674b9cd556bSLois Curfman McInnes +  mapping - the local to global mapping context
675bb25748dSBarry Smith .  N - number of integers
676b9cd556bSLois Curfman McInnes -  in - input indices in local numbering
677bb25748dSBarry Smith 
678bb25748dSBarry Smith    Output Parameter:
679bb25748dSBarry Smith .  out - indices in global numbering
680bb25748dSBarry Smith 
681b9cd556bSLois Curfman McInnes    Notes:
682b9cd556bSLois Curfman McInnes    The in and out array parameters may be identical.
683d4bb536fSBarry Smith 
684a997ad1aSLois Curfman McInnes    Level: advanced
685a997ad1aSLois Curfman McInnes 
68645b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
6870752156aSBarry Smith           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
688d4bb536fSBarry Smith           AOPetscToApplication(), ISGlobalToLocalMappingApply()
689bb25748dSBarry Smith 
690afcb2eb5SJed Brown @*/
691afcb2eb5SJed Brown PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
692afcb2eb5SJed Brown {
693cbc1caf0SMatthew G. Knepley   PetscInt i,bs,Nmax;
69445b6f7e9SBarry Smith 
69545b6f7e9SBarry Smith   PetscFunctionBegin;
696cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
697cbc1caf0SMatthew G. Knepley   bs   = mapping->bs;
698cbc1caf0SMatthew G. Knepley   Nmax = bs*mapping->n;
69945b6f7e9SBarry Smith   if (bs == 1) {
700cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
70145b6f7e9SBarry Smith     for (i=0; i<N; i++) {
70245b6f7e9SBarry Smith       if (in[i] < 0) {
70345b6f7e9SBarry Smith         out[i] = in[i];
70445b6f7e9SBarry Smith         continue;
70545b6f7e9SBarry Smith       }
706e24637baSBarry 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);
70745b6f7e9SBarry Smith       out[i] = idx[in[i]];
70845b6f7e9SBarry Smith     }
70945b6f7e9SBarry Smith   } else {
710cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
71145b6f7e9SBarry Smith     for (i=0; i<N; i++) {
71245b6f7e9SBarry Smith       if (in[i] < 0) {
71345b6f7e9SBarry Smith         out[i] = in[i];
71445b6f7e9SBarry Smith         continue;
71545b6f7e9SBarry Smith       }
716e24637baSBarry 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);
71745b6f7e9SBarry Smith       out[i] = idx[in[i]/bs]*bs + (in[i] % bs);
71845b6f7e9SBarry Smith     }
71945b6f7e9SBarry Smith   }
72045b6f7e9SBarry Smith   PetscFunctionReturn(0);
72145b6f7e9SBarry Smith }
72245b6f7e9SBarry Smith 
72345b6f7e9SBarry Smith /*@
7246006e8d2SBarry Smith    ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local block numbering and converts them to the global block numbering
72545b6f7e9SBarry Smith 
72645b6f7e9SBarry Smith    Not collective
72745b6f7e9SBarry Smith 
72845b6f7e9SBarry Smith    Input Parameters:
72945b6f7e9SBarry Smith +  mapping - the local to global mapping context
73045b6f7e9SBarry Smith .  N - number of integers
7316006e8d2SBarry Smith -  in - input indices in local block numbering
73245b6f7e9SBarry Smith 
73345b6f7e9SBarry Smith    Output Parameter:
7346006e8d2SBarry Smith .  out - indices in global block numbering
73545b6f7e9SBarry Smith 
73645b6f7e9SBarry Smith    Notes:
73745b6f7e9SBarry Smith    The in and out array parameters may be identical.
73845b6f7e9SBarry Smith 
7396006e8d2SBarry Smith    Example:
7406006e8d2SBarry 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
7416006e8d2SBarry Smith      (the first block) would produce 0 and the mapping applied to 1 (the second block) would produce 3.
7426006e8d2SBarry Smith 
74345b6f7e9SBarry Smith    Level: advanced
74445b6f7e9SBarry Smith 
74545b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
74645b6f7e9SBarry Smith           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
74745b6f7e9SBarry Smith           AOPetscToApplication(), ISGlobalToLocalMappingApply()
74845b6f7e9SBarry Smith 
74945b6f7e9SBarry Smith @*/
75045b6f7e9SBarry Smith PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
75145b6f7e9SBarry Smith {
752cbc1caf0SMatthew G. Knepley 
753cbc1caf0SMatthew G. Knepley   PetscFunctionBegin;
754cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
755cbc1caf0SMatthew G. Knepley   {
756afcb2eb5SJed Brown     PetscInt i,Nmax = mapping->n;
757afcb2eb5SJed Brown     const PetscInt *idx = mapping->indices;
758d4bb536fSBarry Smith 
759afcb2eb5SJed Brown     for (i=0; i<N; i++) {
760afcb2eb5SJed Brown       if (in[i] < 0) {
761afcb2eb5SJed Brown         out[i] = in[i];
762afcb2eb5SJed Brown         continue;
763afcb2eb5SJed Brown       }
764e24637baSBarry 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);
765afcb2eb5SJed Brown       out[i] = idx[in[i]];
766afcb2eb5SJed Brown     }
767cbc1caf0SMatthew G. Knepley   }
768afcb2eb5SJed Brown   PetscFunctionReturn(0);
769afcb2eb5SJed Brown }
770d4bb536fSBarry Smith 
7717e99dc12SLawrence Mitchell /*@
772a997ad1aSLois Curfman McInnes     ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers
773a997ad1aSLois Curfman McInnes     specified with a global numbering.
774d4bb536fSBarry Smith 
775b9cd556bSLois Curfman McInnes     Not collective
776b9cd556bSLois Curfman McInnes 
777d4bb536fSBarry Smith     Input Parameters:
778b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
779d4bb536fSBarry Smith .   type - IS_GTOLM_MASK - replaces global indices with no local value with -1
780d4bb536fSBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
781d4bb536fSBarry Smith .   n - number of global indices to map
782b9cd556bSLois Curfman McInnes -   idx - global indices to map
783d4bb536fSBarry Smith 
784d4bb536fSBarry Smith     Output Parameters:
785b9cd556bSLois Curfman McInnes +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
786b9cd556bSLois Curfman McInnes -   idxout - local index of each global index, one must pass in an array long enough
787e182c471SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApply() with
7880298fd71SBarry Smith              idxout == NULL to determine the required length (returned in nout)
789e182c471SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApply()
790e182c471SBarry Smith              a second time to set the values.
791d4bb536fSBarry Smith 
792b9cd556bSLois Curfman McInnes     Notes:
7930298fd71SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
794d4bb536fSBarry Smith 
7959a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
796413f72f0SBarry 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.
797413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
7980f5bd95cSBarry Smith 
799a997ad1aSLois Curfman McInnes     Level: advanced
800a997ad1aSLois Curfman McInnes 
80132fd6b96SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
80232fd6b96SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
80332fd6b96SBarry Smith 
8049d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),
805413f72f0SBarry Smith           ISLocalToGlobalMappingDestroy()
806d4bb536fSBarry Smith @*/
807413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
808d4bb536fSBarry Smith {
8099d90f715SBarry Smith   PetscErrorCode ierr;
8109d90f715SBarry Smith 
8119d90f715SBarry Smith   PetscFunctionBegin;
8129d90f715SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
813413f72f0SBarry Smith   if (!mapping->data) {
814413f72f0SBarry Smith     ierr = ISGlobalToLocalMappingSetUp(mapping);CHKERRQ(ierr);
8159d90f715SBarry Smith   }
816413f72f0SBarry Smith   ierr = (*mapping->ops->globaltolocalmappingapply)(mapping,type,n,idx,nout,idxout);CHKERRQ(ierr);
8179d90f715SBarry Smith   PetscFunctionReturn(0);
8189d90f715SBarry Smith }
8199d90f715SBarry Smith 
820d4fe737eSStefano Zampini /*@
821d4fe737eSStefano Zampini     ISGlobalToLocalMappingApplyIS - Creates from an IS in the global numbering
822d4fe737eSStefano Zampini     a new index set using the local numbering defined in an ISLocalToGlobalMapping
823d4fe737eSStefano Zampini     context.
824d4fe737eSStefano Zampini 
825d4fe737eSStefano Zampini     Not collective
826d4fe737eSStefano Zampini 
827d4fe737eSStefano Zampini     Input Parameters:
828d4fe737eSStefano Zampini +   mapping - mapping between local and global numbering
8292785b321SStefano Zampini .   type - IS_GTOLM_MASK - replaces global indices with no local value with -1
8302785b321SStefano Zampini            IS_GTOLM_DROP - drops the indices with no local value from the output list
831d4fe737eSStefano Zampini -   is - index set in global numbering
832d4fe737eSStefano Zampini 
833d4fe737eSStefano Zampini     Output Parameters:
834d4fe737eSStefano Zampini .   newis - index set in local numbering
835d4fe737eSStefano Zampini 
8364cb36875SStefano Zampini     Notes:
8374cb36875SStefano Zampini     The output IS will be sequential, as it encodes a purely local operation
8384cb36875SStefano Zampini 
839d4fe737eSStefano Zampini     Level: advanced
840d4fe737eSStefano Zampini 
841d4fe737eSStefano Zampini .seealso: ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(),
842d4fe737eSStefano Zampini           ISLocalToGlobalMappingDestroy()
843d4fe737eSStefano Zampini @*/
844413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,IS is,IS *newis)
845d4fe737eSStefano Zampini {
846d4fe737eSStefano Zampini   PetscErrorCode ierr;
847d4fe737eSStefano Zampini   PetscInt       n,nout,*idxout;
848d4fe737eSStefano Zampini   const PetscInt *idxin;
849d4fe737eSStefano Zampini 
850d4fe737eSStefano Zampini   PetscFunctionBegin;
851d4fe737eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
852d4fe737eSStefano Zampini   PetscValidHeaderSpecific(is,IS_CLASSID,3);
853d4fe737eSStefano Zampini   PetscValidPointer(newis,4);
854d4fe737eSStefano Zampini 
855d4fe737eSStefano Zampini   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
856d4fe737eSStefano Zampini   ierr = ISGetIndices(is,&idxin);CHKERRQ(ierr);
857d4fe737eSStefano Zampini   if (type == IS_GTOLM_MASK) {
858d4fe737eSStefano Zampini     ierr = PetscMalloc1(n,&idxout);CHKERRQ(ierr);
859d4fe737eSStefano Zampini   } else {
860d4fe737eSStefano Zampini     ierr = ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,NULL);CHKERRQ(ierr);
861d4fe737eSStefano Zampini     ierr = PetscMalloc1(nout,&idxout);CHKERRQ(ierr);
862d4fe737eSStefano Zampini   }
863d4fe737eSStefano Zampini   ierr = ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,idxout);CHKERRQ(ierr);
864d4fe737eSStefano Zampini   ierr = ISRestoreIndices(is,&idxin);CHKERRQ(ierr);
8654cb36875SStefano Zampini   ierr = ISCreateGeneral(PETSC_COMM_SELF,nout,idxout,PETSC_OWN_POINTER,newis);CHKERRQ(ierr);
866d4fe737eSStefano Zampini   PetscFunctionReturn(0);
867d4fe737eSStefano Zampini }
868d4fe737eSStefano Zampini 
8699d90f715SBarry Smith /*@
8709d90f715SBarry Smith     ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers
8719d90f715SBarry Smith     specified with a block global numbering.
8729d90f715SBarry Smith 
8739d90f715SBarry Smith     Not collective
8749d90f715SBarry Smith 
8759d90f715SBarry Smith     Input Parameters:
8769d90f715SBarry Smith +   mapping - mapping between local and global numbering
8779d90f715SBarry Smith .   type - IS_GTOLM_MASK - replaces global indices with no local value with -1
8789d90f715SBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
8799d90f715SBarry Smith .   n - number of global indices to map
8809d90f715SBarry Smith -   idx - global indices to map
8819d90f715SBarry Smith 
8829d90f715SBarry Smith     Output Parameters:
8839d90f715SBarry Smith +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
8849d90f715SBarry Smith -   idxout - local index of each global index, one must pass in an array long enough
8859d90f715SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApplyBlock() with
8869d90f715SBarry Smith              idxout == NULL to determine the required length (returned in nout)
8879d90f715SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApplyBlock()
8889d90f715SBarry Smith              a second time to set the values.
8899d90f715SBarry Smith 
8909d90f715SBarry Smith     Notes:
8919d90f715SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
8929d90f715SBarry Smith 
8939a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
894413f72f0SBarry 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.
895413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
8969d90f715SBarry Smith 
8979d90f715SBarry Smith     Level: advanced
8989d90f715SBarry Smith 
8999d90f715SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
9009d90f715SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
9019d90f715SBarry Smith 
9029d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(),
903413f72f0SBarry Smith           ISLocalToGlobalMappingDestroy()
9049d90f715SBarry Smith @*/
905413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,
9069d90f715SBarry Smith                                                  PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
9079d90f715SBarry Smith {
9086849ba73SBarry Smith   PetscErrorCode ierr;
909d4bb536fSBarry Smith 
9103a40ed3dSBarry Smith   PetscFunctionBegin;
9110700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
912413f72f0SBarry Smith   if (!mapping->data) {
913413f72f0SBarry Smith     ierr = ISGlobalToLocalMappingSetUp(mapping);CHKERRQ(ierr);
914d4bb536fSBarry Smith   }
915413f72f0SBarry Smith   ierr = (*mapping->ops->globaltolocalmappingapplyblock)(mapping,type,n,idx,nout,idxout);CHKERRQ(ierr);
9163a40ed3dSBarry Smith   PetscFunctionReturn(0);
917d4bb536fSBarry Smith }
91890f02eecSBarry Smith 
919413f72f0SBarry Smith 
92089d82c54SBarry Smith /*@C
9216a818285SBarry Smith     ISLocalToGlobalMappingGetBlockInfo - Gets the neighbor information for each processor and
92289d82c54SBarry Smith      each index shared by more than one processor
92389d82c54SBarry Smith 
92489d82c54SBarry Smith     Collective on ISLocalToGlobalMapping
92589d82c54SBarry Smith 
92689d82c54SBarry Smith     Input Parameters:
92789d82c54SBarry Smith .   mapping - the mapping from local to global indexing
92889d82c54SBarry Smith 
92989d82c54SBarry Smith     Output Parameter:
93089d82c54SBarry Smith +   nproc - number of processors that are connected to this one
93189d82c54SBarry Smith .   proc - neighboring processors
93207b52d57SBarry Smith .   numproc - number of indices for each subdomain (processor)
9333463a7baSJed Brown -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
93489d82c54SBarry Smith 
93589d82c54SBarry Smith     Level: advanced
93689d82c54SBarry Smith 
9372cfcea29SBarry Smith     Fortran Usage:
9382cfcea29SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
9392cfcea29SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
9402cfcea29SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
9412cfcea29SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
9422cfcea29SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
9432cfcea29SBarry Smith 
9442cfcea29SBarry Smith 
94507b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
94607b52d57SBarry Smith           ISLocalToGlobalMappingRestoreInfo()
94789d82c54SBarry Smith @*/
9486a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
94989d82c54SBarry Smith {
9506849ba73SBarry Smith   PetscErrorCode ierr;
951268a049cSStefano Zampini 
952268a049cSStefano Zampini   PetscFunctionBegin;
953268a049cSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
954268a049cSStefano Zampini   if (mapping->info_cached) {
955268a049cSStefano Zampini     *nproc    = mapping->info_nproc;
956268a049cSStefano Zampini     *procs    = mapping->info_procs;
957268a049cSStefano Zampini     *numprocs = mapping->info_numprocs;
958268a049cSStefano Zampini     *indices  = mapping->info_indices;
959268a049cSStefano Zampini   } else {
960268a049cSStefano Zampini     ierr = ISLocalToGlobalMappingGetBlockInfo_Private(mapping,nproc,procs,numprocs,indices);CHKERRQ(ierr);
961268a049cSStefano Zampini   }
962268a049cSStefano Zampini   PetscFunctionReturn(0);
963268a049cSStefano Zampini }
964268a049cSStefano Zampini 
965268a049cSStefano Zampini static PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
966268a049cSStefano Zampini {
967268a049cSStefano Zampini   PetscErrorCode ierr;
96897f1f81fSBarry Smith   PetscMPIInt    size,rank,tag1,tag2,tag3,*len,*source,imdex;
96932dcc486SBarry Smith   PetscInt       i,n = mapping->n,Ng,ng,max = 0,*lindices = mapping->indices;
97032dcc486SBarry Smith   PetscInt       *nprocs,*owner,nsends,*sends,j,*starts,nmax,nrecvs,*recvs,proc;
97197f1f81fSBarry Smith   PetscInt       cnt,scale,*ownedsenders,*nownedsenders,rstart,nowned;
97232dcc486SBarry Smith   PetscInt       node,nownedm,nt,*sends2,nsends2,*starts2,*lens2,*dest,nrecvs2,*starts3,*recvs2,k,*bprocs,*tmp;
97332dcc486SBarry Smith   PetscInt       first_procs,first_numprocs,*first_indices;
97489d82c54SBarry Smith   MPI_Request    *recv_waits,*send_waits;
97530dcb7c9SBarry Smith   MPI_Status     recv_status,*send_status,*recv_statuses;
976ce94432eSBarry Smith   MPI_Comm       comm;
977ace3abfcSBarry Smith   PetscBool      debug = PETSC_FALSE;
97889d82c54SBarry Smith 
97989d82c54SBarry Smith   PetscFunctionBegin;
980ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)mapping,&comm);CHKERRQ(ierr);
98124cf384cSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
98224cf384cSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
98324cf384cSBarry Smith   if (size == 1) {
98424cf384cSBarry Smith     *nproc         = 0;
9850298fd71SBarry Smith     *procs         = NULL;
98695dccacaSBarry Smith     ierr           = PetscNew(numprocs);CHKERRQ(ierr);
9871e2105dcSBarry Smith     (*numprocs)[0] = 0;
98895dccacaSBarry Smith     ierr           = PetscNew(indices);CHKERRQ(ierr);
9890298fd71SBarry Smith     (*indices)[0]  = NULL;
990268a049cSStefano Zampini     /* save info for reuse */
991268a049cSStefano Zampini     mapping->info_nproc = *nproc;
992268a049cSStefano Zampini     mapping->info_procs = *procs;
993268a049cSStefano Zampini     mapping->info_numprocs = *numprocs;
994268a049cSStefano Zampini     mapping->info_indices = *indices;
995268a049cSStefano Zampini     mapping->info_cached = PETSC_TRUE;
99624cf384cSBarry Smith     PetscFunctionReturn(0);
99724cf384cSBarry Smith   }
99824cf384cSBarry Smith 
999c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)mapping)->options,NULL,"-islocaltoglobalmappinggetinfo_debug",&debug,NULL);CHKERRQ(ierr);
100007b52d57SBarry Smith 
10013677ff5aSBarry Smith   /*
10026a818285SBarry Smith     Notes on ISLocalToGlobalMappingGetBlockInfo
10033677ff5aSBarry Smith 
10043677ff5aSBarry Smith     globally owned node - the nodes that have been assigned to this processor in global
10053677ff5aSBarry Smith            numbering, just for this routine.
10063677ff5aSBarry Smith 
10073677ff5aSBarry Smith     nontrivial globally owned node - node assigned to this processor that is on a subdomain
10083677ff5aSBarry Smith            boundary (i.e. is has more than one local owner)
10093677ff5aSBarry Smith 
10103677ff5aSBarry Smith     locally owned node - node that exists on this processors subdomain
10113677ff5aSBarry Smith 
10123677ff5aSBarry Smith     nontrivial locally owned node - node that is not in the interior (i.e. has more than one
10133677ff5aSBarry Smith            local subdomain
10143677ff5aSBarry Smith   */
101524cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag1);CHKERRQ(ierr);
101624cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag2);CHKERRQ(ierr);
101724cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag3);CHKERRQ(ierr);
101889d82c54SBarry Smith 
101989d82c54SBarry Smith   for (i=0; i<n; i++) {
102089d82c54SBarry Smith     if (lindices[i] > max) max = lindices[i];
102189d82c54SBarry Smith   }
1022b2566f29SBarry Smith   ierr   = MPIU_Allreduce(&max,&Ng,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr);
102378058e43SBarry Smith   Ng++;
102489d82c54SBarry Smith   ierr   = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
102589d82c54SBarry Smith   ierr   = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1026bc8ff85bSBarry Smith   scale  = Ng/size + 1;
1027a2e34c3dSBarry Smith   ng     = scale; if (rank == size-1) ng = Ng - scale*(size-1); ng = PetscMax(1,ng);
1028caba0dd0SBarry Smith   rstart = scale*rank;
102989d82c54SBarry Smith 
103089d82c54SBarry Smith   /* determine ownership ranges of global indices */
1031785e854fSJed Brown   ierr = PetscMalloc1(2*size,&nprocs);CHKERRQ(ierr);
1032580bdb30SBarry Smith   ierr = PetscArrayzero(nprocs,2*size);CHKERRQ(ierr);
103389d82c54SBarry Smith 
103489d82c54SBarry Smith   /* determine owners of each local node  */
1035785e854fSJed Brown   ierr = PetscMalloc1(n,&owner);CHKERRQ(ierr);
103689d82c54SBarry Smith   for (i=0; i<n; i++) {
10373677ff5aSBarry Smith     proc             = lindices[i]/scale; /* processor that globally owns this index */
103827c402fcSBarry Smith     nprocs[2*proc+1] = 1;                 /* processor globally owns at least one of ours */
10393677ff5aSBarry Smith     owner[i]         = proc;
104027c402fcSBarry Smith     nprocs[2*proc]++;                     /* count of how many that processor globally owns of ours */
104189d82c54SBarry Smith   }
104227c402fcSBarry Smith   nsends = 0; for (i=0; i<size; i++) nsends += nprocs[2*i+1];
10437904a332SBarry Smith   ierr = PetscInfo1(mapping,"Number of global owners for my local data %D\n",nsends);CHKERRQ(ierr);
104489d82c54SBarry Smith 
104589d82c54SBarry Smith   /* inform other processors of number of messages and max length*/
104627c402fcSBarry Smith   ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr);
10477904a332SBarry Smith   ierr = PetscInfo1(mapping,"Number of local owners for my global data %D\n",nrecvs);CHKERRQ(ierr);
104889d82c54SBarry Smith 
104989d82c54SBarry Smith   /* post receives for owned rows */
1050785e854fSJed Brown   ierr = PetscMalloc1((2*nrecvs+1)*(nmax+1),&recvs);CHKERRQ(ierr);
1051854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&recv_waits);CHKERRQ(ierr);
105289d82c54SBarry Smith   for (i=0; i<nrecvs; i++) {
105332dcc486SBarry Smith     ierr = MPI_Irecv(recvs+2*nmax*i,2*nmax,MPIU_INT,MPI_ANY_SOURCE,tag1,comm,recv_waits+i);CHKERRQ(ierr);
105489d82c54SBarry Smith   }
105589d82c54SBarry Smith 
105689d82c54SBarry Smith   /* pack messages containing lists of local nodes to owners */
1057854ce69bSBarry Smith   ierr      = PetscMalloc1(2*n+1,&sends);CHKERRQ(ierr);
1058854ce69bSBarry Smith   ierr      = PetscMalloc1(size+1,&starts);CHKERRQ(ierr);
105989d82c54SBarry Smith   starts[0] = 0;
1060f6e5521dSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
106189d82c54SBarry Smith   for (i=0; i<n; i++) {
106289d82c54SBarry Smith     sends[starts[owner[i]]++] = lindices[i];
106330dcb7c9SBarry Smith     sends[starts[owner[i]]++] = i;
106489d82c54SBarry Smith   }
106589d82c54SBarry Smith   ierr = PetscFree(owner);CHKERRQ(ierr);
106689d82c54SBarry Smith   starts[0] = 0;
1067f6e5521dSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
106889d82c54SBarry Smith 
106989d82c54SBarry Smith   /* send the messages */
1070854ce69bSBarry Smith   ierr = PetscMalloc1(nsends+1,&send_waits);CHKERRQ(ierr);
1071854ce69bSBarry Smith   ierr = PetscMalloc1(nsends+1,&dest);CHKERRQ(ierr);
107289d82c54SBarry Smith   cnt = 0;
107389d82c54SBarry Smith   for (i=0; i<size; i++) {
107427c402fcSBarry Smith     if (nprocs[2*i]) {
107532dcc486SBarry Smith       ierr      = MPI_Isend(sends+starts[i],2*nprocs[2*i],MPIU_INT,i,tag1,comm,send_waits+cnt);CHKERRQ(ierr);
107630dcb7c9SBarry Smith       dest[cnt] = i;
107789d82c54SBarry Smith       cnt++;
107889d82c54SBarry Smith     }
107989d82c54SBarry Smith   }
108089d82c54SBarry Smith   ierr = PetscFree(starts);CHKERRQ(ierr);
108189d82c54SBarry Smith 
108289d82c54SBarry Smith   /* wait on receives */
1083854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&source);CHKERRQ(ierr);
1084854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&len);CHKERRQ(ierr);
108589d82c54SBarry Smith   cnt  = nrecvs;
1086580bdb30SBarry Smith   ierr = PetscCalloc1(ng+1,&nownedsenders);CHKERRQ(ierr);
108789d82c54SBarry Smith   while (cnt) {
108889d82c54SBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr);
108989d82c54SBarry Smith     /* unpack receives into our local space */
109032dcc486SBarry Smith     ierr          = MPI_Get_count(&recv_status,MPIU_INT,&len[imdex]);CHKERRQ(ierr);
109189d82c54SBarry Smith     source[imdex] = recv_status.MPI_SOURCE;
109230dcb7c9SBarry Smith     len[imdex]    = len[imdex]/2;
1093caba0dd0SBarry Smith     /* count how many local owners for each of my global owned indices */
109430dcb7c9SBarry Smith     for (i=0; i<len[imdex]; i++) nownedsenders[recvs[2*imdex*nmax+2*i]-rstart]++;
109589d82c54SBarry Smith     cnt--;
109689d82c54SBarry Smith   }
109789d82c54SBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
109889d82c54SBarry Smith 
109930dcb7c9SBarry Smith   /* count how many globally owned indices are on an edge multiplied by how many processors own them. */
1100bc8ff85bSBarry Smith   nowned  = 0;
1101bc8ff85bSBarry Smith   nownedm = 0;
1102bc8ff85bSBarry Smith   for (i=0; i<ng; i++) {
1103bc8ff85bSBarry Smith     if (nownedsenders[i] > 1) {nownedm += nownedsenders[i]; nowned++;}
1104bc8ff85bSBarry Smith   }
1105bc8ff85bSBarry Smith 
1106bc8ff85bSBarry Smith   /* create single array to contain rank of all local owners of each globally owned index */
1107854ce69bSBarry Smith   ierr      = PetscMalloc1(nownedm+1,&ownedsenders);CHKERRQ(ierr);
1108854ce69bSBarry Smith   ierr      = PetscMalloc1(ng+1,&starts);CHKERRQ(ierr);
1109bc8ff85bSBarry Smith   starts[0] = 0;
1110bc8ff85bSBarry Smith   for (i=1; i<ng; i++) {
1111bc8ff85bSBarry Smith     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
1112bc8ff85bSBarry Smith     else starts[i] = starts[i-1];
1113bc8ff85bSBarry Smith   }
1114bc8ff85bSBarry Smith 
111530dcb7c9SBarry Smith   /* for each nontrival globally owned node list all arriving processors */
1116bc8ff85bSBarry Smith   for (i=0; i<nrecvs; i++) {
1117bc8ff85bSBarry Smith     for (j=0; j<len[i]; j++) {
111830dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
1119f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) ownedsenders[starts[node]++] = source[i];
1120bc8ff85bSBarry Smith     }
1121bc8ff85bSBarry Smith   }
1122bc8ff85bSBarry Smith 
112307b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
112430dcb7c9SBarry Smith     starts[0] = 0;
112530dcb7c9SBarry Smith     for (i=1; i<ng; i++) {
112630dcb7c9SBarry Smith       if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
112730dcb7c9SBarry Smith       else starts[i] = starts[i-1];
112830dcb7c9SBarry Smith     }
112930dcb7c9SBarry Smith     for (i=0; i<ng; i++) {
113030dcb7c9SBarry Smith       if (nownedsenders[i] > 1) {
11317904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"[%d] global node %D local owner processors: ",rank,i+rstart);CHKERRQ(ierr);
113230dcb7c9SBarry Smith         for (j=0; j<nownedsenders[i]; j++) {
11337904a332SBarry Smith           ierr = PetscSynchronizedPrintf(comm,"%D ",ownedsenders[starts[i]+j]);CHKERRQ(ierr);
113430dcb7c9SBarry Smith         }
113530dcb7c9SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
113630dcb7c9SBarry Smith       }
113730dcb7c9SBarry Smith     }
11380ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
113907b52d57SBarry Smith   } /* -----------------------------------  */
114030dcb7c9SBarry Smith 
11413677ff5aSBarry Smith   /* wait on original sends */
11423a96401aSBarry Smith   if (nsends) {
1143785e854fSJed Brown     ierr = PetscMalloc1(nsends,&send_status);CHKERRQ(ierr);
11443a96401aSBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr);
11453a96401aSBarry Smith     ierr = PetscFree(send_status);CHKERRQ(ierr);
11463a96401aSBarry Smith   }
114789d82c54SBarry Smith   ierr = PetscFree(send_waits);CHKERRQ(ierr);
11483a96401aSBarry Smith   ierr = PetscFree(sends);CHKERRQ(ierr);
11493677ff5aSBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
11503677ff5aSBarry Smith 
11513677ff5aSBarry Smith   /* pack messages to send back to local owners */
115230dcb7c9SBarry Smith   starts[0] = 0;
115330dcb7c9SBarry Smith   for (i=1; i<ng; i++) {
115430dcb7c9SBarry Smith     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
115530dcb7c9SBarry Smith     else starts[i] = starts[i-1];
115630dcb7c9SBarry Smith   }
115730dcb7c9SBarry Smith   nsends2 = nrecvs;
1158854ce69bSBarry Smith   ierr    = PetscMalloc1(nsends2+1,&nprocs);CHKERRQ(ierr); /* length of each message */
115930dcb7c9SBarry Smith   for (i=0; i<nrecvs; i++) {
116030dcb7c9SBarry Smith     nprocs[i] = 1;
116130dcb7c9SBarry Smith     for (j=0; j<len[i]; j++) {
116230dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
1163f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) nprocs[i] += 2 + nownedsenders[node];
116430dcb7c9SBarry Smith     }
116530dcb7c9SBarry Smith   }
1166f6e5521dSKarl Rupp   nt = 0;
1167f6e5521dSKarl Rupp   for (i=0; i<nsends2; i++) nt += nprocs[i];
1168f6e5521dSKarl Rupp 
1169854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,&sends2);CHKERRQ(ierr);
1170854ce69bSBarry Smith   ierr = PetscMalloc1(nsends2+1,&starts2);CHKERRQ(ierr);
1171f6e5521dSKarl Rupp 
1172f6e5521dSKarl Rupp   starts2[0] = 0;
1173f6e5521dSKarl Rupp   for (i=1; i<nsends2; i++) starts2[i] = starts2[i-1] + nprocs[i-1];
117430dcb7c9SBarry Smith   /*
117530dcb7c9SBarry Smith      Each message is 1 + nprocs[i] long, and consists of
117630dcb7c9SBarry Smith        (0) the number of nodes being sent back
117730dcb7c9SBarry Smith        (1) the local node number,
117830dcb7c9SBarry Smith        (2) the number of processors sharing it,
117930dcb7c9SBarry Smith        (3) the processors sharing it
118030dcb7c9SBarry Smith   */
118130dcb7c9SBarry Smith   for (i=0; i<nsends2; i++) {
118230dcb7c9SBarry Smith     cnt = 1;
118330dcb7c9SBarry Smith     sends2[starts2[i]] = 0;
118430dcb7c9SBarry Smith     for (j=0; j<len[i]; j++) {
118530dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
118630dcb7c9SBarry Smith       if (nownedsenders[node] > 1) {
118730dcb7c9SBarry Smith         sends2[starts2[i]]++;
118830dcb7c9SBarry Smith         sends2[starts2[i]+cnt++] = recvs[2*i*nmax+2*j+1];
118930dcb7c9SBarry Smith         sends2[starts2[i]+cnt++] = nownedsenders[node];
1190580bdb30SBarry Smith         ierr = PetscArraycpy(&sends2[starts2[i]+cnt],&ownedsenders[starts[node]],nownedsenders[node]);CHKERRQ(ierr);
119130dcb7c9SBarry Smith         cnt += nownedsenders[node];
119230dcb7c9SBarry Smith       }
119330dcb7c9SBarry Smith     }
119430dcb7c9SBarry Smith   }
119530dcb7c9SBarry Smith 
119630dcb7c9SBarry Smith   /* receive the message lengths */
119730dcb7c9SBarry Smith   nrecvs2 = nsends;
1198854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&lens2);CHKERRQ(ierr);
1199854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&starts3);CHKERRQ(ierr);
1200854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&recv_waits);CHKERRQ(ierr);
120130dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
1202d44834fbSBarry Smith     ierr = MPI_Irecv(&lens2[i],1,MPIU_INT,dest[i],tag2,comm,recv_waits+i);CHKERRQ(ierr);
120330dcb7c9SBarry Smith   }
1204d44834fbSBarry Smith 
12058a8e0b3aSBarry Smith   /* send the message lengths */
12068a8e0b3aSBarry Smith   for (i=0; i<nsends2; i++) {
12078a8e0b3aSBarry Smith     ierr = MPI_Send(&nprocs[i],1,MPIU_INT,source[i],tag2,comm);CHKERRQ(ierr);
12088a8e0b3aSBarry Smith   }
12098a8e0b3aSBarry Smith 
1210d44834fbSBarry Smith   /* wait on receives of lens */
12110c468ba9SBarry Smith   if (nrecvs2) {
1212785e854fSJed Brown     ierr = PetscMalloc1(nrecvs2,&recv_statuses);CHKERRQ(ierr);
1213d44834fbSBarry Smith     ierr = MPI_Waitall(nrecvs2,recv_waits,recv_statuses);CHKERRQ(ierr);
1214d44834fbSBarry Smith     ierr = PetscFree(recv_statuses);CHKERRQ(ierr);
12150c468ba9SBarry Smith   }
1216a2ea699eSBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
1217d44834fbSBarry Smith 
121830dcb7c9SBarry Smith   starts3[0] = 0;
1219d44834fbSBarry Smith   nt         = 0;
122030dcb7c9SBarry Smith   for (i=0; i<nrecvs2-1; i++) {
122130dcb7c9SBarry Smith     starts3[i+1] = starts3[i] + lens2[i];
1222d44834fbSBarry Smith     nt          += lens2[i];
122330dcb7c9SBarry Smith   }
122476466f69SStefano Zampini   if (nrecvs2) nt += lens2[nrecvs2-1];
1225d44834fbSBarry Smith 
1226854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,&recvs2);CHKERRQ(ierr);
1227854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs2+1,&recv_waits);CHKERRQ(ierr);
122852b72c4aSBarry Smith   for (i=0; i<nrecvs2; i++) {
122932dcc486SBarry Smith     ierr = MPI_Irecv(recvs2+starts3[i],lens2[i],MPIU_INT,dest[i],tag3,comm,recv_waits+i);CHKERRQ(ierr);
123030dcb7c9SBarry Smith   }
123130dcb7c9SBarry Smith 
123230dcb7c9SBarry Smith   /* send the messages */
1233854ce69bSBarry Smith   ierr = PetscMalloc1(nsends2+1,&send_waits);CHKERRQ(ierr);
123430dcb7c9SBarry Smith   for (i=0; i<nsends2; i++) {
123532dcc486SBarry Smith     ierr = MPI_Isend(sends2+starts2[i],nprocs[i],MPIU_INT,source[i],tag3,comm,send_waits+i);CHKERRQ(ierr);
123630dcb7c9SBarry Smith   }
123730dcb7c9SBarry Smith 
123830dcb7c9SBarry Smith   /* wait on receives */
12390c468ba9SBarry Smith   if (nrecvs2) {
1240785e854fSJed Brown     ierr = PetscMalloc1(nrecvs2,&recv_statuses);CHKERRQ(ierr);
124130dcb7c9SBarry Smith     ierr = MPI_Waitall(nrecvs2,recv_waits,recv_statuses);CHKERRQ(ierr);
124230dcb7c9SBarry Smith     ierr = PetscFree(recv_statuses);CHKERRQ(ierr);
12430c468ba9SBarry Smith   }
124430dcb7c9SBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
124530dcb7c9SBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
124630dcb7c9SBarry Smith 
124707b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
124830dcb7c9SBarry Smith     cnt = 0;
124930dcb7c9SBarry Smith     for (i=0; i<nrecvs2; i++) {
125030dcb7c9SBarry Smith       nt = recvs2[cnt++];
125130dcb7c9SBarry Smith       for (j=0; j<nt; j++) {
12527904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"[%d] local node %D number of subdomains %D: ",rank,recvs2[cnt],recvs2[cnt+1]);CHKERRQ(ierr);
125330dcb7c9SBarry Smith         for (k=0; k<recvs2[cnt+1]; k++) {
12547904a332SBarry Smith           ierr = PetscSynchronizedPrintf(comm,"%D ",recvs2[cnt+2+k]);CHKERRQ(ierr);
125530dcb7c9SBarry Smith         }
125630dcb7c9SBarry Smith         cnt += 2 + recvs2[cnt+1];
125730dcb7c9SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
125830dcb7c9SBarry Smith       }
125930dcb7c9SBarry Smith     }
12600ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
126107b52d57SBarry Smith   } /* -----------------------------------  */
126230dcb7c9SBarry Smith 
126330dcb7c9SBarry Smith   /* count number subdomains for each local node */
1264580bdb30SBarry Smith   ierr = PetscCalloc1(size,&nprocs);CHKERRQ(ierr);
126530dcb7c9SBarry Smith   cnt  = 0;
126630dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
126730dcb7c9SBarry Smith     nt = recvs2[cnt++];
126830dcb7c9SBarry Smith     for (j=0; j<nt; j++) {
1269f6e5521dSKarl Rupp       for (k=0; k<recvs2[cnt+1]; k++) nprocs[recvs2[cnt+2+k]]++;
127030dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt+1];
127130dcb7c9SBarry Smith     }
127230dcb7c9SBarry Smith   }
127330dcb7c9SBarry Smith   nt = 0; for (i=0; i<size; i++) nt += (nprocs[i] > 0);
127430dcb7c9SBarry Smith   *nproc    = nt;
1275854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,procs);CHKERRQ(ierr);
1276854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,numprocs);CHKERRQ(ierr);
1277854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,indices);CHKERRQ(ierr);
12780298fd71SBarry Smith   for (i=0;i<nt+1;i++) (*indices)[i]=NULL;
1279785e854fSJed Brown   ierr = PetscMalloc1(size,&bprocs);CHKERRQ(ierr);
128030dcb7c9SBarry Smith   cnt  = 0;
128130dcb7c9SBarry Smith   for (i=0; i<size; i++) {
128230dcb7c9SBarry Smith     if (nprocs[i] > 0) {
128330dcb7c9SBarry Smith       bprocs[i]        = cnt;
128430dcb7c9SBarry Smith       (*procs)[cnt]    = i;
128530dcb7c9SBarry Smith       (*numprocs)[cnt] = nprocs[i];
1286785e854fSJed Brown       ierr             = PetscMalloc1(nprocs[i],&(*indices)[cnt]);CHKERRQ(ierr);
128730dcb7c9SBarry Smith       cnt++;
128830dcb7c9SBarry Smith     }
128930dcb7c9SBarry Smith   }
129030dcb7c9SBarry Smith 
129130dcb7c9SBarry Smith   /* make the list of subdomains for each nontrivial local node */
1292580bdb30SBarry Smith   ierr = PetscArrayzero(*numprocs,nt);CHKERRQ(ierr);
129330dcb7c9SBarry Smith   cnt  = 0;
129430dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
129530dcb7c9SBarry Smith     nt = recvs2[cnt++];
129630dcb7c9SBarry Smith     for (j=0; j<nt; j++) {
1297f6e5521dSKarl Rupp       for (k=0; k<recvs2[cnt+1]; k++) (*indices)[bprocs[recvs2[cnt+2+k]]][(*numprocs)[bprocs[recvs2[cnt+2+k]]]++] = recvs2[cnt];
129830dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt+1];
129930dcb7c9SBarry Smith     }
130030dcb7c9SBarry Smith   }
130130dcb7c9SBarry Smith   ierr = PetscFree(bprocs);CHKERRQ(ierr);
130207b52d57SBarry Smith   ierr = PetscFree(recvs2);CHKERRQ(ierr);
130330dcb7c9SBarry Smith 
130407b52d57SBarry Smith   /* sort the node indexing by their global numbers */
130507b52d57SBarry Smith   nt = *nproc;
130607b52d57SBarry Smith   for (i=0; i<nt; i++) {
1307854ce69bSBarry Smith     ierr = PetscMalloc1((*numprocs)[i],&tmp);CHKERRQ(ierr);
1308f6e5521dSKarl Rupp     for (j=0; j<(*numprocs)[i]; j++) tmp[j] = lindices[(*indices)[i][j]];
130907b52d57SBarry Smith     ierr = PetscSortIntWithArray((*numprocs)[i],tmp,(*indices)[i]);CHKERRQ(ierr);
131007b52d57SBarry Smith     ierr = PetscFree(tmp);CHKERRQ(ierr);
131107b52d57SBarry Smith   }
131207b52d57SBarry Smith 
131307b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
131430dcb7c9SBarry Smith     nt = *nproc;
131530dcb7c9SBarry Smith     for (i=0; i<nt; i++) {
13167904a332SBarry Smith       ierr = PetscSynchronizedPrintf(comm,"[%d] subdomain %D number of indices %D: ",rank,(*procs)[i],(*numprocs)[i]);CHKERRQ(ierr);
131730dcb7c9SBarry Smith       for (j=0; j<(*numprocs)[i]; j++) {
13187904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"%D ",(*indices)[i][j]);CHKERRQ(ierr);
131930dcb7c9SBarry Smith       }
132030dcb7c9SBarry Smith       ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
132130dcb7c9SBarry Smith     }
13220ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
132307b52d57SBarry Smith   } /* -----------------------------------  */
132430dcb7c9SBarry Smith 
132530dcb7c9SBarry Smith   /* wait on sends */
132630dcb7c9SBarry Smith   if (nsends2) {
1327785e854fSJed Brown     ierr = PetscMalloc1(nsends2,&send_status);CHKERRQ(ierr);
132830dcb7c9SBarry Smith     ierr = MPI_Waitall(nsends2,send_waits,send_status);CHKERRQ(ierr);
132930dcb7c9SBarry Smith     ierr = PetscFree(send_status);CHKERRQ(ierr);
133030dcb7c9SBarry Smith   }
133130dcb7c9SBarry Smith 
133230dcb7c9SBarry Smith   ierr = PetscFree(starts3);CHKERRQ(ierr);
133330dcb7c9SBarry Smith   ierr = PetscFree(dest);CHKERRQ(ierr);
133430dcb7c9SBarry Smith   ierr = PetscFree(send_waits);CHKERRQ(ierr);
13353677ff5aSBarry Smith 
1336bc8ff85bSBarry Smith   ierr = PetscFree(nownedsenders);CHKERRQ(ierr);
1337bc8ff85bSBarry Smith   ierr = PetscFree(ownedsenders);CHKERRQ(ierr);
1338bc8ff85bSBarry Smith   ierr = PetscFree(starts);CHKERRQ(ierr);
133930dcb7c9SBarry Smith   ierr = PetscFree(starts2);CHKERRQ(ierr);
134030dcb7c9SBarry Smith   ierr = PetscFree(lens2);CHKERRQ(ierr);
134189d82c54SBarry Smith 
134289d82c54SBarry Smith   ierr = PetscFree(source);CHKERRQ(ierr);
134397f1f81fSBarry Smith   ierr = PetscFree(len);CHKERRQ(ierr);
134489d82c54SBarry Smith   ierr = PetscFree(recvs);CHKERRQ(ierr);
13453a96401aSBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
134630dcb7c9SBarry Smith   ierr = PetscFree(sends2);CHKERRQ(ierr);
134724cf384cSBarry Smith 
134824cf384cSBarry Smith   /* put the information about myself as the first entry in the list */
134924cf384cSBarry Smith   first_procs    = (*procs)[0];
135024cf384cSBarry Smith   first_numprocs = (*numprocs)[0];
135124cf384cSBarry Smith   first_indices  = (*indices)[0];
135224cf384cSBarry Smith   for (i=0; i<*nproc; i++) {
135324cf384cSBarry Smith     if ((*procs)[i] == rank) {
135424cf384cSBarry Smith       (*procs)[0]    = (*procs)[i];
135524cf384cSBarry Smith       (*numprocs)[0] = (*numprocs)[i];
135624cf384cSBarry Smith       (*indices)[0]  = (*indices)[i];
135724cf384cSBarry Smith       (*procs)[i]    = first_procs;
135824cf384cSBarry Smith       (*numprocs)[i] = first_numprocs;
135924cf384cSBarry Smith       (*indices)[i]  = first_indices;
136024cf384cSBarry Smith       break;
136124cf384cSBarry Smith     }
136224cf384cSBarry Smith   }
1363268a049cSStefano Zampini 
1364268a049cSStefano Zampini   /* save info for reuse */
1365268a049cSStefano Zampini   mapping->info_nproc = *nproc;
1366268a049cSStefano Zampini   mapping->info_procs = *procs;
1367268a049cSStefano Zampini   mapping->info_numprocs = *numprocs;
1368268a049cSStefano Zampini   mapping->info_indices = *indices;
1369268a049cSStefano Zampini   mapping->info_cached = PETSC_TRUE;
137089d82c54SBarry Smith   PetscFunctionReturn(0);
137189d82c54SBarry Smith }
137289d82c54SBarry Smith 
13736a818285SBarry Smith /*@C
13746a818285SBarry Smith     ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by ISLocalToGlobalMappingGetBlockInfo()
13756a818285SBarry Smith 
13766a818285SBarry Smith     Collective on ISLocalToGlobalMapping
13776a818285SBarry Smith 
13786a818285SBarry Smith     Input Parameters:
13796a818285SBarry Smith .   mapping - the mapping from local to global indexing
13806a818285SBarry Smith 
13816a818285SBarry Smith     Output Parameter:
13826a818285SBarry Smith +   nproc - number of processors that are connected to this one
13836a818285SBarry Smith .   proc - neighboring processors
13846a818285SBarry Smith .   numproc - number of indices for each processor
13856a818285SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
13866a818285SBarry Smith 
13876a818285SBarry Smith     Level: advanced
13886a818285SBarry Smith 
13896a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
13906a818285SBarry Smith           ISLocalToGlobalMappingGetInfo()
13916a818285SBarry Smith @*/
13926a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
13936a818285SBarry Smith {
13946a818285SBarry Smith   PetscErrorCode ierr;
13956a818285SBarry Smith 
13966a818285SBarry Smith   PetscFunctionBegin;
1397cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1398268a049cSStefano Zampini   if (mapping->info_free) {
13996a818285SBarry Smith     ierr = PetscFree(*numprocs);CHKERRQ(ierr);
14006a818285SBarry Smith     if (*indices) {
1401268a049cSStefano Zampini       PetscInt i;
1402268a049cSStefano Zampini 
14036a818285SBarry Smith       ierr = PetscFree((*indices)[0]);CHKERRQ(ierr);
14046a818285SBarry Smith       for (i=1; i<*nproc; i++) {
14056a818285SBarry Smith         ierr = PetscFree((*indices)[i]);CHKERRQ(ierr);
14066a818285SBarry Smith       }
14076a818285SBarry Smith       ierr = PetscFree(*indices);CHKERRQ(ierr);
14086a818285SBarry Smith     }
1409268a049cSStefano Zampini   }
1410268a049cSStefano Zampini   *nproc    = 0;
1411268a049cSStefano Zampini   *procs    = NULL;
1412268a049cSStefano Zampini   *numprocs = NULL;
1413268a049cSStefano Zampini   *indices  = NULL;
14146a818285SBarry Smith   PetscFunctionReturn(0);
14156a818285SBarry Smith }
14166a818285SBarry Smith 
14176a818285SBarry Smith /*@C
14186a818285SBarry Smith     ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each processor and
14196a818285SBarry Smith      each index shared by more than one processor
14206a818285SBarry Smith 
14216a818285SBarry Smith     Collective on ISLocalToGlobalMapping
14226a818285SBarry Smith 
14236a818285SBarry Smith     Input Parameters:
14246a818285SBarry Smith .   mapping - the mapping from local to global indexing
14256a818285SBarry Smith 
14266a818285SBarry Smith     Output Parameter:
14276a818285SBarry Smith +   nproc - number of processors that are connected to this one
14286a818285SBarry Smith .   proc - neighboring processors
14296a818285SBarry Smith .   numproc - number of indices for each subdomain (processor)
14306a818285SBarry Smith -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
14316a818285SBarry Smith 
14326a818285SBarry Smith     Level: advanced
14336a818285SBarry Smith 
14341bd0b88eSStefano Zampini     Notes: The user needs to call ISLocalToGlobalMappingRestoreInfo when the data is no longer needed.
14351bd0b88eSStefano Zampini 
14366a818285SBarry Smith     Fortran Usage:
14376a818285SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
14386a818285SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
14396a818285SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
14406a818285SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
14416a818285SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
14426a818285SBarry Smith 
14436a818285SBarry Smith 
14446a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
14456a818285SBarry Smith           ISLocalToGlobalMappingRestoreInfo()
14466a818285SBarry Smith @*/
14476a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
14486a818285SBarry Smith {
14496a818285SBarry Smith   PetscErrorCode ierr;
1450268a049cSStefano Zampini   PetscInt       **bindices = NULL,*bnumprocs = NULL,bs = mapping->bs,i,j,k;
14516a818285SBarry Smith 
14526a818285SBarry Smith   PetscFunctionBegin;
14536a818285SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1454268a049cSStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockInfo(mapping,nproc,procs,&bnumprocs,&bindices);CHKERRQ(ierr);
1455268a049cSStefano Zampini   if (bs > 1) { /* we need to expand the cached info */
1456732f65e3SBarry Smith     ierr = PetscCalloc1(*nproc,&*indices);CHKERRQ(ierr);
1457268a049cSStefano Zampini     ierr = PetscCalloc1(*nproc,&*numprocs);CHKERRQ(ierr);
14586a818285SBarry Smith     for (i=0; i<*nproc; i++) {
1459268a049cSStefano Zampini       ierr = PetscMalloc1(bs*bnumprocs[i],&(*indices)[i]);CHKERRQ(ierr);
1460268a049cSStefano Zampini       for (j=0; j<bnumprocs[i]; j++) {
14616a818285SBarry Smith         for (k=0; k<bs; k++) {
14626a818285SBarry Smith           (*indices)[i][j*bs+k] = bs*bindices[i][j] + k;
14636a818285SBarry Smith         }
14646a818285SBarry Smith       }
1465268a049cSStefano Zampini       (*numprocs)[i] = bnumprocs[i]*bs;
14666a818285SBarry Smith     }
1467268a049cSStefano Zampini     mapping->info_free = PETSC_TRUE;
1468268a049cSStefano Zampini   } else {
1469268a049cSStefano Zampini     *numprocs = bnumprocs;
1470268a049cSStefano Zampini     *indices  = bindices;
14716a818285SBarry Smith   }
14726a818285SBarry Smith   PetscFunctionReturn(0);
14736a818285SBarry Smith }
14746a818285SBarry Smith 
147507b52d57SBarry Smith /*@C
147607b52d57SBarry Smith     ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by ISLocalToGlobalMappingGetInfo()
147789d82c54SBarry Smith 
147807b52d57SBarry Smith     Collective on ISLocalToGlobalMapping
147907b52d57SBarry Smith 
148007b52d57SBarry Smith     Input Parameters:
148107b52d57SBarry Smith .   mapping - the mapping from local to global indexing
148207b52d57SBarry Smith 
148307b52d57SBarry Smith     Output Parameter:
148407b52d57SBarry Smith +   nproc - number of processors that are connected to this one
148507b52d57SBarry Smith .   proc - neighboring processors
148607b52d57SBarry Smith .   numproc - number of indices for each processor
148707b52d57SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
148807b52d57SBarry Smith 
148907b52d57SBarry Smith     Level: advanced
149007b52d57SBarry Smith 
149107b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
149207b52d57SBarry Smith           ISLocalToGlobalMappingGetInfo()
149307b52d57SBarry Smith @*/
14947087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
149507b52d57SBarry Smith {
14966849ba73SBarry Smith   PetscErrorCode ierr;
149707b52d57SBarry Smith 
149807b52d57SBarry Smith   PetscFunctionBegin;
14996a818285SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockInfo(mapping,nproc,procs,numprocs,indices);CHKERRQ(ierr);
150007b52d57SBarry Smith   PetscFunctionReturn(0);
150107b52d57SBarry Smith }
150286994e45SJed Brown 
150386994e45SJed Brown /*@C
15041bd0b88eSStefano Zampini     ISLocalToGlobalMappingGetNodeInfo - Gets the neighbor information for each node
15051bd0b88eSStefano Zampini 
15061bd0b88eSStefano Zampini     Collective on ISLocalToGlobalMapping
15071bd0b88eSStefano Zampini 
15081bd0b88eSStefano Zampini     Input Parameters:
15091bd0b88eSStefano Zampini .   mapping - the mapping from local to global indexing
15101bd0b88eSStefano Zampini 
15111bd0b88eSStefano Zampini     Output Parameter:
15121bd0b88eSStefano Zampini +   nnodes - number of local nodes (same ISLocalToGlobalMappingGetSize())
15131bd0b88eSStefano Zampini .   count - number of neighboring processors per node
15141bd0b88eSStefano Zampini -   indices - indices of processes sharing the node (sorted)
15151bd0b88eSStefano Zampini 
15161bd0b88eSStefano Zampini     Level: advanced
15171bd0b88eSStefano Zampini 
15181bd0b88eSStefano Zampini     Notes: The user needs to call ISLocalToGlobalMappingRestoreInfo when the data is no longer needed.
15191bd0b88eSStefano Zampini 
15201bd0b88eSStefano Zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
15211bd0b88eSStefano Zampini           ISLocalToGlobalMappingGetInfo(), ISLocalToGlobalMappingRestoreNodeInfo()
15221bd0b88eSStefano Zampini @*/
15231bd0b88eSStefano Zampini PetscErrorCode  ISLocalToGlobalMappingGetNodeInfo(ISLocalToGlobalMapping mapping,PetscInt *nnodes,PetscInt *count[],PetscInt **indices[])
15241bd0b88eSStefano Zampini {
15251bd0b88eSStefano Zampini   PetscInt       n;
15261bd0b88eSStefano Zampini   PetscErrorCode ierr;
15271bd0b88eSStefano Zampini 
15281bd0b88eSStefano Zampini   PetscFunctionBegin;
15291bd0b88eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
15301bd0b88eSStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(mapping,&n);CHKERRQ(ierr);
15311bd0b88eSStefano Zampini   if (!mapping->info_nodec) {
15321bd0b88eSStefano Zampini     PetscInt i,m,n_neigh,*neigh,*n_shared,**shared;
15331bd0b88eSStefano Zampini 
1534*071fcb05SBarry Smith     ierr = PetscMalloc2(n+1,&mapping->info_nodec,n,&mapping->info_nodei);CHKERRQ(ierr);
15351bd0b88eSStefano Zampini     ierr = ISLocalToGlobalMappingGetInfo(mapping,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
1536*071fcb05SBarry Smith     for (i=0;i<n;i++) { mapping->info_nodec[i] = 1;}
1537*071fcb05SBarry Smith     m = n;
1538*071fcb05SBarry Smith     mapping->info_nodec[n] = 0;
15391bd0b88eSStefano Zampini     for (i=1;i<n_neigh;i++) {
15401bd0b88eSStefano Zampini       PetscInt j;
15411bd0b88eSStefano Zampini 
15421bd0b88eSStefano Zampini       m += n_shared[i];
15431bd0b88eSStefano Zampini       for (j=0;j<n_shared[i];j++) mapping->info_nodec[shared[i][j]] += 1;
15441bd0b88eSStefano Zampini     }
15451bd0b88eSStefano Zampini     if (n) { ierr = PetscMalloc1(m,&mapping->info_nodei[0]);CHKERRQ(ierr); }
15461bd0b88eSStefano Zampini     for (i=1;i<n;i++) mapping->info_nodei[i] = mapping->info_nodei[i-1] + mapping->info_nodec[i-1];
1547580bdb30SBarry Smith     ierr = PetscArrayzero(mapping->info_nodec,n);CHKERRQ(ierr);
15481bd0b88eSStefano Zampini     for (i=0;i<n;i++) { mapping->info_nodec[i] = 1; mapping->info_nodei[i][0] = neigh[0]; }
15491bd0b88eSStefano Zampini     for (i=1;i<n_neigh;i++) {
15501bd0b88eSStefano Zampini       PetscInt j;
15511bd0b88eSStefano Zampini 
15521bd0b88eSStefano Zampini       for (j=0;j<n_shared[i];j++) {
15531bd0b88eSStefano Zampini         PetscInt k = shared[i][j];
15541bd0b88eSStefano Zampini 
15551bd0b88eSStefano Zampini         mapping->info_nodei[k][mapping->info_nodec[k]] = neigh[i];
15561bd0b88eSStefano Zampini         mapping->info_nodec[k] += 1;
15571bd0b88eSStefano Zampini       }
15581bd0b88eSStefano Zampini     }
15591bd0b88eSStefano Zampini     for (i=0;i<n;i++) { ierr = PetscSortRemoveDupsInt(&mapping->info_nodec[i],mapping->info_nodei[i]);CHKERRQ(ierr); }
15601bd0b88eSStefano Zampini     ierr = ISLocalToGlobalMappingRestoreInfo(mapping,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
15611bd0b88eSStefano Zampini   }
15621bd0b88eSStefano Zampini   if (nnodes)  *nnodes  = n;
15631bd0b88eSStefano Zampini   if (count)   *count   = mapping->info_nodec;
15641bd0b88eSStefano Zampini   if (indices) *indices = mapping->info_nodei;
15651bd0b88eSStefano Zampini   PetscFunctionReturn(0);
15661bd0b88eSStefano Zampini }
15671bd0b88eSStefano Zampini 
15681bd0b88eSStefano Zampini /*@C
15691bd0b88eSStefano Zampini     ISLocalToGlobalMappingRestoreNodeInfo - Frees the memory allocated by ISLocalToGlobalMappingGetNodeInfo()
15701bd0b88eSStefano Zampini 
15711bd0b88eSStefano Zampini     Collective on ISLocalToGlobalMapping
15721bd0b88eSStefano Zampini 
15731bd0b88eSStefano Zampini     Input Parameters:
15741bd0b88eSStefano Zampini .   mapping - the mapping from local to global indexing
15751bd0b88eSStefano Zampini 
15761bd0b88eSStefano Zampini     Output Parameter:
15771bd0b88eSStefano Zampini +   nnodes - number of local nodes
15781bd0b88eSStefano Zampini .   count - number of neighboring processors per node
15791bd0b88eSStefano Zampini -   indices - indices of processes sharing the node (sorted)
15801bd0b88eSStefano Zampini 
15811bd0b88eSStefano Zampini     Level: advanced
15821bd0b88eSStefano Zampini 
15831bd0b88eSStefano Zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
15841bd0b88eSStefano Zampini           ISLocalToGlobalMappingGetInfo()
15851bd0b88eSStefano Zampini @*/
15861bd0b88eSStefano Zampini PetscErrorCode  ISLocalToGlobalMappingRestoreNodeInfo(ISLocalToGlobalMapping mapping,PetscInt *nnodes,PetscInt *count[],PetscInt **indices[])
15871bd0b88eSStefano Zampini {
15881bd0b88eSStefano Zampini   PetscFunctionBegin;
15891bd0b88eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
15901bd0b88eSStefano Zampini   if (nnodes)  *nnodes  = 0;
15911bd0b88eSStefano Zampini   if (count)   *count   = NULL;
15921bd0b88eSStefano Zampini   if (indices) *indices = NULL;
15931bd0b88eSStefano Zampini   PetscFunctionReturn(0);
15941bd0b88eSStefano Zampini }
15951bd0b88eSStefano Zampini 
15961bd0b88eSStefano Zampini /*@C
1597107e9a97SBarry Smith    ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped
159886994e45SJed Brown 
159986994e45SJed Brown    Not Collective
160086994e45SJed Brown 
160186994e45SJed Brown    Input Arguments:
160286994e45SJed Brown . ltog - local to global mapping
160386994e45SJed Brown 
160486994e45SJed Brown    Output Arguments:
1605565245c5SBarry Smith . array - array of indices, the length of this array may be obtained with ISLocalToGlobalMappingGetSize()
160686994e45SJed Brown 
160786994e45SJed Brown    Level: advanced
160886994e45SJed Brown 
160995452b02SPatrick Sanan    Notes:
161095452b02SPatrick Sanan     ISLocalToGlobalMappingGetSize() returns the length the this array
1611107e9a97SBarry Smith 
1612107e9a97SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreIndices(), ISLocalToGlobalMappingGetBlockIndices(), ISLocalToGlobalMappingRestoreBlockIndices()
161386994e45SJed Brown @*/
16147087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
161586994e45SJed Brown {
161686994e45SJed Brown   PetscFunctionBegin;
161786994e45SJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
161886994e45SJed Brown   PetscValidPointer(array,2);
161945b6f7e9SBarry Smith   if (ltog->bs == 1) {
162086994e45SJed Brown     *array = ltog->indices;
162145b6f7e9SBarry Smith   } else {
162245b6f7e9SBarry Smith     PetscInt       *jj,k,i,j,n = ltog->n, bs = ltog->bs;
162345b6f7e9SBarry Smith     const PetscInt *ii;
162445b6f7e9SBarry Smith     PetscErrorCode ierr;
162545b6f7e9SBarry Smith 
162645b6f7e9SBarry Smith     ierr = PetscMalloc1(bs*n,&jj);CHKERRQ(ierr);
162745b6f7e9SBarry Smith     *array = jj;
162845b6f7e9SBarry Smith     k    = 0;
162945b6f7e9SBarry Smith     ii   = ltog->indices;
163045b6f7e9SBarry Smith     for (i=0; i<n; i++)
163145b6f7e9SBarry Smith       for (j=0; j<bs; j++)
163245b6f7e9SBarry Smith         jj[k++] = bs*ii[i] + j;
163345b6f7e9SBarry Smith   }
163486994e45SJed Brown   PetscFunctionReturn(0);
163586994e45SJed Brown }
163686994e45SJed Brown 
163786994e45SJed Brown /*@C
1638193a2b41SJulian Andrej    ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with ISLocalToGlobalMappingGetIndices()
163986994e45SJed Brown 
164086994e45SJed Brown    Not Collective
164186994e45SJed Brown 
164286994e45SJed Brown    Input Arguments:
164386994e45SJed Brown + ltog - local to global mapping
164486994e45SJed Brown - array - array of indices
164586994e45SJed Brown 
164686994e45SJed Brown    Level: advanced
164786994e45SJed Brown 
164886994e45SJed Brown .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices()
164986994e45SJed Brown @*/
16507087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
165186994e45SJed Brown {
165286994e45SJed Brown   PetscFunctionBegin;
165386994e45SJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
165486994e45SJed Brown   PetscValidPointer(array,2);
165545b6f7e9SBarry Smith   if (ltog->bs == 1 && *array != ltog->indices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
165645b6f7e9SBarry Smith 
165745b6f7e9SBarry Smith   if (ltog->bs > 1) {
165845b6f7e9SBarry Smith     PetscErrorCode ierr;
165945b6f7e9SBarry Smith     ierr = PetscFree(*(void**)array);CHKERRQ(ierr);
166045b6f7e9SBarry Smith   }
166145b6f7e9SBarry Smith   PetscFunctionReturn(0);
166245b6f7e9SBarry Smith }
166345b6f7e9SBarry Smith 
166445b6f7e9SBarry Smith /*@C
166545b6f7e9SBarry Smith    ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block
166645b6f7e9SBarry Smith 
166745b6f7e9SBarry Smith    Not Collective
166845b6f7e9SBarry Smith 
166945b6f7e9SBarry Smith    Input Arguments:
167045b6f7e9SBarry Smith . ltog - local to global mapping
167145b6f7e9SBarry Smith 
167245b6f7e9SBarry Smith    Output Arguments:
167345b6f7e9SBarry Smith . array - array of indices
167445b6f7e9SBarry Smith 
167545b6f7e9SBarry Smith    Level: advanced
167645b6f7e9SBarry Smith 
167745b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreBlockIndices()
167845b6f7e9SBarry Smith @*/
167945b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
168045b6f7e9SBarry Smith {
168145b6f7e9SBarry Smith   PetscFunctionBegin;
168245b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
168345b6f7e9SBarry Smith   PetscValidPointer(array,2);
168445b6f7e9SBarry Smith   *array = ltog->indices;
168545b6f7e9SBarry Smith   PetscFunctionReturn(0);
168645b6f7e9SBarry Smith }
168745b6f7e9SBarry Smith 
168845b6f7e9SBarry Smith /*@C
168945b6f7e9SBarry Smith    ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with ISLocalToGlobalMappingGetBlockIndices()
169045b6f7e9SBarry Smith 
169145b6f7e9SBarry Smith    Not Collective
169245b6f7e9SBarry Smith 
169345b6f7e9SBarry Smith    Input Arguments:
169445b6f7e9SBarry Smith + ltog - local to global mapping
169545b6f7e9SBarry Smith - array - array of indices
169645b6f7e9SBarry Smith 
169745b6f7e9SBarry Smith    Level: advanced
169845b6f7e9SBarry Smith 
169945b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices()
170045b6f7e9SBarry Smith @*/
170145b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
170245b6f7e9SBarry Smith {
170345b6f7e9SBarry Smith   PetscFunctionBegin;
170445b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
170545b6f7e9SBarry Smith   PetscValidPointer(array,2);
170686994e45SJed Brown   if (*array != ltog->indices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
17070298fd71SBarry Smith   *array = NULL;
170886994e45SJed Brown   PetscFunctionReturn(0);
170986994e45SJed Brown }
1710f7efa3c7SJed Brown 
1711f7efa3c7SJed Brown /*@C
1712f7efa3c7SJed Brown    ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings
1713f7efa3c7SJed Brown 
1714f7efa3c7SJed Brown    Not Collective
1715f7efa3c7SJed Brown 
1716f7efa3c7SJed Brown    Input Arguments:
1717f7efa3c7SJed Brown + comm - communicator for the new mapping, must contain the communicator of every mapping to concatenate
1718f7efa3c7SJed Brown . n - number of mappings to concatenate
1719f7efa3c7SJed Brown - ltogs - local to global mappings
1720f7efa3c7SJed Brown 
1721f7efa3c7SJed Brown    Output Arguments:
1722f7efa3c7SJed Brown . ltogcat - new mapping
1723f7efa3c7SJed Brown 
17249d90f715SBarry Smith    Note: this currently always returns a mapping with block size of 1
17259d90f715SBarry Smith 
17269d90f715SBarry Smith    Developer Note: If all the input mapping have the same block size we could easily handle that as a special case
17279d90f715SBarry Smith 
1728f7efa3c7SJed Brown    Level: advanced
1729f7efa3c7SJed Brown 
1730f7efa3c7SJed Brown .seealso: ISLocalToGlobalMappingCreate()
1731f7efa3c7SJed Brown @*/
1732f7efa3c7SJed Brown PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm,PetscInt n,const ISLocalToGlobalMapping ltogs[],ISLocalToGlobalMapping *ltogcat)
1733f7efa3c7SJed Brown {
1734f7efa3c7SJed Brown   PetscInt       i,cnt,m,*idx;
1735f7efa3c7SJed Brown   PetscErrorCode ierr;
1736f7efa3c7SJed Brown 
1737f7efa3c7SJed Brown   PetscFunctionBegin;
1738f7efa3c7SJed Brown   if (n < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have a non-negative number of mappings, given %D",n);
1739f7efa3c7SJed Brown   if (n > 0) PetscValidPointer(ltogs,3);
1740f7efa3c7SJed Brown   for (i=0; i<n; i++) PetscValidHeaderSpecific(ltogs[i],IS_LTOGM_CLASSID,3);
1741f7efa3c7SJed Brown   PetscValidPointer(ltogcat,4);
1742f7efa3c7SJed Brown   for (cnt=0,i=0; i<n; i++) {
1743f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
1744f7efa3c7SJed Brown     cnt += m;
1745f7efa3c7SJed Brown   }
1746785e854fSJed Brown   ierr = PetscMalloc1(cnt,&idx);CHKERRQ(ierr);
1747f7efa3c7SJed Brown   for (cnt=0,i=0; i<n; i++) {
1748f7efa3c7SJed Brown     const PetscInt *subidx;
1749f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
1750f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetIndices(ltogs[i],&subidx);CHKERRQ(ierr);
1751580bdb30SBarry Smith     ierr = PetscArraycpy(&idx[cnt],subidx,m);CHKERRQ(ierr);
1752f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingRestoreIndices(ltogs[i],&subidx);CHKERRQ(ierr);
1753f7efa3c7SJed Brown     cnt += m;
1754f7efa3c7SJed Brown   }
1755f0413b6fSBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,1,cnt,idx,PETSC_OWN_POINTER,ltogcat);CHKERRQ(ierr);
1756f7efa3c7SJed Brown   PetscFunctionReturn(0);
1757f7efa3c7SJed Brown }
175804a59952SBarry Smith 
1759413f72f0SBarry Smith /*MC
1760413f72f0SBarry Smith       ISLOCALTOGLOBALMAPPINGBASIC - basic implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is
1761413f72f0SBarry Smith                                     used this is good for only small and moderate size problems.
1762413f72f0SBarry Smith 
1763413f72f0SBarry Smith    Options Database Keys:
1764413f72f0SBarry Smith +   -islocaltoglobalmapping_type basic - select this method
1765413f72f0SBarry Smith 
1766413f72f0SBarry Smith    Level: beginner
1767413f72f0SBarry Smith 
1768413f72f0SBarry Smith .seealso:  ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType(), ISLOCALTOGLOBALMAPPINGHASH
1769413f72f0SBarry Smith M*/
1770413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Basic(ISLocalToGlobalMapping ltog)
1771413f72f0SBarry Smith {
1772413f72f0SBarry Smith   PetscFunctionBegin;
1773413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Basic;
1774413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Basic;
1775413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Basic;
1776413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Basic;
1777413f72f0SBarry Smith   PetscFunctionReturn(0);
1778413f72f0SBarry Smith }
1779413f72f0SBarry Smith 
1780413f72f0SBarry Smith /*MC
1781413f72f0SBarry Smith       ISLOCALTOGLOBALMAPPINGHASH - hash implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is
1782ed56e8eaSBarry Smith                                     used this is good for large memory problems.
1783413f72f0SBarry Smith 
1784413f72f0SBarry Smith    Options Database Keys:
1785413f72f0SBarry Smith +   -islocaltoglobalmapping_type hash - select this method
1786413f72f0SBarry Smith 
1787ed56e8eaSBarry Smith 
178895452b02SPatrick Sanan    Notes:
178995452b02SPatrick Sanan     This is selected automatically for large problems if the user does not set the type.
1790ed56e8eaSBarry Smith 
1791413f72f0SBarry Smith    Level: beginner
1792413f72f0SBarry Smith 
1793413f72f0SBarry Smith .seealso:  ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType(), ISLOCALTOGLOBALMAPPINGHASH
1794413f72f0SBarry Smith M*/
1795413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Hash(ISLocalToGlobalMapping ltog)
1796413f72f0SBarry Smith {
1797413f72f0SBarry Smith   PetscFunctionBegin;
1798413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Hash;
1799413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Hash;
1800413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Hash;
1801413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Hash;
1802413f72f0SBarry Smith   PetscFunctionReturn(0);
1803413f72f0SBarry Smith }
1804413f72f0SBarry Smith 
1805413f72f0SBarry Smith 
1806413f72f0SBarry Smith /*@C
1807413f72f0SBarry Smith     ISLocalToGlobalMappingRegister -  Adds a method for applying a global to local mapping with an ISLocalToGlobalMapping
1808413f72f0SBarry Smith 
1809413f72f0SBarry Smith    Not Collective
1810413f72f0SBarry Smith 
1811413f72f0SBarry Smith    Input Parameters:
1812413f72f0SBarry Smith +  sname - name of a new method
1813413f72f0SBarry Smith -  routine_create - routine to create method context
1814413f72f0SBarry Smith 
1815413f72f0SBarry Smith    Notes:
1816ed56e8eaSBarry Smith    ISLocalToGlobalMappingRegister() may be called multiple times to add several user-defined mappings.
1817413f72f0SBarry Smith 
1818413f72f0SBarry Smith    Sample usage:
1819413f72f0SBarry Smith .vb
1820413f72f0SBarry Smith    ISLocalToGlobalMappingRegister("my_mapper",MyCreate);
1821413f72f0SBarry Smith .ve
1822413f72f0SBarry Smith 
1823ed56e8eaSBarry Smith    Then, your mapping can be chosen with the procedural interface via
1824413f72f0SBarry Smith $     ISLocalToGlobalMappingSetType(ltog,"my_mapper")
1825413f72f0SBarry Smith    or at runtime via the option
1826ed56e8eaSBarry Smith $     -islocaltoglobalmapping_type my_mapper
1827413f72f0SBarry Smith 
1828413f72f0SBarry Smith    Level: advanced
1829413f72f0SBarry Smith 
1830413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingRegisterAll(), ISLocalToGlobalMappingRegisterDestroy(), ISLOCALTOGLOBALMAPPINGBASIC, ISLOCALTOGLOBALMAPPINGHASH
1831413f72f0SBarry Smith 
1832413f72f0SBarry Smith @*/
1833413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRegister(const char sname[],PetscErrorCode (*function)(ISLocalToGlobalMapping))
1834413f72f0SBarry Smith {
1835413f72f0SBarry Smith   PetscErrorCode ierr;
1836413f72f0SBarry Smith 
1837413f72f0SBarry Smith   PetscFunctionBegin;
18381d36bdfdSBarry Smith   ierr = ISInitializePackage();CHKERRQ(ierr);
1839413f72f0SBarry Smith   ierr = PetscFunctionListAdd(&ISLocalToGlobalMappingList,sname,function);CHKERRQ(ierr);
1840413f72f0SBarry Smith   PetscFunctionReturn(0);
1841413f72f0SBarry Smith }
1842413f72f0SBarry Smith 
1843413f72f0SBarry Smith /*@C
1844ed56e8eaSBarry Smith    ISLocalToGlobalMappingSetType - Builds ISLocalToGlobalMapping for a particular global to local mapping approach.
1845413f72f0SBarry Smith 
1846413f72f0SBarry Smith    Logically Collective on ISLocalToGlobalMapping
1847413f72f0SBarry Smith 
1848413f72f0SBarry Smith    Input Parameters:
1849413f72f0SBarry Smith +  ltog - the ISLocalToGlobalMapping object
1850413f72f0SBarry Smith -  type - a known method
1851413f72f0SBarry Smith 
1852413f72f0SBarry Smith    Options Database Key:
1853ed56e8eaSBarry Smith .  -islocaltoglobalmapping_type  <method> - Sets the method; use -help for a list
1854413f72f0SBarry Smith     of available methods (for instance, basic or hash)
1855413f72f0SBarry Smith 
1856413f72f0SBarry Smith    Notes:
1857413f72f0SBarry Smith    See "petsc/include/petscis.h" for available methods
1858413f72f0SBarry Smith 
1859413f72f0SBarry Smith   Normally, it is best to use the ISLocalToGlobalMappingSetFromOptions() command and
1860413f72f0SBarry Smith   then set the ISLocalToGlobalMapping type from the options database rather than by using
1861413f72f0SBarry Smith   this routine.
1862413f72f0SBarry Smith 
1863413f72f0SBarry Smith   Level: intermediate
1864413f72f0SBarry Smith 
1865413f72f0SBarry Smith   Developer Note: ISLocalToGlobalMappingRegister() is used to add new types to ISLocalToGlobalMappingList from which they
1866413f72f0SBarry Smith   are accessed by ISLocalToGlobalMappingSetType().
1867413f72f0SBarry Smith 
1868413f72f0SBarry Smith .seealso: PCSetType(), ISLocalToGlobalMappingType, ISLocalToGlobalMappingRegister(), ISLocalToGlobalMappingCreate()
1869413f72f0SBarry Smith 
1870413f72f0SBarry Smith @*/
1871413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType type)
1872413f72f0SBarry Smith {
1873413f72f0SBarry Smith   PetscErrorCode ierr,(*r)(ISLocalToGlobalMapping);
1874413f72f0SBarry Smith   PetscBool      match;
1875413f72f0SBarry Smith 
1876413f72f0SBarry Smith   PetscFunctionBegin;
1877413f72f0SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
1878413f72f0SBarry Smith   PetscValidCharPointer(type,2);
1879413f72f0SBarry Smith 
1880413f72f0SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)ltog,type,&match);CHKERRQ(ierr);
1881413f72f0SBarry Smith   if (match) PetscFunctionReturn(0);
1882413f72f0SBarry Smith 
1883413f72f0SBarry Smith   ierr =  PetscFunctionListFind(ISLocalToGlobalMappingList,type,&r);CHKERRQ(ierr);
1884413f72f0SBarry Smith   if (!r) SETERRQ1(PetscObjectComm((PetscObject)ltog),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested ISLocalToGlobalMapping type %s",type);
1885413f72f0SBarry Smith   /* Destroy the previous private LTOG context */
1886413f72f0SBarry Smith   if (ltog->ops->destroy) {
1887413f72f0SBarry Smith     ierr              = (*ltog->ops->destroy)(ltog);CHKERRQ(ierr);
1888413f72f0SBarry Smith     ltog->ops->destroy = NULL;
1889413f72f0SBarry Smith   }
1890413f72f0SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)ltog,type);CHKERRQ(ierr);
1891413f72f0SBarry Smith   ierr = (*r)(ltog);CHKERRQ(ierr);
1892413f72f0SBarry Smith   PetscFunctionReturn(0);
1893413f72f0SBarry Smith }
1894413f72f0SBarry Smith 
1895413f72f0SBarry Smith PetscBool ISLocalToGlobalMappingRegisterAllCalled = PETSC_FALSE;
1896413f72f0SBarry Smith 
1897413f72f0SBarry Smith /*@C
1898413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAll - Registers all of the local to global mapping components in the IS package.
1899413f72f0SBarry Smith 
1900413f72f0SBarry Smith   Not Collective
1901413f72f0SBarry Smith 
1902413f72f0SBarry Smith   Level: advanced
1903413f72f0SBarry Smith 
1904413f72f0SBarry Smith .seealso:  ISRegister(),  ISLocalToGlobalRegister()
1905413f72f0SBarry Smith @*/
1906413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRegisterAll(void)
1907413f72f0SBarry Smith {
1908413f72f0SBarry Smith   PetscErrorCode ierr;
1909413f72f0SBarry Smith 
1910413f72f0SBarry Smith   PetscFunctionBegin;
1911413f72f0SBarry Smith   if (ISLocalToGlobalMappingRegisterAllCalled) PetscFunctionReturn(0);
1912413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAllCalled = PETSC_TRUE;
1913413f72f0SBarry Smith 
1914413f72f0SBarry Smith   ierr = ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGBASIC, ISLocalToGlobalMappingCreate_Basic);CHKERRQ(ierr);
1915413f72f0SBarry Smith   ierr = ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGHASH, ISLocalToGlobalMappingCreate_Hash);CHKERRQ(ierr);
1916413f72f0SBarry Smith   PetscFunctionReturn(0);
1917413f72f0SBarry Smith }
191804a59952SBarry Smith 
1919