xref: /petsc/src/vec/is/utils/isltog.c (revision 820f2d463ad7e8160c4f507e09c88f621629a79d)
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 
186528b96dSMatthew G. Knepley /*@C
196528b96dSMatthew G. Knepley   ISGetPointRange - Returns a description of the points in an IS suitable for traversal
20413f72f0SBarry Smith 
216528b96dSMatthew G. Knepley   Not collective
226528b96dSMatthew G. Knepley 
236528b96dSMatthew G. Knepley   Input Parameter:
246528b96dSMatthew G. Knepley . pointIS - The IS object
256528b96dSMatthew G. Knepley 
266528b96dSMatthew G. Knepley   Output Parameters:
276528b96dSMatthew G. Knepley + pStart - The first index, see notes
286528b96dSMatthew G. Knepley . pEnd   - One past the last index, see notes
296528b96dSMatthew G. Knepley - points - The indices, see notes
306528b96dSMatthew G. Knepley 
316528b96dSMatthew G. Knepley   Notes:
326528b96dSMatthew G. Knepley   If the IS contains contiguous indices in an ISSTRIDE, then the indices are contained in [pStart, pEnd) and points = NULL. Otherwise, pStart = 0, pEnd = numIndices, and points is an array of the indices. This supports the following pattern
336528b96dSMatthew G. Knepley $ ISGetPointRange(is, &pStart, &pEnd, &points);
346528b96dSMatthew G. Knepley $ for (p = pStart; p < pEnd; ++p) {
356528b96dSMatthew G. Knepley $   const PetscInt point = points ? points[p] : p;
366528b96dSMatthew G. Knepley $ }
376528b96dSMatthew G. Knepley $ ISRestorePointRange(is, &pstart, &pEnd, &points);
386528b96dSMatthew G. Knepley 
396528b96dSMatthew G. Knepley   Level: intermediate
406528b96dSMatthew G. Knepley 
416528b96dSMatthew G. Knepley .seealso: ISRestorePointRange(), ISGetPointSubrange(), ISGetIndices(), ISCreateStride()
426528b96dSMatthew G. Knepley @*/
439305a4c7SMatthew G. Knepley PetscErrorCode ISGetPointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points)
449305a4c7SMatthew G. Knepley {
459305a4c7SMatthew G. Knepley   PetscInt       numCells, step = 1;
469305a4c7SMatthew G. Knepley   PetscBool      isStride;
479305a4c7SMatthew G. Knepley   PetscErrorCode ierr;
489305a4c7SMatthew G. Knepley 
499305a4c7SMatthew G. Knepley   PetscFunctionBeginHot;
509305a4c7SMatthew G. Knepley   *pStart = 0;
519305a4c7SMatthew G. Knepley   *points = NULL;
529305a4c7SMatthew G. Knepley   ierr = ISGetLocalSize(pointIS, &numCells);CHKERRQ(ierr);
539305a4c7SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) pointIS, ISSTRIDE, &isStride);CHKERRQ(ierr);
549305a4c7SMatthew G. Knepley   if (isStride) {ierr = ISStrideGetInfo(pointIS, pStart, &step);CHKERRQ(ierr);}
559305a4c7SMatthew G. Knepley   *pEnd   = *pStart + numCells;
569305a4c7SMatthew G. Knepley   if (!isStride || step != 1) {ierr = ISGetIndices(pointIS, points);CHKERRQ(ierr);}
579305a4c7SMatthew G. Knepley   PetscFunctionReturn(0);
589305a4c7SMatthew G. Knepley }
599305a4c7SMatthew G. Knepley 
606528b96dSMatthew G. Knepley /*@C
616528b96dSMatthew G. Knepley   ISRestorePointRange - Destroys the traversal description
626528b96dSMatthew G. Knepley 
636528b96dSMatthew G. Knepley   Not collective
646528b96dSMatthew G. Knepley 
656528b96dSMatthew G. Knepley   Input Parameters:
666528b96dSMatthew G. Knepley + pointIS - The IS object
676528b96dSMatthew G. Knepley . pStart  - The first index, from ISGetPointRange()
686528b96dSMatthew G. Knepley . pEnd    - One past the last index, from ISGetPointRange()
696528b96dSMatthew G. Knepley - points  - The indices, from ISGetPointRange()
706528b96dSMatthew G. Knepley 
716528b96dSMatthew G. Knepley   Notes:
726528b96dSMatthew G. Knepley   If the IS contains contiguous indices in an ISSTRIDE, then the indices are contained in [pStart, pEnd) and points = NULL. Otherwise, pStart = 0, pEnd = numIndices, and points is an array of the indices. This supports the following pattern
736528b96dSMatthew G. Knepley $ ISGetPointRange(is, &pStart, &pEnd, &points);
746528b96dSMatthew G. Knepley $ for (p = pStart; p < pEnd; ++p) {
756528b96dSMatthew G. Knepley $   const PetscInt point = points ? points[p] : p;
766528b96dSMatthew G. Knepley $ }
776528b96dSMatthew G. Knepley $ ISRestorePointRange(is, &pstart, &pEnd, &points);
786528b96dSMatthew G. Knepley 
796528b96dSMatthew G. Knepley   Level: intermediate
806528b96dSMatthew G. Knepley 
816528b96dSMatthew G. Knepley .seealso: ISGetPointRange(), ISGetPointSubrange(), ISGetIndices(), ISCreateStride()
826528b96dSMatthew G. Knepley @*/
839305a4c7SMatthew G. Knepley PetscErrorCode ISRestorePointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points)
849305a4c7SMatthew G. Knepley {
859305a4c7SMatthew G. Knepley   PetscInt       step = 1;
869305a4c7SMatthew G. Knepley   PetscBool      isStride;
879305a4c7SMatthew G. Knepley   PetscErrorCode ierr;
889305a4c7SMatthew G. Knepley 
899305a4c7SMatthew G. Knepley   PetscFunctionBeginHot;
909305a4c7SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) pointIS, ISSTRIDE, &isStride);CHKERRQ(ierr);
919305a4c7SMatthew G. Knepley   if (isStride) {ierr = ISStrideGetInfo(pointIS, pStart, &step);CHKERRQ(ierr);}
929305a4c7SMatthew G. Knepley   if (!isStride || step != 1) {ierr = ISGetIndices(pointIS, points);CHKERRQ(ierr);}
939305a4c7SMatthew G. Knepley   PetscFunctionReturn(0);
949305a4c7SMatthew G. Knepley }
959305a4c7SMatthew G. Knepley 
966528b96dSMatthew G. Knepley /*@C
976528b96dSMatthew G. Knepley   ISGetPointSubrange - Configures the input IS to be a subrange for the traversal information given
986528b96dSMatthew G. Knepley 
996528b96dSMatthew G. Knepley   Not collective
1006528b96dSMatthew G. Knepley 
1016528b96dSMatthew G. Knepley   Input Parameters:
1026528b96dSMatthew G. Knepley + subpointIS - The IS object to be configured
1036528b96dSMatthew G. Knepley . pStar   t  - The first index of the subrange
1046528b96dSMatthew G. Knepley . pEnd       - One past the last index for the subrange
1056528b96dSMatthew G. Knepley - points     - The indices for the entire range, from ISGetPointRange()
1066528b96dSMatthew G. Knepley 
1076528b96dSMatthew G. Knepley   Output Parameters:
1086528b96dSMatthew G. Knepley . subpointIS - The IS object now configured to be a subrange
1096528b96dSMatthew G. Knepley 
1106528b96dSMatthew G. Knepley   Notes:
1116528b96dSMatthew G. Knepley   The input IS will now respond properly to calls to ISGetPointRange() and return the subrange.
1126528b96dSMatthew G. Knepley 
1136528b96dSMatthew G. Knepley   Level: intermediate
1146528b96dSMatthew G. Knepley 
1156528b96dSMatthew G. Knepley .seealso: ISGetPointRange(), ISRestorePointRange(), ISGetIndices(), ISCreateStride()
1166528b96dSMatthew G. Knepley @*/
1179305a4c7SMatthew G. Knepley PetscErrorCode ISGetPointSubrange(IS subpointIS, PetscInt pStart, PetscInt pEnd, const PetscInt *points)
1189305a4c7SMatthew G. Knepley {
1199305a4c7SMatthew G. Knepley   PetscErrorCode ierr;
1209305a4c7SMatthew G. Knepley 
1219305a4c7SMatthew G. Knepley   PetscFunctionBeginHot;
1229305a4c7SMatthew G. Knepley   if (points) {
1239305a4c7SMatthew G. Knepley     ierr = ISSetType(subpointIS, ISGENERAL);CHKERRQ(ierr);
1249305a4c7SMatthew G. Knepley     ierr = ISGeneralSetIndices(subpointIS, pEnd-pStart, &points[pStart], PETSC_USE_POINTER);CHKERRQ(ierr);
1259305a4c7SMatthew G. Knepley   } else {
1269305a4c7SMatthew G. Knepley     ierr = ISSetType(subpointIS, ISSTRIDE);CHKERRQ(ierr);
1279305a4c7SMatthew G. Knepley     ierr = ISStrideSetStride(subpointIS, pEnd-pStart, pStart, 1);CHKERRQ(ierr);
1289305a4c7SMatthew G. Knepley   }
1299305a4c7SMatthew G. Knepley   PetscFunctionReturn(0);
1309305a4c7SMatthew G. Knepley }
1319305a4c7SMatthew G. Knepley 
132413f72f0SBarry Smith /* -----------------------------------------------------------------------------------------*/
133413f72f0SBarry Smith 
134413f72f0SBarry Smith /*
135413f72f0SBarry Smith     Creates the global mapping information in the ISLocalToGlobalMapping structure
136413f72f0SBarry Smith 
137413f72f0SBarry Smith     If the user has not selected how to handle the global to local mapping then use HASH for "large" problems
138413f72f0SBarry Smith */
139413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp(ISLocalToGlobalMapping mapping)
140413f72f0SBarry Smith {
141413f72f0SBarry Smith   PetscInt       i,*idx = mapping->indices,n = mapping->n,end,start;
142413f72f0SBarry Smith   PetscErrorCode ierr;
143413f72f0SBarry Smith 
144413f72f0SBarry Smith   PetscFunctionBegin;
145413f72f0SBarry Smith   if (mapping->data) PetscFunctionReturn(0);
146413f72f0SBarry Smith   end   = 0;
147413f72f0SBarry Smith   start = PETSC_MAX_INT;
148413f72f0SBarry Smith 
149413f72f0SBarry Smith   for (i=0; i<n; i++) {
150413f72f0SBarry Smith     if (idx[i] < 0) continue;
151413f72f0SBarry Smith     if (idx[i] < start) start = idx[i];
152413f72f0SBarry Smith     if (idx[i] > end)   end   = idx[i];
153413f72f0SBarry Smith   }
154413f72f0SBarry Smith   if (start > end) {start = 0; end = -1;}
155413f72f0SBarry Smith   mapping->globalstart = start;
156413f72f0SBarry Smith   mapping->globalend   = end;
157413f72f0SBarry Smith   if (!((PetscObject)mapping)->type_name) {
158413f72f0SBarry Smith     if ((end - start) > PetscMax(4*n,1000000)) {
1597f79407eSBarry Smith       ierr = ISLocalToGlobalMappingSetType(mapping,ISLOCALTOGLOBALMAPPINGHASH);CHKERRQ(ierr);
160413f72f0SBarry Smith     } else {
1617f79407eSBarry Smith       ierr = ISLocalToGlobalMappingSetType(mapping,ISLOCALTOGLOBALMAPPINGBASIC);CHKERRQ(ierr);
162413f72f0SBarry Smith     }
163413f72f0SBarry Smith   }
164413f72f0SBarry Smith   ierr = (*mapping->ops->globaltolocalmappingsetup)(mapping);CHKERRQ(ierr);
165413f72f0SBarry Smith   PetscFunctionReturn(0);
166413f72f0SBarry Smith }
167413f72f0SBarry Smith 
168413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp_Basic(ISLocalToGlobalMapping mapping)
169413f72f0SBarry Smith {
170413f72f0SBarry Smith   PetscErrorCode              ierr;
171413f72f0SBarry Smith   PetscInt                    i,*idx = mapping->indices,n = mapping->n,end,start,*globals;
172413f72f0SBarry Smith   ISLocalToGlobalMapping_Basic *map;
173413f72f0SBarry Smith 
174413f72f0SBarry Smith   PetscFunctionBegin;
175413f72f0SBarry Smith   start            = mapping->globalstart;
176413f72f0SBarry Smith   end              = mapping->globalend;
177413f72f0SBarry Smith   ierr             = PetscNew(&map);CHKERRQ(ierr);
178413f72f0SBarry Smith   ierr             = PetscMalloc1(end-start+2,&globals);CHKERRQ(ierr);
179413f72f0SBarry Smith   map->globals     = globals;
180413f72f0SBarry Smith   for (i=0; i<end-start+1; i++) globals[i] = -1;
181413f72f0SBarry Smith   for (i=0; i<n; i++) {
182413f72f0SBarry Smith     if (idx[i] < 0) continue;
183413f72f0SBarry Smith     globals[idx[i] - start] = i;
184413f72f0SBarry Smith   }
185413f72f0SBarry Smith   mapping->data = (void*)map;
186413f72f0SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)mapping,(end-start+1)*sizeof(PetscInt));CHKERRQ(ierr);
187413f72f0SBarry Smith   PetscFunctionReturn(0);
188413f72f0SBarry Smith }
189413f72f0SBarry Smith 
190413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp_Hash(ISLocalToGlobalMapping mapping)
191413f72f0SBarry Smith {
192413f72f0SBarry Smith   PetscErrorCode              ierr;
193413f72f0SBarry Smith   PetscInt                    i,*idx = mapping->indices,n = mapping->n;
194413f72f0SBarry Smith   ISLocalToGlobalMapping_Hash *map;
195413f72f0SBarry Smith 
196413f72f0SBarry Smith   PetscFunctionBegin;
197413f72f0SBarry Smith   ierr = PetscNew(&map);CHKERRQ(ierr);
198e8f14785SLisandro Dalcin   ierr = PetscHMapICreate(&map->globalht);CHKERRQ(ierr);
199413f72f0SBarry Smith   for (i=0; i<n; i++) {
200413f72f0SBarry Smith     if (idx[i] < 0) continue;
201e8f14785SLisandro Dalcin     ierr = PetscHMapISet(map->globalht,idx[i],i);CHKERRQ(ierr);
202413f72f0SBarry Smith   }
203413f72f0SBarry Smith   mapping->data = (void*)map;
204413f72f0SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)mapping,2*n*sizeof(PetscInt));CHKERRQ(ierr);
205413f72f0SBarry Smith   PetscFunctionReturn(0);
206413f72f0SBarry Smith }
207413f72f0SBarry Smith 
208413f72f0SBarry Smith static PetscErrorCode ISLocalToGlobalMappingDestroy_Basic(ISLocalToGlobalMapping mapping)
209413f72f0SBarry Smith {
210413f72f0SBarry Smith   PetscErrorCode              ierr;
211413f72f0SBarry Smith   ISLocalToGlobalMapping_Basic *map  = (ISLocalToGlobalMapping_Basic *)mapping->data;
212413f72f0SBarry Smith 
213413f72f0SBarry Smith   PetscFunctionBegin;
214413f72f0SBarry Smith   if (!map) PetscFunctionReturn(0);
215413f72f0SBarry Smith   ierr = PetscFree(map->globals);CHKERRQ(ierr);
216413f72f0SBarry Smith   ierr = PetscFree(mapping->data);CHKERRQ(ierr);
217413f72f0SBarry Smith   PetscFunctionReturn(0);
218413f72f0SBarry Smith }
219413f72f0SBarry Smith 
220413f72f0SBarry Smith static PetscErrorCode ISLocalToGlobalMappingDestroy_Hash(ISLocalToGlobalMapping mapping)
221413f72f0SBarry Smith {
222413f72f0SBarry Smith   PetscErrorCode              ierr;
223413f72f0SBarry Smith   ISLocalToGlobalMapping_Hash *map  = (ISLocalToGlobalMapping_Hash*)mapping->data;
224413f72f0SBarry Smith 
225413f72f0SBarry Smith   PetscFunctionBegin;
226413f72f0SBarry Smith   if (!map) PetscFunctionReturn(0);
227e8f14785SLisandro Dalcin   ierr = PetscHMapIDestroy(&map->globalht);CHKERRQ(ierr);
228413f72f0SBarry Smith   ierr = PetscFree(mapping->data);CHKERRQ(ierr);
229413f72f0SBarry Smith   PetscFunctionReturn(0);
230413f72f0SBarry Smith }
231413f72f0SBarry Smith 
232413f72f0SBarry Smith #define GTOLTYPE _Basic
233413f72f0SBarry Smith #define GTOLNAME _Basic
234541bf97eSAdrian Croucher #define GTOLBS mapping->bs
235413f72f0SBarry Smith #define GTOL(g, local) do {                  \
236413f72f0SBarry Smith     local = map->globals[g/bs - start];      \
2370040bde1SJunchao Zhang     if (local >= 0) local = bs*local + (g % bs); \
238413f72f0SBarry Smith   } while (0)
239541bf97eSAdrian Croucher 
240413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
241413f72f0SBarry Smith 
242413f72f0SBarry Smith #define GTOLTYPE _Basic
243413f72f0SBarry Smith #define GTOLNAME Block_Basic
244541bf97eSAdrian Croucher #define GTOLBS 1
245413f72f0SBarry Smith #define GTOL(g, local) do {                  \
246413f72f0SBarry Smith     local = map->globals[g - start];         \
247413f72f0SBarry Smith   } while (0)
248413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
249413f72f0SBarry Smith 
250413f72f0SBarry Smith #define GTOLTYPE _Hash
251413f72f0SBarry Smith #define GTOLNAME _Hash
252541bf97eSAdrian Croucher #define GTOLBS mapping->bs
253413f72f0SBarry Smith #define GTOL(g, local) do {                         \
254e8f14785SLisandro Dalcin     (void)PetscHMapIGet(map->globalht,g/bs,&local); \
2550040bde1SJunchao Zhang     if (local >= 0) local = bs*local + (g % bs);   \
256413f72f0SBarry Smith    } while (0)
257413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
258413f72f0SBarry Smith 
259413f72f0SBarry Smith #define GTOLTYPE _Hash
260413f72f0SBarry Smith #define GTOLNAME Block_Hash
261541bf97eSAdrian Croucher #define GTOLBS 1
262413f72f0SBarry Smith #define GTOL(g, local) do {                         \
263e8f14785SLisandro Dalcin     (void)PetscHMapIGet(map->globalht,g,&local);    \
264413f72f0SBarry Smith   } while (0)
265413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
266413f72f0SBarry Smith 
2676658fb44Sstefano_zampini /*@
2686658fb44Sstefano_zampini     ISLocalToGlobalMappingDuplicate - Duplicates the local to global mapping object
2696658fb44Sstefano_zampini 
2706658fb44Sstefano_zampini     Not Collective
2716658fb44Sstefano_zampini 
2726658fb44Sstefano_zampini     Input Parameter:
2736658fb44Sstefano_zampini .   ltog - local to global mapping
2746658fb44Sstefano_zampini 
2756658fb44Sstefano_zampini     Output Parameter:
2766658fb44Sstefano_zampini .   nltog - the duplicated local to global mapping
2776658fb44Sstefano_zampini 
2786658fb44Sstefano_zampini     Level: advanced
2796658fb44Sstefano_zampini 
2806658fb44Sstefano_zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
2816658fb44Sstefano_zampini @*/
2826658fb44Sstefano_zampini PetscErrorCode  ISLocalToGlobalMappingDuplicate(ISLocalToGlobalMapping ltog,ISLocalToGlobalMapping* nltog)
2836658fb44Sstefano_zampini {
2846658fb44Sstefano_zampini   PetscErrorCode ierr;
2856658fb44Sstefano_zampini 
2866658fb44Sstefano_zampini   PetscFunctionBegin;
2876658fb44Sstefano_zampini   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
2886658fb44Sstefano_zampini   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)ltog),ltog->bs,ltog->n,ltog->indices,PETSC_COPY_VALUES,nltog);CHKERRQ(ierr);
2896658fb44Sstefano_zampini   PetscFunctionReturn(0);
2906658fb44Sstefano_zampini }
2916658fb44Sstefano_zampini 
292565245c5SBarry Smith /*@
293107e9a97SBarry Smith     ISLocalToGlobalMappingGetSize - Gets the local size of a local to global mapping
2943b9aefa3SBarry Smith 
2953b9aefa3SBarry Smith     Not Collective
2963b9aefa3SBarry Smith 
2973b9aefa3SBarry Smith     Input Parameter:
2983b9aefa3SBarry Smith .   ltog - local to global mapping
2993b9aefa3SBarry Smith 
3003b9aefa3SBarry Smith     Output Parameter:
301107e9a97SBarry Smith .   n - the number of entries in the local mapping, ISLocalToGlobalMappingGetIndices() returns an array of this length
3023b9aefa3SBarry Smith 
3033b9aefa3SBarry Smith     Level: advanced
3043b9aefa3SBarry Smith 
3053b9aefa3SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
3063b9aefa3SBarry Smith @*/
3077087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping mapping,PetscInt *n)
3083b9aefa3SBarry Smith {
3093b9aefa3SBarry Smith   PetscFunctionBegin;
3100700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
3114482741eSBarry Smith   PetscValidIntPointer(n,2);
312107e9a97SBarry Smith   *n = mapping->bs*mapping->n;
3133b9aefa3SBarry Smith   PetscFunctionReturn(0);
3143b9aefa3SBarry Smith }
3153b9aefa3SBarry Smith 
3165a5d4f66SBarry Smith /*@C
317fe2efc57SMark    ISLocalToGlobalMappingViewFromOptions - View from Options
318fe2efc57SMark 
319fe2efc57SMark    Collective on ISLocalToGlobalMapping
320fe2efc57SMark 
321fe2efc57SMark    Input Parameters:
322fe2efc57SMark +  A - the local to global mapping object
323736c3998SJose E. Roman .  obj - Optional object
324736c3998SJose E. Roman -  name - command line option
325fe2efc57SMark 
326fe2efc57SMark    Level: intermediate
327fe2efc57SMark .seealso:  ISLocalToGlobalMapping, ISLocalToGlobalMappingView, PetscObjectViewFromOptions(), ISLocalToGlobalMappingCreate()
328fe2efc57SMark @*/
329fe2efc57SMark PetscErrorCode  ISLocalToGlobalMappingViewFromOptions(ISLocalToGlobalMapping A,PetscObject obj,const char name[])
330fe2efc57SMark {
331fe2efc57SMark   PetscErrorCode ierr;
332fe2efc57SMark 
333fe2efc57SMark   PetscFunctionBegin;
334fe2efc57SMark   PetscValidHeaderSpecific(A,IS_LTOGM_CLASSID,1);
335fe2efc57SMark   ierr = PetscObjectViewFromOptions((PetscObject)A,obj,name);CHKERRQ(ierr);
336fe2efc57SMark   PetscFunctionReturn(0);
337fe2efc57SMark }
338fe2efc57SMark 
339fe2efc57SMark /*@C
3405a5d4f66SBarry Smith     ISLocalToGlobalMappingView - View a local to global mapping
3415a5d4f66SBarry Smith 
342b9cd556bSLois Curfman McInnes     Not Collective
343b9cd556bSLois Curfman McInnes 
3445a5d4f66SBarry Smith     Input Parameters:
3453b9aefa3SBarry Smith +   ltog - local to global mapping
3463b9aefa3SBarry Smith -   viewer - viewer
3475a5d4f66SBarry Smith 
348a997ad1aSLois Curfman McInnes     Level: advanced
349a997ad1aSLois Curfman McInnes 
3505a5d4f66SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate()
3515a5d4f66SBarry Smith @*/
3527087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingView(ISLocalToGlobalMapping mapping,PetscViewer viewer)
3535a5d4f66SBarry Smith {
35432dcc486SBarry Smith   PetscInt       i;
35532dcc486SBarry Smith   PetscMPIInt    rank;
356ace3abfcSBarry Smith   PetscBool      iascii;
3576849ba73SBarry Smith   PetscErrorCode ierr;
3585a5d4f66SBarry Smith 
3595a5d4f66SBarry Smith   PetscFunctionBegin;
3600700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
3613050cee2SBarry Smith   if (!viewer) {
362ce94432eSBarry Smith     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mapping),&viewer);CHKERRQ(ierr);
3633050cee2SBarry Smith   }
3640700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
3655a5d4f66SBarry Smith 
366ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mapping),&rank);CHKERRMPI(ierr);
367251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
36832077d6dSBarry Smith   if (iascii) {
36998c3331eSBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)mapping,viewer);CHKERRQ(ierr);
3701575c14dSBarry Smith     ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
3715a5d4f66SBarry Smith     for (i=0; i<mapping->n; i++) {
3727904a332SBarry Smith       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D %D\n",rank,i,mapping->indices[i]);CHKERRQ(ierr);
3736831982aSBarry Smith     }
374b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
3751575c14dSBarry Smith     ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
3761575c14dSBarry Smith   }
3775a5d4f66SBarry Smith   PetscFunctionReturn(0);
3785a5d4f66SBarry Smith }
3795a5d4f66SBarry Smith 
3801f428162SBarry Smith /*@
3812bdab257SBarry Smith     ISLocalToGlobalMappingCreateIS - Creates a mapping between a local (0 to n)
3822bdab257SBarry Smith     ordering and a global parallel ordering.
3832bdab257SBarry Smith 
3840f5bd95cSBarry Smith     Not collective
385b9cd556bSLois Curfman McInnes 
386a997ad1aSLois Curfman McInnes     Input Parameter:
3878c03b21aSDmitry Karpeev .   is - index set containing the global numbers for each local number
3882bdab257SBarry Smith 
389a997ad1aSLois Curfman McInnes     Output Parameter:
3902bdab257SBarry Smith .   mapping - new mapping data structure
3912bdab257SBarry Smith 
39295452b02SPatrick Sanan     Notes:
39395452b02SPatrick Sanan     the block size of the IS determines the block size of the mapping
394a997ad1aSLois Curfman McInnes     Level: advanced
395a997ad1aSLois Curfman McInnes 
3967e99dc12SLawrence Mitchell .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetFromOptions()
3972bdab257SBarry Smith @*/
3987087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingCreateIS(IS is,ISLocalToGlobalMapping *mapping)
3992bdab257SBarry Smith {
4006849ba73SBarry Smith   PetscErrorCode ierr;
4013bbf0e92SBarry Smith   PetscInt       n,bs;
4025d0c19d7SBarry Smith   const PetscInt *indices;
4032bdab257SBarry Smith   MPI_Comm       comm;
4043bbf0e92SBarry Smith   PetscBool      isblock;
4053a40ed3dSBarry Smith 
4063a40ed3dSBarry Smith   PetscFunctionBegin;
4070700a824SBarry Smith   PetscValidHeaderSpecific(is,IS_CLASSID,1);
4084482741eSBarry Smith   PetscValidPointer(mapping,2);
4092bdab257SBarry Smith 
4102bdab257SBarry Smith   ierr = PetscObjectGetComm((PetscObject)is,&comm);CHKERRQ(ierr);
4113b9aefa3SBarry Smith   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
4123bbf0e92SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)is,ISBLOCK,&isblock);CHKERRQ(ierr);
4136006e8d2SBarry Smith   if (!isblock) {
414f0413b6fSBarry Smith     ierr = ISGetIndices(is,&indices);CHKERRQ(ierr);
415f0413b6fSBarry Smith     ierr = ISLocalToGlobalMappingCreate(comm,1,n,indices,PETSC_COPY_VALUES,mapping);CHKERRQ(ierr);
4162bdab257SBarry Smith     ierr = ISRestoreIndices(is,&indices);CHKERRQ(ierr);
4176006e8d2SBarry Smith   } else {
4186006e8d2SBarry Smith     ierr = ISGetBlockSize(is,&bs);CHKERRQ(ierr);
419f0413b6fSBarry Smith     ierr = ISBlockGetIndices(is,&indices);CHKERRQ(ierr);
42028bc9809SBarry Smith     ierr = ISLocalToGlobalMappingCreate(comm,bs,n/bs,indices,PETSC_COPY_VALUES,mapping);CHKERRQ(ierr);
421f0413b6fSBarry Smith     ierr = ISBlockRestoreIndices(is,&indices);CHKERRQ(ierr);
4226006e8d2SBarry Smith   }
4233a40ed3dSBarry Smith   PetscFunctionReturn(0);
4242bdab257SBarry Smith }
4255a5d4f66SBarry Smith 
426a4d96a55SJed Brown /*@C
427a4d96a55SJed Brown     ISLocalToGlobalMappingCreateSF - Creates a mapping between a local (0 to n)
428a4d96a55SJed Brown     ordering and a global parallel ordering.
429a4d96a55SJed Brown 
430a4d96a55SJed Brown     Collective
431a4d96a55SJed Brown 
432a4d96a55SJed Brown     Input Parameter:
433a4d96a55SJed Brown +   sf - star forest mapping contiguous local indices to (rank, offset)
434a4d96a55SJed Brown -   start - first global index on this process
435a4d96a55SJed Brown 
436a4d96a55SJed Brown     Output Parameter:
437a4d96a55SJed Brown .   mapping - new mapping data structure
438a4d96a55SJed Brown 
439a4d96a55SJed Brown     Level: advanced
440a4d96a55SJed Brown 
4417e99dc12SLawrence Mitchell .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingSetFromOptions()
442a4d96a55SJed Brown @*/
443a4d96a55SJed Brown PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF sf,PetscInt start,ISLocalToGlobalMapping *mapping)
444a4d96a55SJed Brown {
445a4d96a55SJed Brown   PetscErrorCode ierr;
446a4d96a55SJed Brown   PetscInt       i,maxlocal,nroots,nleaves,*globals,*ltog;
447a4d96a55SJed Brown   const PetscInt *ilocal;
448a4d96a55SJed Brown   MPI_Comm       comm;
449a4d96a55SJed Brown 
450a4d96a55SJed Brown   PetscFunctionBegin;
451a4d96a55SJed Brown   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
452a4d96a55SJed Brown   PetscValidPointer(mapping,3);
453a4d96a55SJed Brown 
454a4d96a55SJed Brown   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
4550298fd71SBarry Smith   ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr);
456f6e5521dSKarl Rupp   if (ilocal) {
457f6e5521dSKarl Rupp     for (i=0,maxlocal=0; i<nleaves; i++) maxlocal = PetscMax(maxlocal,ilocal[i]+1);
458f6e5521dSKarl Rupp   }
459a4d96a55SJed Brown   else maxlocal = nleaves;
460785e854fSJed Brown   ierr = PetscMalloc1(nroots,&globals);CHKERRQ(ierr);
461785e854fSJed Brown   ierr = PetscMalloc1(maxlocal,&ltog);CHKERRQ(ierr);
462a4d96a55SJed Brown   for (i=0; i<nroots; i++) globals[i] = start + i;
463a4d96a55SJed Brown   for (i=0; i<maxlocal; i++) ltog[i] = -1;
464ad227feaSJunchao Zhang   ierr = PetscSFBcastBegin(sf,MPIU_INT,globals,ltog,MPI_REPLACE);CHKERRQ(ierr);
465ad227feaSJunchao Zhang   ierr = PetscSFBcastEnd(sf,MPIU_INT,globals,ltog,MPI_REPLACE);CHKERRQ(ierr);
466f0413b6fSBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,1,maxlocal,ltog,PETSC_OWN_POINTER,mapping);CHKERRQ(ierr);
467a4d96a55SJed Brown   ierr = PetscFree(globals);CHKERRQ(ierr);
468a4d96a55SJed Brown   PetscFunctionReturn(0);
469a4d96a55SJed Brown }
470b46b645bSBarry Smith 
47163fa5c83Sstefano_zampini /*@
47263fa5c83Sstefano_zampini     ISLocalToGlobalMappingSetBlockSize - Sets the blocksize of the mapping
47363fa5c83Sstefano_zampini 
47463fa5c83Sstefano_zampini     Not collective
47563fa5c83Sstefano_zampini 
47663fa5c83Sstefano_zampini     Input Parameters:
477a2b725a8SWilliam Gropp +   mapping - mapping data structure
478a2b725a8SWilliam Gropp -   bs - the blocksize
47963fa5c83Sstefano_zampini 
48063fa5c83Sstefano_zampini     Level: advanced
48163fa5c83Sstefano_zampini 
48263fa5c83Sstefano_zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS()
48363fa5c83Sstefano_zampini @*/
48463fa5c83Sstefano_zampini PetscErrorCode  ISLocalToGlobalMappingSetBlockSize(ISLocalToGlobalMapping mapping,PetscInt bs)
48563fa5c83Sstefano_zampini {
486a59f3c4dSstefano_zampini   PetscInt       *nid;
487a59f3c4dSstefano_zampini   const PetscInt *oid;
488a59f3c4dSstefano_zampini   PetscInt       i,cn,on,obs,nn;
48963fa5c83Sstefano_zampini   PetscErrorCode ierr;
49063fa5c83Sstefano_zampini 
49163fa5c83Sstefano_zampini   PetscFunctionBegin;
49263fa5c83Sstefano_zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
49363fa5c83Sstefano_zampini   if (bs < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid block size %D",bs);
49463fa5c83Sstefano_zampini   if (bs == mapping->bs) PetscFunctionReturn(0);
49563fa5c83Sstefano_zampini   on  = mapping->n;
49663fa5c83Sstefano_zampini   obs = mapping->bs;
49763fa5c83Sstefano_zampini   oid = mapping->indices;
49863fa5c83Sstefano_zampini   nn  = (on*obs)/bs;
49963fa5c83Sstefano_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);
500a59f3c4dSstefano_zampini 
50163fa5c83Sstefano_zampini   ierr = PetscMalloc1(nn,&nid);CHKERRQ(ierr);
502a59f3c4dSstefano_zampini   ierr = ISLocalToGlobalMappingGetIndices(mapping,&oid);CHKERRQ(ierr);
503a59f3c4dSstefano_zampini   for (i=0;i<nn;i++) {
504a59f3c4dSstefano_zampini     PetscInt j;
505a59f3c4dSstefano_zampini     for (j=0,cn=0;j<bs-1;j++) {
506a59f3c4dSstefano_zampini       if (oid[i*bs+j] < 0) { cn++; continue; }
507a59f3c4dSstefano_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]);
508a59f3c4dSstefano_zampini     }
509a59f3c4dSstefano_zampini     if (oid[i*bs+j] < 0) cn++;
5108b7cb0e6Sstefano_zampini     if (cn) {
511a59f3c4dSstefano_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);
512a59f3c4dSstefano_zampini       nid[i] = -1;
5138b7cb0e6Sstefano_zampini     } else {
514a59f3c4dSstefano_zampini       nid[i] = oid[i*bs]/bs;
51563fa5c83Sstefano_zampini     }
51663fa5c83Sstefano_zampini   }
517a59f3c4dSstefano_zampini   ierr = ISLocalToGlobalMappingRestoreIndices(mapping,&oid);CHKERRQ(ierr);
518a59f3c4dSstefano_zampini 
51963fa5c83Sstefano_zampini   mapping->n           = nn;
52063fa5c83Sstefano_zampini   mapping->bs          = bs;
52163fa5c83Sstefano_zampini   ierr                 = PetscFree(mapping->indices);CHKERRQ(ierr);
52263fa5c83Sstefano_zampini   mapping->indices     = nid;
523c9345713Sstefano_zampini   mapping->globalstart = 0;
524c9345713Sstefano_zampini   mapping->globalend   = 0;
5251bd0b88eSStefano Zampini 
5261bd0b88eSStefano Zampini   /* reset the cached information */
5271bd0b88eSStefano Zampini   ierr = PetscFree(mapping->info_procs);CHKERRQ(ierr);
5281bd0b88eSStefano Zampini   ierr = PetscFree(mapping->info_numprocs);CHKERRQ(ierr);
5291bd0b88eSStefano Zampini   if (mapping->info_indices) {
5301bd0b88eSStefano Zampini     PetscInt i;
5311bd0b88eSStefano Zampini 
5321bd0b88eSStefano Zampini     ierr = PetscFree((mapping->info_indices)[0]);CHKERRQ(ierr);
5331bd0b88eSStefano Zampini     for (i=1; i<mapping->info_nproc; i++) {
5341bd0b88eSStefano Zampini       ierr = PetscFree(mapping->info_indices[i]);CHKERRQ(ierr);
5351bd0b88eSStefano Zampini     }
5361bd0b88eSStefano Zampini     ierr = PetscFree(mapping->info_indices);CHKERRQ(ierr);
5371bd0b88eSStefano Zampini   }
5381bd0b88eSStefano Zampini   mapping->info_cached = PETSC_FALSE;
5391bd0b88eSStefano Zampini 
540413f72f0SBarry Smith   if (mapping->ops->destroy) {
541413f72f0SBarry Smith     ierr = (*mapping->ops->destroy)(mapping);CHKERRQ(ierr);
542413f72f0SBarry Smith   }
54363fa5c83Sstefano_zampini   PetscFunctionReturn(0);
54463fa5c83Sstefano_zampini }
54563fa5c83Sstefano_zampini 
54645b6f7e9SBarry Smith /*@
54745b6f7e9SBarry Smith     ISLocalToGlobalMappingGetBlockSize - Gets the blocksize of the mapping
54845b6f7e9SBarry Smith     ordering and a global parallel ordering.
54945b6f7e9SBarry Smith 
55045b6f7e9SBarry Smith     Not Collective
55145b6f7e9SBarry Smith 
55245b6f7e9SBarry Smith     Input Parameters:
55345b6f7e9SBarry Smith .   mapping - mapping data structure
55445b6f7e9SBarry Smith 
55545b6f7e9SBarry Smith     Output Parameter:
55645b6f7e9SBarry Smith .   bs - the blocksize
55745b6f7e9SBarry Smith 
55845b6f7e9SBarry Smith     Level: advanced
55945b6f7e9SBarry Smith 
56045b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS()
56145b6f7e9SBarry Smith @*/
56245b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping mapping,PetscInt *bs)
56345b6f7e9SBarry Smith {
56445b6f7e9SBarry Smith   PetscFunctionBegin;
565cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
56645b6f7e9SBarry Smith   *bs = mapping->bs;
56745b6f7e9SBarry Smith   PetscFunctionReturn(0);
56845b6f7e9SBarry Smith }
56945b6f7e9SBarry Smith 
570ba5bb76aSSatish Balay /*@
57190f02eecSBarry Smith     ISLocalToGlobalMappingCreate - Creates a mapping between a local (0 to n)
57290f02eecSBarry Smith     ordering and a global parallel ordering.
5732362add9SBarry Smith 
57489d82c54SBarry Smith     Not Collective, but communicator may have more than one process
575b9cd556bSLois Curfman McInnes 
5762362add9SBarry Smith     Input Parameters:
57789d82c54SBarry Smith +   comm - MPI communicator
578f0413b6fSBarry Smith .   bs - the block size
57928bc9809SBarry Smith .   n - the number of local elements divided by the block size, or equivalently the number of block indices
58028bc9809SBarry 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
581d5ad8652SBarry Smith -   mode - see PetscCopyMode
5822362add9SBarry Smith 
583a997ad1aSLois Curfman McInnes     Output Parameter:
58490f02eecSBarry Smith .   mapping - new mapping data structure
5852362add9SBarry Smith 
58695452b02SPatrick Sanan     Notes:
58795452b02SPatrick Sanan     There is one integer value in indices per block and it represents the actual indices bs*idx + j, where j=0,..,bs-1
588413f72f0SBarry Smith 
5899a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
590413f72f0SBarry 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.
591413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
592413f72f0SBarry Smith 
593a997ad1aSLois Curfman McInnes     Level: advanced
594a997ad1aSLois Curfman McInnes 
595413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingSetFromOptions(), ISLOCALTOGLOBALMAPPINGBASIC, ISLOCALTOGLOBALMAPPINGHASH
596413f72f0SBarry Smith           ISLocalToGlobalMappingSetType(), ISLocalToGlobalMappingType
5972362add9SBarry Smith @*/
59860c7cefcSBarry Smith PetscErrorCode  ISLocalToGlobalMappingCreate(MPI_Comm comm,PetscInt bs,PetscInt n,const PetscInt indices[],PetscCopyMode mode,ISLocalToGlobalMapping *mapping)
5992362add9SBarry Smith {
6006849ba73SBarry Smith   PetscErrorCode ierr;
60132dcc486SBarry Smith   PetscInt       *in;
602b46b645bSBarry Smith 
603b46b645bSBarry Smith   PetscFunctionBegin;
60473911063SBarry Smith   if (n) PetscValidIntPointer(indices,3);
6054482741eSBarry Smith   PetscValidPointer(mapping,4);
606b46b645bSBarry Smith 
6070298fd71SBarry Smith   *mapping = NULL;
608607a6623SBarry Smith   ierr = ISInitializePackage();CHKERRQ(ierr);
6092362add9SBarry Smith 
61073107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(*mapping,IS_LTOGM_CLASSID,"ISLocalToGlobalMapping","Local to global mapping","IS",
61160c7cefcSBarry Smith                            comm,ISLocalToGlobalMappingDestroy,ISLocalToGlobalMappingView);CHKERRQ(ierr);
612d4bb536fSBarry Smith   (*mapping)->n             = n;
613f0413b6fSBarry Smith   (*mapping)->bs            = bs;
614268a049cSStefano Zampini   (*mapping)->info_cached   = PETSC_FALSE;
615268a049cSStefano Zampini   (*mapping)->info_free     = PETSC_FALSE;
616268a049cSStefano Zampini   (*mapping)->info_procs    = NULL;
617268a049cSStefano Zampini   (*mapping)->info_numprocs = NULL;
618268a049cSStefano Zampini   (*mapping)->info_indices  = NULL;
6191bd0b88eSStefano Zampini   (*mapping)->info_nodec    = NULL;
6201bd0b88eSStefano Zampini   (*mapping)->info_nodei    = NULL;
621413f72f0SBarry Smith 
622413f72f0SBarry Smith   (*mapping)->ops->globaltolocalmappingapply      = NULL;
623413f72f0SBarry Smith   (*mapping)->ops->globaltolocalmappingapplyblock = NULL;
624413f72f0SBarry Smith   (*mapping)->ops->destroy                        = NULL;
625d5ad8652SBarry Smith   if (mode == PETSC_COPY_VALUES) {
626785e854fSJed Brown     ierr = PetscMalloc1(n,&in);CHKERRQ(ierr);
627580bdb30SBarry Smith     ierr = PetscArraycpy(in,indices,n);CHKERRQ(ierr);
628d5ad8652SBarry Smith     (*mapping)->indices = in;
6296389a1a1SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt));CHKERRQ(ierr);
6306389a1a1SBarry Smith   } else if (mode == PETSC_OWN_POINTER) {
6316389a1a1SBarry Smith     (*mapping)->indices = (PetscInt*)indices;
6326389a1a1SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt));CHKERRQ(ierr);
6336389a1a1SBarry Smith   }
63460c7cefcSBarry Smith   else SETERRQ(comm,PETSC_ERR_SUP,"Cannot currently use PETSC_USE_POINTER");
6353a40ed3dSBarry Smith   PetscFunctionReturn(0);
6362362add9SBarry Smith }
6372362add9SBarry Smith 
638413f72f0SBarry Smith PetscFunctionList ISLocalToGlobalMappingList = NULL;
639413f72f0SBarry Smith 
640413f72f0SBarry Smith 
64190f02eecSBarry Smith /*@
6427e99dc12SLawrence Mitchell    ISLocalToGlobalMappingSetFromOptions - Set mapping options from the options database.
6437e99dc12SLawrence Mitchell 
6447e99dc12SLawrence Mitchell    Not collective
6457e99dc12SLawrence Mitchell 
6467e99dc12SLawrence Mitchell    Input Parameters:
6477e99dc12SLawrence Mitchell .  mapping - mapping data structure
6487e99dc12SLawrence Mitchell 
6497e99dc12SLawrence Mitchell    Level: advanced
6507e99dc12SLawrence Mitchell 
6517e99dc12SLawrence Mitchell @*/
6527e99dc12SLawrence Mitchell PetscErrorCode ISLocalToGlobalMappingSetFromOptions(ISLocalToGlobalMapping mapping)
6537e99dc12SLawrence Mitchell {
6547e99dc12SLawrence Mitchell   PetscErrorCode             ierr;
655413f72f0SBarry Smith   char                       type[256];
656413f72f0SBarry Smith   ISLocalToGlobalMappingType defaulttype = "Not set";
6577e99dc12SLawrence Mitchell   PetscBool                  flg;
6587e99dc12SLawrence Mitchell 
6597e99dc12SLawrence Mitchell   PetscFunctionBegin;
6607e99dc12SLawrence Mitchell   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
661413f72f0SBarry Smith   ierr = ISLocalToGlobalMappingRegisterAll();CHKERRQ(ierr);
6627e99dc12SLawrence Mitchell   ierr = PetscObjectOptionsBegin((PetscObject)mapping);CHKERRQ(ierr);
663413f72f0SBarry Smith   ierr = PetscOptionsFList("-islocaltoglobalmapping_type","ISLocalToGlobalMapping method","ISLocalToGlobalMappingSetType",ISLocalToGlobalMappingList,(char*)(((PetscObject)mapping)->type_name) ? ((PetscObject)mapping)->type_name : defaulttype,type,256,&flg);CHKERRQ(ierr);
664413f72f0SBarry Smith   if (flg) {
665413f72f0SBarry Smith     ierr = ISLocalToGlobalMappingSetType(mapping,type);CHKERRQ(ierr);
666413f72f0SBarry Smith   }
6677e99dc12SLawrence Mitchell   ierr = PetscOptionsEnd();CHKERRQ(ierr);
6687e99dc12SLawrence Mitchell   PetscFunctionReturn(0);
6697e99dc12SLawrence Mitchell }
6707e99dc12SLawrence Mitchell 
6717e99dc12SLawrence Mitchell /*@
67290f02eecSBarry Smith    ISLocalToGlobalMappingDestroy - Destroys a mapping between a local (0 to n)
67390f02eecSBarry Smith    ordering and a global parallel ordering.
67490f02eecSBarry Smith 
6750f5bd95cSBarry Smith    Note Collective
676b9cd556bSLois Curfman McInnes 
67790f02eecSBarry Smith    Input Parameters:
67890f02eecSBarry Smith .  mapping - mapping data structure
67990f02eecSBarry Smith 
680a997ad1aSLois Curfman McInnes    Level: advanced
681a997ad1aSLois Curfman McInnes 
6823acfe500SLois Curfman McInnes .seealso: ISLocalToGlobalMappingCreate()
68390f02eecSBarry Smith @*/
6846bf464f9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping *mapping)
68590f02eecSBarry Smith {
686dfbe8321SBarry Smith   PetscErrorCode ierr;
6875fd66863SKarl Rupp 
6883a40ed3dSBarry Smith   PetscFunctionBegin;
6896bf464f9SBarry Smith   if (!*mapping) PetscFunctionReturn(0);
6906bf464f9SBarry Smith   PetscValidHeaderSpecific((*mapping),IS_LTOGM_CLASSID,1);
6914c8fdceaSLisandro Dalcin   if (--((PetscObject)(*mapping))->refct > 0) {*mapping = NULL;PetscFunctionReturn(0);}
6926bf464f9SBarry Smith   ierr = PetscFree((*mapping)->indices);CHKERRQ(ierr);
693268a049cSStefano Zampini   ierr = PetscFree((*mapping)->info_procs);CHKERRQ(ierr);
694268a049cSStefano Zampini   ierr = PetscFree((*mapping)->info_numprocs);CHKERRQ(ierr);
695268a049cSStefano Zampini   if ((*mapping)->info_indices) {
696268a049cSStefano Zampini     PetscInt i;
697268a049cSStefano Zampini 
698268a049cSStefano Zampini     ierr = PetscFree(((*mapping)->info_indices)[0]);CHKERRQ(ierr);
699268a049cSStefano Zampini     for (i=1; i<(*mapping)->info_nproc; i++) {
700268a049cSStefano Zampini       ierr = PetscFree(((*mapping)->info_indices)[i]);CHKERRQ(ierr);
701268a049cSStefano Zampini     }
702268a049cSStefano Zampini     ierr = PetscFree((*mapping)->info_indices);CHKERRQ(ierr);
703268a049cSStefano Zampini   }
7041bd0b88eSStefano Zampini   if ((*mapping)->info_nodei) {
7051bd0b88eSStefano Zampini     ierr = PetscFree(((*mapping)->info_nodei)[0]);CHKERRQ(ierr);
7061bd0b88eSStefano Zampini   }
707071fcb05SBarry Smith   ierr = PetscFree2((*mapping)->info_nodec,(*mapping)->info_nodei);CHKERRQ(ierr);
708413f72f0SBarry Smith   if ((*mapping)->ops->destroy) {
709413f72f0SBarry Smith     ierr = (*(*mapping)->ops->destroy)(*mapping);CHKERRQ(ierr);
710413f72f0SBarry Smith   }
711d38fa0fbSBarry Smith   ierr     = PetscHeaderDestroy(mapping);CHKERRQ(ierr);
7124c8fdceaSLisandro Dalcin   *mapping = NULL;
7133a40ed3dSBarry Smith   PetscFunctionReturn(0);
71490f02eecSBarry Smith }
71590f02eecSBarry Smith 
71690f02eecSBarry Smith /*@
7173acfe500SLois Curfman McInnes     ISLocalToGlobalMappingApplyIS - Creates from an IS in the local numbering
7183acfe500SLois Curfman McInnes     a new index set using the global numbering defined in an ISLocalToGlobalMapping
7193acfe500SLois Curfman McInnes     context.
72090f02eecSBarry Smith 
7214cb36875SStefano Zampini     Collective on is
722b9cd556bSLois Curfman McInnes 
72390f02eecSBarry Smith     Input Parameters:
724b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
725b9cd556bSLois Curfman McInnes -   is - index set in local numbering
72690f02eecSBarry Smith 
72790f02eecSBarry Smith     Output Parameters:
72890f02eecSBarry Smith .   newis - index set in global numbering
72990f02eecSBarry Smith 
7304cb36875SStefano Zampini     Notes:
7314cb36875SStefano Zampini     The output IS will have the same communicator of the input IS.
7324cb36875SStefano Zampini 
733a997ad1aSLois Curfman McInnes     Level: advanced
734a997ad1aSLois Curfman McInnes 
73590f02eecSBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),
736d4bb536fSBarry Smith           ISLocalToGlobalMappingDestroy(), ISGlobalToLocalMappingApply()
73790f02eecSBarry Smith @*/
7387087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping,IS is,IS *newis)
73990f02eecSBarry Smith {
7406849ba73SBarry Smith   PetscErrorCode ierr;
741e24637baSBarry Smith   PetscInt       n,*idxout;
7425d0c19d7SBarry Smith   const PetscInt *idxin;
7433a40ed3dSBarry Smith 
7443a40ed3dSBarry Smith   PetscFunctionBegin;
7450700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
7460700a824SBarry Smith   PetscValidHeaderSpecific(is,IS_CLASSID,2);
7474482741eSBarry Smith   PetscValidPointer(newis,3);
74890f02eecSBarry Smith 
7493b9aefa3SBarry Smith   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
75090f02eecSBarry Smith   ierr = ISGetIndices(is,&idxin);CHKERRQ(ierr);
751785e854fSJed Brown   ierr = PetscMalloc1(n,&idxout);CHKERRQ(ierr);
752e24637baSBarry Smith   ierr = ISLocalToGlobalMappingApply(mapping,n,idxin,idxout);CHKERRQ(ierr);
7533b9aefa3SBarry Smith   ierr = ISRestoreIndices(is,&idxin);CHKERRQ(ierr);
754543f3098SMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),n,idxout,PETSC_OWN_POINTER,newis);CHKERRQ(ierr);
7553a40ed3dSBarry Smith   PetscFunctionReturn(0);
75690f02eecSBarry Smith }
75790f02eecSBarry Smith 
758b89cb25eSSatish Balay /*@
7593acfe500SLois Curfman McInnes    ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering
7603acfe500SLois Curfman McInnes    and converts them to the global numbering.
76190f02eecSBarry Smith 
762b9cd556bSLois Curfman McInnes    Not collective
763b9cd556bSLois Curfman McInnes 
764bb25748dSBarry Smith    Input Parameters:
765b9cd556bSLois Curfman McInnes +  mapping - the local to global mapping context
766bb25748dSBarry Smith .  N - number of integers
767b9cd556bSLois Curfman McInnes -  in - input indices in local numbering
768bb25748dSBarry Smith 
769bb25748dSBarry Smith    Output Parameter:
770bb25748dSBarry Smith .  out - indices in global numbering
771bb25748dSBarry Smith 
772b9cd556bSLois Curfman McInnes    Notes:
773b9cd556bSLois Curfman McInnes    The in and out array parameters may be identical.
774d4bb536fSBarry Smith 
775a997ad1aSLois Curfman McInnes    Level: advanced
776a997ad1aSLois Curfman McInnes 
77745b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
7780752156aSBarry Smith           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
779d4bb536fSBarry Smith           AOPetscToApplication(), ISGlobalToLocalMappingApply()
780bb25748dSBarry Smith 
781afcb2eb5SJed Brown @*/
782afcb2eb5SJed Brown PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
783afcb2eb5SJed Brown {
784cbc1caf0SMatthew G. Knepley   PetscInt i,bs,Nmax;
78545b6f7e9SBarry Smith 
78645b6f7e9SBarry Smith   PetscFunctionBegin;
787cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
788cbc1caf0SMatthew G. Knepley   bs   = mapping->bs;
789cbc1caf0SMatthew G. Knepley   Nmax = bs*mapping->n;
79045b6f7e9SBarry Smith   if (bs == 1) {
791cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
79245b6f7e9SBarry Smith     for (i=0; i<N; i++) {
79345b6f7e9SBarry Smith       if (in[i] < 0) {
79445b6f7e9SBarry Smith         out[i] = in[i];
79545b6f7e9SBarry Smith         continue;
79645b6f7e9SBarry Smith       }
797e24637baSBarry 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);
79845b6f7e9SBarry Smith       out[i] = idx[in[i]];
79945b6f7e9SBarry Smith     }
80045b6f7e9SBarry Smith   } else {
801cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
80245b6f7e9SBarry Smith     for (i=0; i<N; i++) {
80345b6f7e9SBarry Smith       if (in[i] < 0) {
80445b6f7e9SBarry Smith         out[i] = in[i];
80545b6f7e9SBarry Smith         continue;
80645b6f7e9SBarry Smith       }
807e24637baSBarry 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);
80845b6f7e9SBarry Smith       out[i] = idx[in[i]/bs]*bs + (in[i] % bs);
80945b6f7e9SBarry Smith     }
81045b6f7e9SBarry Smith   }
81145b6f7e9SBarry Smith   PetscFunctionReturn(0);
81245b6f7e9SBarry Smith }
81345b6f7e9SBarry Smith 
81445b6f7e9SBarry Smith /*@
8156006e8d2SBarry Smith    ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local block numbering and converts them to the global block numbering
81645b6f7e9SBarry Smith 
81745b6f7e9SBarry Smith    Not collective
81845b6f7e9SBarry Smith 
81945b6f7e9SBarry Smith    Input Parameters:
82045b6f7e9SBarry Smith +  mapping - the local to global mapping context
82145b6f7e9SBarry Smith .  N - number of integers
8226006e8d2SBarry Smith -  in - input indices in local block numbering
82345b6f7e9SBarry Smith 
82445b6f7e9SBarry Smith    Output Parameter:
8256006e8d2SBarry Smith .  out - indices in global block numbering
82645b6f7e9SBarry Smith 
82745b6f7e9SBarry Smith    Notes:
82845b6f7e9SBarry Smith    The in and out array parameters may be identical.
82945b6f7e9SBarry Smith 
8306006e8d2SBarry Smith    Example:
8316006e8d2SBarry 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
8326006e8d2SBarry Smith      (the first block) would produce 0 and the mapping applied to 1 (the second block) would produce 3.
8336006e8d2SBarry Smith 
83445b6f7e9SBarry Smith    Level: advanced
83545b6f7e9SBarry Smith 
83645b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
83745b6f7e9SBarry Smith           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
83845b6f7e9SBarry Smith           AOPetscToApplication(), ISGlobalToLocalMappingApply()
83945b6f7e9SBarry Smith 
84045b6f7e9SBarry Smith @*/
84145b6f7e9SBarry Smith PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
84245b6f7e9SBarry Smith {
843cbc1caf0SMatthew G. Knepley 
844cbc1caf0SMatthew G. Knepley   PetscFunctionBegin;
845cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
846cbc1caf0SMatthew G. Knepley   {
847afcb2eb5SJed Brown     PetscInt i,Nmax = mapping->n;
848afcb2eb5SJed Brown     const PetscInt *idx = mapping->indices;
849d4bb536fSBarry Smith 
850afcb2eb5SJed Brown     for (i=0; i<N; i++) {
851afcb2eb5SJed Brown       if (in[i] < 0) {
852afcb2eb5SJed Brown         out[i] = in[i];
853afcb2eb5SJed Brown         continue;
854afcb2eb5SJed Brown       }
855e24637baSBarry 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);
856afcb2eb5SJed Brown       out[i] = idx[in[i]];
857afcb2eb5SJed Brown     }
858cbc1caf0SMatthew G. Knepley   }
859afcb2eb5SJed Brown   PetscFunctionReturn(0);
860afcb2eb5SJed Brown }
861d4bb536fSBarry Smith 
8627e99dc12SLawrence Mitchell /*@
863a997ad1aSLois Curfman McInnes     ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers
864a997ad1aSLois Curfman McInnes     specified with a global numbering.
865d4bb536fSBarry Smith 
866b9cd556bSLois Curfman McInnes     Not collective
867b9cd556bSLois Curfman McInnes 
868d4bb536fSBarry Smith     Input Parameters:
869b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
8700040bde1SJunchao Zhang .   type - IS_GTOLM_MASK - maps global indices with no local value to -1 in the output list (i.e., mask them)
871d4bb536fSBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
872d4bb536fSBarry Smith .   n - number of global indices to map
873b9cd556bSLois Curfman McInnes -   idx - global indices to map
874d4bb536fSBarry Smith 
875d4bb536fSBarry Smith     Output Parameters:
876b9cd556bSLois Curfman McInnes +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
877b9cd556bSLois Curfman McInnes -   idxout - local index of each global index, one must pass in an array long enough
878e182c471SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApply() with
8790298fd71SBarry Smith              idxout == NULL to determine the required length (returned in nout)
880e182c471SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApply()
881e182c471SBarry Smith              a second time to set the values.
882d4bb536fSBarry Smith 
883b9cd556bSLois Curfman McInnes     Notes:
8840298fd71SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
885d4bb536fSBarry Smith 
8869a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
887413f72f0SBarry 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.
888413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
8890f5bd95cSBarry Smith 
890a997ad1aSLois Curfman McInnes     Level: advanced
891a997ad1aSLois Curfman McInnes 
89232fd6b96SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
89332fd6b96SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
89432fd6b96SBarry Smith 
8959d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),
896413f72f0SBarry Smith           ISLocalToGlobalMappingDestroy()
897d4bb536fSBarry Smith @*/
898413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
899d4bb536fSBarry Smith {
9009d90f715SBarry Smith   PetscErrorCode ierr;
9019d90f715SBarry Smith 
9029d90f715SBarry Smith   PetscFunctionBegin;
9039d90f715SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
904413f72f0SBarry Smith   if (!mapping->data) {
905413f72f0SBarry Smith     ierr = ISGlobalToLocalMappingSetUp(mapping);CHKERRQ(ierr);
9069d90f715SBarry Smith   }
907413f72f0SBarry Smith   ierr = (*mapping->ops->globaltolocalmappingapply)(mapping,type,n,idx,nout,idxout);CHKERRQ(ierr);
9089d90f715SBarry Smith   PetscFunctionReturn(0);
9099d90f715SBarry Smith }
9109d90f715SBarry Smith 
911d4fe737eSStefano Zampini /*@
912d4fe737eSStefano Zampini     ISGlobalToLocalMappingApplyIS - Creates from an IS in the global numbering
913d4fe737eSStefano Zampini     a new index set using the local numbering defined in an ISLocalToGlobalMapping
914d4fe737eSStefano Zampini     context.
915d4fe737eSStefano Zampini 
916d4fe737eSStefano Zampini     Not collective
917d4fe737eSStefano Zampini 
918d4fe737eSStefano Zampini     Input Parameters:
919d4fe737eSStefano Zampini +   mapping - mapping between local and global numbering
9200040bde1SJunchao Zhang .   type - IS_GTOLM_MASK - maps global indices with no local value to -1 in the output list (i.e., mask them)
9212785b321SStefano Zampini            IS_GTOLM_DROP - drops the indices with no local value from the output list
922d4fe737eSStefano Zampini -   is - index set in global numbering
923d4fe737eSStefano Zampini 
924d4fe737eSStefano Zampini     Output Parameters:
925d4fe737eSStefano Zampini .   newis - index set in local numbering
926d4fe737eSStefano Zampini 
9274cb36875SStefano Zampini     Notes:
9284cb36875SStefano Zampini     The output IS will be sequential, as it encodes a purely local operation
9294cb36875SStefano Zampini 
930d4fe737eSStefano Zampini     Level: advanced
931d4fe737eSStefano Zampini 
932d4fe737eSStefano Zampini .seealso: ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(),
933d4fe737eSStefano Zampini           ISLocalToGlobalMappingDestroy()
934d4fe737eSStefano Zampini @*/
935413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,IS is,IS *newis)
936d4fe737eSStefano Zampini {
937d4fe737eSStefano Zampini   PetscErrorCode ierr;
938d4fe737eSStefano Zampini   PetscInt       n,nout,*idxout;
939d4fe737eSStefano Zampini   const PetscInt *idxin;
940d4fe737eSStefano Zampini 
941d4fe737eSStefano Zampini   PetscFunctionBegin;
942d4fe737eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
943d4fe737eSStefano Zampini   PetscValidHeaderSpecific(is,IS_CLASSID,3);
944d4fe737eSStefano Zampini   PetscValidPointer(newis,4);
945d4fe737eSStefano Zampini 
946d4fe737eSStefano Zampini   ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr);
947d4fe737eSStefano Zampini   ierr = ISGetIndices(is,&idxin);CHKERRQ(ierr);
948d4fe737eSStefano Zampini   if (type == IS_GTOLM_MASK) {
949d4fe737eSStefano Zampini     ierr = PetscMalloc1(n,&idxout);CHKERRQ(ierr);
950d4fe737eSStefano Zampini   } else {
951d4fe737eSStefano Zampini     ierr = ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,NULL);CHKERRQ(ierr);
952d4fe737eSStefano Zampini     ierr = PetscMalloc1(nout,&idxout);CHKERRQ(ierr);
953d4fe737eSStefano Zampini   }
954d4fe737eSStefano Zampini   ierr = ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,idxout);CHKERRQ(ierr);
955d4fe737eSStefano Zampini   ierr = ISRestoreIndices(is,&idxin);CHKERRQ(ierr);
9564cb36875SStefano Zampini   ierr = ISCreateGeneral(PETSC_COMM_SELF,nout,idxout,PETSC_OWN_POINTER,newis);CHKERRQ(ierr);
957d4fe737eSStefano Zampini   PetscFunctionReturn(0);
958d4fe737eSStefano Zampini }
959d4fe737eSStefano Zampini 
9609d90f715SBarry Smith /*@
9619d90f715SBarry Smith     ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers
9629d90f715SBarry Smith     specified with a block global numbering.
9639d90f715SBarry Smith 
9649d90f715SBarry Smith     Not collective
9659d90f715SBarry Smith 
9669d90f715SBarry Smith     Input Parameters:
9679d90f715SBarry Smith +   mapping - mapping between local and global numbering
9680040bde1SJunchao Zhang .   type - IS_GTOLM_MASK - maps global indices with no local value to -1 in the output list (i.e., mask them)
9699d90f715SBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
9709d90f715SBarry Smith .   n - number of global indices to map
9719d90f715SBarry Smith -   idx - global indices to map
9729d90f715SBarry Smith 
9739d90f715SBarry Smith     Output Parameters:
9749d90f715SBarry Smith +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
9759d90f715SBarry Smith -   idxout - local index of each global index, one must pass in an array long enough
9769d90f715SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApplyBlock() with
9779d90f715SBarry Smith              idxout == NULL to determine the required length (returned in nout)
9789d90f715SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApplyBlock()
9799d90f715SBarry Smith              a second time to set the values.
9809d90f715SBarry Smith 
9819d90f715SBarry Smith     Notes:
9829d90f715SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
9839d90f715SBarry Smith 
9849a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
985413f72f0SBarry 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.
986413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
9879d90f715SBarry Smith 
9889d90f715SBarry Smith     Level: advanced
9899d90f715SBarry Smith 
9909d90f715SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
9919d90f715SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
9929d90f715SBarry Smith 
9939d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(),
994413f72f0SBarry Smith           ISLocalToGlobalMappingDestroy()
9959d90f715SBarry Smith @*/
996413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,
9979d90f715SBarry Smith                                                  PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
9989d90f715SBarry Smith {
9996849ba73SBarry Smith   PetscErrorCode ierr;
1000d4bb536fSBarry Smith 
10013a40ed3dSBarry Smith   PetscFunctionBegin;
10020700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1003413f72f0SBarry Smith   if (!mapping->data) {
1004413f72f0SBarry Smith     ierr = ISGlobalToLocalMappingSetUp(mapping);CHKERRQ(ierr);
1005d4bb536fSBarry Smith   }
1006413f72f0SBarry Smith   ierr = (*mapping->ops->globaltolocalmappingapplyblock)(mapping,type,n,idx,nout,idxout);CHKERRQ(ierr);
10073a40ed3dSBarry Smith   PetscFunctionReturn(0);
1008d4bb536fSBarry Smith }
100990f02eecSBarry Smith 
1010413f72f0SBarry Smith 
101189d82c54SBarry Smith /*@C
10126a818285SBarry Smith     ISLocalToGlobalMappingGetBlockInfo - Gets the neighbor information for each processor and
101389d82c54SBarry Smith      each index shared by more than one processor
101489d82c54SBarry Smith 
101589d82c54SBarry Smith     Collective on ISLocalToGlobalMapping
101689d82c54SBarry Smith 
101789d82c54SBarry Smith     Input Parameters:
101889d82c54SBarry Smith .   mapping - the mapping from local to global indexing
101989d82c54SBarry Smith 
102089d82c54SBarry Smith     Output Parameter:
102189d82c54SBarry Smith +   nproc - number of processors that are connected to this one
102289d82c54SBarry Smith .   proc - neighboring processors
102307b52d57SBarry Smith .   numproc - number of indices for each subdomain (processor)
10243463a7baSJed Brown -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
102589d82c54SBarry Smith 
102689d82c54SBarry Smith     Level: advanced
102789d82c54SBarry Smith 
10282cfcea29SBarry Smith     Fortran Usage:
10292cfcea29SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
10302cfcea29SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
10312cfcea29SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
10322cfcea29SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
10332cfcea29SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
10342cfcea29SBarry Smith 
10352cfcea29SBarry Smith 
103607b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
103707b52d57SBarry Smith           ISLocalToGlobalMappingRestoreInfo()
103889d82c54SBarry Smith @*/
10396a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
104089d82c54SBarry Smith {
10416849ba73SBarry Smith   PetscErrorCode ierr;
1042268a049cSStefano Zampini 
1043268a049cSStefano Zampini   PetscFunctionBegin;
1044268a049cSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1045268a049cSStefano Zampini   if (mapping->info_cached) {
1046268a049cSStefano Zampini     *nproc    = mapping->info_nproc;
1047268a049cSStefano Zampini     *procs    = mapping->info_procs;
1048268a049cSStefano Zampini     *numprocs = mapping->info_numprocs;
1049268a049cSStefano Zampini     *indices  = mapping->info_indices;
1050268a049cSStefano Zampini   } else {
1051268a049cSStefano Zampini     ierr = ISLocalToGlobalMappingGetBlockInfo_Private(mapping,nproc,procs,numprocs,indices);CHKERRQ(ierr);
1052268a049cSStefano Zampini   }
1053268a049cSStefano Zampini   PetscFunctionReturn(0);
1054268a049cSStefano Zampini }
1055268a049cSStefano Zampini 
1056268a049cSStefano Zampini static PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
1057268a049cSStefano Zampini {
1058268a049cSStefano Zampini   PetscErrorCode ierr;
105997f1f81fSBarry Smith   PetscMPIInt    size,rank,tag1,tag2,tag3,*len,*source,imdex;
106032dcc486SBarry Smith   PetscInt       i,n = mapping->n,Ng,ng,max = 0,*lindices = mapping->indices;
106132dcc486SBarry Smith   PetscInt       *nprocs,*owner,nsends,*sends,j,*starts,nmax,nrecvs,*recvs,proc;
106297f1f81fSBarry Smith   PetscInt       cnt,scale,*ownedsenders,*nownedsenders,rstart,nowned;
106332dcc486SBarry Smith   PetscInt       node,nownedm,nt,*sends2,nsends2,*starts2,*lens2,*dest,nrecvs2,*starts3,*recvs2,k,*bprocs,*tmp;
106432dcc486SBarry Smith   PetscInt       first_procs,first_numprocs,*first_indices;
106589d82c54SBarry Smith   MPI_Request    *recv_waits,*send_waits;
106630dcb7c9SBarry Smith   MPI_Status     recv_status,*send_status,*recv_statuses;
1067ce94432eSBarry Smith   MPI_Comm       comm;
1068ace3abfcSBarry Smith   PetscBool      debug = PETSC_FALSE;
106989d82c54SBarry Smith 
107089d82c54SBarry Smith   PetscFunctionBegin;
1071ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)mapping,&comm);CHKERRQ(ierr);
1072ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
1073ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
107424cf384cSBarry Smith   if (size == 1) {
107524cf384cSBarry Smith     *nproc         = 0;
10760298fd71SBarry Smith     *procs         = NULL;
107795dccacaSBarry Smith     ierr           = PetscNew(numprocs);CHKERRQ(ierr);
10781e2105dcSBarry Smith     (*numprocs)[0] = 0;
107995dccacaSBarry Smith     ierr           = PetscNew(indices);CHKERRQ(ierr);
10800298fd71SBarry Smith     (*indices)[0]  = NULL;
1081268a049cSStefano Zampini     /* save info for reuse */
1082268a049cSStefano Zampini     mapping->info_nproc = *nproc;
1083268a049cSStefano Zampini     mapping->info_procs = *procs;
1084268a049cSStefano Zampini     mapping->info_numprocs = *numprocs;
1085268a049cSStefano Zampini     mapping->info_indices = *indices;
1086268a049cSStefano Zampini     mapping->info_cached = PETSC_TRUE;
108724cf384cSBarry Smith     PetscFunctionReturn(0);
108824cf384cSBarry Smith   }
108924cf384cSBarry Smith 
1090c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)mapping)->options,NULL,"-islocaltoglobalmappinggetinfo_debug",&debug,NULL);CHKERRQ(ierr);
109107b52d57SBarry Smith 
10923677ff5aSBarry Smith   /*
10936a818285SBarry Smith     Notes on ISLocalToGlobalMappingGetBlockInfo
10943677ff5aSBarry Smith 
10953677ff5aSBarry Smith     globally owned node - the nodes that have been assigned to this processor in global
10963677ff5aSBarry Smith            numbering, just for this routine.
10973677ff5aSBarry Smith 
10983677ff5aSBarry Smith     nontrivial globally owned node - node assigned to this processor that is on a subdomain
10993677ff5aSBarry Smith            boundary (i.e. is has more than one local owner)
11003677ff5aSBarry Smith 
11013677ff5aSBarry Smith     locally owned node - node that exists on this processors subdomain
11023677ff5aSBarry Smith 
11033677ff5aSBarry Smith     nontrivial locally owned node - node that is not in the interior (i.e. has more than one
11043677ff5aSBarry Smith            local subdomain
11053677ff5aSBarry Smith   */
110624cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag1);CHKERRQ(ierr);
110724cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag2);CHKERRQ(ierr);
110824cf384cSBarry Smith   ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag3);CHKERRQ(ierr);
110989d82c54SBarry Smith 
111089d82c54SBarry Smith   for (i=0; i<n; i++) {
111189d82c54SBarry Smith     if (lindices[i] > max) max = lindices[i];
111289d82c54SBarry Smith   }
1113*820f2d46SBarry Smith   ierr   = MPIU_Allreduce(&max,&Ng,1,MPIU_INT,MPI_MAX,comm);CHKERRMPI(ierr);
111478058e43SBarry Smith   Ng++;
111555b25c41SPierre Jolivet   ierr   = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
111655b25c41SPierre Jolivet   ierr   = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
1117bc8ff85bSBarry Smith   scale  = Ng/size + 1;
1118a2e34c3dSBarry Smith   ng     = scale; if (rank == size-1) ng = Ng - scale*(size-1); ng = PetscMax(1,ng);
1119caba0dd0SBarry Smith   rstart = scale*rank;
112089d82c54SBarry Smith 
112189d82c54SBarry Smith   /* determine ownership ranges of global indices */
1122785e854fSJed Brown   ierr = PetscMalloc1(2*size,&nprocs);CHKERRQ(ierr);
1123580bdb30SBarry Smith   ierr = PetscArrayzero(nprocs,2*size);CHKERRQ(ierr);
112489d82c54SBarry Smith 
112589d82c54SBarry Smith   /* determine owners of each local node  */
1126785e854fSJed Brown   ierr = PetscMalloc1(n,&owner);CHKERRQ(ierr);
112789d82c54SBarry Smith   for (i=0; i<n; i++) {
11283677ff5aSBarry Smith     proc             = lindices[i]/scale; /* processor that globally owns this index */
112927c402fcSBarry Smith     nprocs[2*proc+1] = 1;                 /* processor globally owns at least one of ours */
11303677ff5aSBarry Smith     owner[i]         = proc;
113127c402fcSBarry Smith     nprocs[2*proc]++;                     /* count of how many that processor globally owns of ours */
113289d82c54SBarry Smith   }
113327c402fcSBarry Smith   nsends = 0; for (i=0; i<size; i++) nsends += nprocs[2*i+1];
11347904a332SBarry Smith   ierr = PetscInfo1(mapping,"Number of global owners for my local data %D\n",nsends);CHKERRQ(ierr);
113589d82c54SBarry Smith 
113689d82c54SBarry Smith   /* inform other processors of number of messages and max length*/
113727c402fcSBarry Smith   ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr);
11387904a332SBarry Smith   ierr = PetscInfo1(mapping,"Number of local owners for my global data %D\n",nrecvs);CHKERRQ(ierr);
113989d82c54SBarry Smith 
114089d82c54SBarry Smith   /* post receives for owned rows */
1141785e854fSJed Brown   ierr = PetscMalloc1((2*nrecvs+1)*(nmax+1),&recvs);CHKERRQ(ierr);
1142854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&recv_waits);CHKERRQ(ierr);
114389d82c54SBarry Smith   for (i=0; i<nrecvs; i++) {
1144ffc4695bSBarry Smith     ierr = MPI_Irecv(recvs+2*nmax*i,2*nmax,MPIU_INT,MPI_ANY_SOURCE,tag1,comm,recv_waits+i);CHKERRMPI(ierr);
114589d82c54SBarry Smith   }
114689d82c54SBarry Smith 
114789d82c54SBarry Smith   /* pack messages containing lists of local nodes to owners */
1148854ce69bSBarry Smith   ierr      = PetscMalloc1(2*n+1,&sends);CHKERRQ(ierr);
1149854ce69bSBarry Smith   ierr      = PetscMalloc1(size+1,&starts);CHKERRQ(ierr);
115089d82c54SBarry Smith   starts[0] = 0;
1151f6e5521dSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
115289d82c54SBarry Smith   for (i=0; i<n; i++) {
115389d82c54SBarry Smith     sends[starts[owner[i]]++] = lindices[i];
115430dcb7c9SBarry Smith     sends[starts[owner[i]]++] = i;
115589d82c54SBarry Smith   }
115689d82c54SBarry Smith   ierr = PetscFree(owner);CHKERRQ(ierr);
115789d82c54SBarry Smith   starts[0] = 0;
1158f6e5521dSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
115989d82c54SBarry Smith 
116089d82c54SBarry Smith   /* send the messages */
1161854ce69bSBarry Smith   ierr = PetscMalloc1(nsends+1,&send_waits);CHKERRQ(ierr);
1162854ce69bSBarry Smith   ierr = PetscMalloc1(nsends+1,&dest);CHKERRQ(ierr);
116389d82c54SBarry Smith   cnt = 0;
116489d82c54SBarry Smith   for (i=0; i<size; i++) {
116527c402fcSBarry Smith     if (nprocs[2*i]) {
116655b25c41SPierre Jolivet       ierr      = MPI_Isend(sends+starts[i],2*nprocs[2*i],MPIU_INT,i,tag1,comm,send_waits+cnt);CHKERRMPI(ierr);
116730dcb7c9SBarry Smith       dest[cnt] = i;
116889d82c54SBarry Smith       cnt++;
116989d82c54SBarry Smith     }
117089d82c54SBarry Smith   }
117189d82c54SBarry Smith   ierr = PetscFree(starts);CHKERRQ(ierr);
117289d82c54SBarry Smith 
117389d82c54SBarry Smith   /* wait on receives */
1174854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&source);CHKERRQ(ierr);
1175854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&len);CHKERRQ(ierr);
117689d82c54SBarry Smith   cnt  = nrecvs;
1177580bdb30SBarry Smith   ierr = PetscCalloc1(ng+1,&nownedsenders);CHKERRQ(ierr);
117889d82c54SBarry Smith   while (cnt) {
1179ffc4695bSBarry Smith     ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRMPI(ierr);
118089d82c54SBarry Smith     /* unpack receives into our local space */
118155b25c41SPierre Jolivet     ierr          = MPI_Get_count(&recv_status,MPIU_INT,&len[imdex]);CHKERRMPI(ierr);
118289d82c54SBarry Smith     source[imdex] = recv_status.MPI_SOURCE;
118330dcb7c9SBarry Smith     len[imdex]    = len[imdex]/2;
1184caba0dd0SBarry Smith     /* count how many local owners for each of my global owned indices */
118530dcb7c9SBarry Smith     for (i=0; i<len[imdex]; i++) nownedsenders[recvs[2*imdex*nmax+2*i]-rstart]++;
118689d82c54SBarry Smith     cnt--;
118789d82c54SBarry Smith   }
118889d82c54SBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
118989d82c54SBarry Smith 
119030dcb7c9SBarry Smith   /* count how many globally owned indices are on an edge multiplied by how many processors own them. */
1191bc8ff85bSBarry Smith   nowned  = 0;
1192bc8ff85bSBarry Smith   nownedm = 0;
1193bc8ff85bSBarry Smith   for (i=0; i<ng; i++) {
1194bc8ff85bSBarry Smith     if (nownedsenders[i] > 1) {nownedm += nownedsenders[i]; nowned++;}
1195bc8ff85bSBarry Smith   }
1196bc8ff85bSBarry Smith 
1197bc8ff85bSBarry Smith   /* create single array to contain rank of all local owners of each globally owned index */
1198854ce69bSBarry Smith   ierr      = PetscMalloc1(nownedm+1,&ownedsenders);CHKERRQ(ierr);
1199854ce69bSBarry Smith   ierr      = PetscMalloc1(ng+1,&starts);CHKERRQ(ierr);
1200bc8ff85bSBarry Smith   starts[0] = 0;
1201bc8ff85bSBarry Smith   for (i=1; i<ng; i++) {
1202bc8ff85bSBarry Smith     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
1203bc8ff85bSBarry Smith     else starts[i] = starts[i-1];
1204bc8ff85bSBarry Smith   }
1205bc8ff85bSBarry Smith 
120630dcb7c9SBarry Smith   /* for each nontrival globally owned node list all arriving processors */
1207bc8ff85bSBarry Smith   for (i=0; i<nrecvs; i++) {
1208bc8ff85bSBarry Smith     for (j=0; j<len[i]; j++) {
120930dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
1210f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) ownedsenders[starts[node]++] = source[i];
1211bc8ff85bSBarry Smith     }
1212bc8ff85bSBarry Smith   }
1213bc8ff85bSBarry Smith 
121407b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
121530dcb7c9SBarry Smith     starts[0] = 0;
121630dcb7c9SBarry Smith     for (i=1; i<ng; i++) {
121730dcb7c9SBarry Smith       if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
121830dcb7c9SBarry Smith       else starts[i] = starts[i-1];
121930dcb7c9SBarry Smith     }
122030dcb7c9SBarry Smith     for (i=0; i<ng; i++) {
122130dcb7c9SBarry Smith       if (nownedsenders[i] > 1) {
12227904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"[%d] global node %D local owner processors: ",rank,i+rstart);CHKERRQ(ierr);
122330dcb7c9SBarry Smith         for (j=0; j<nownedsenders[i]; j++) {
12247904a332SBarry Smith           ierr = PetscSynchronizedPrintf(comm,"%D ",ownedsenders[starts[i]+j]);CHKERRQ(ierr);
122530dcb7c9SBarry Smith         }
122630dcb7c9SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
122730dcb7c9SBarry Smith       }
122830dcb7c9SBarry Smith     }
12290ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
123007b52d57SBarry Smith   } /* -----------------------------------  */
123130dcb7c9SBarry Smith 
12323677ff5aSBarry Smith   /* wait on original sends */
12333a96401aSBarry Smith   if (nsends) {
1234785e854fSJed Brown     ierr = PetscMalloc1(nsends,&send_status);CHKERRQ(ierr);
1235ffc4695bSBarry Smith     ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRMPI(ierr);
12363a96401aSBarry Smith     ierr = PetscFree(send_status);CHKERRQ(ierr);
12373a96401aSBarry Smith   }
123889d82c54SBarry Smith   ierr = PetscFree(send_waits);CHKERRQ(ierr);
12393a96401aSBarry Smith   ierr = PetscFree(sends);CHKERRQ(ierr);
12403677ff5aSBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
12413677ff5aSBarry Smith 
12423677ff5aSBarry Smith   /* pack messages to send back to local owners */
124330dcb7c9SBarry Smith   starts[0] = 0;
124430dcb7c9SBarry Smith   for (i=1; i<ng; i++) {
124530dcb7c9SBarry Smith     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
124630dcb7c9SBarry Smith     else starts[i] = starts[i-1];
124730dcb7c9SBarry Smith   }
124830dcb7c9SBarry Smith   nsends2 = nrecvs;
1249854ce69bSBarry Smith   ierr    = PetscMalloc1(nsends2+1,&nprocs);CHKERRQ(ierr); /* length of each message */
125030dcb7c9SBarry Smith   for (i=0; i<nrecvs; i++) {
125130dcb7c9SBarry Smith     nprocs[i] = 1;
125230dcb7c9SBarry Smith     for (j=0; j<len[i]; j++) {
125330dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
1254f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) nprocs[i] += 2 + nownedsenders[node];
125530dcb7c9SBarry Smith     }
125630dcb7c9SBarry Smith   }
1257f6e5521dSKarl Rupp   nt = 0;
1258f6e5521dSKarl Rupp   for (i=0; i<nsends2; i++) nt += nprocs[i];
1259f6e5521dSKarl Rupp 
1260854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,&sends2);CHKERRQ(ierr);
1261854ce69bSBarry Smith   ierr = PetscMalloc1(nsends2+1,&starts2);CHKERRQ(ierr);
1262f6e5521dSKarl Rupp 
1263f6e5521dSKarl Rupp   starts2[0] = 0;
1264f6e5521dSKarl Rupp   for (i=1; i<nsends2; i++) starts2[i] = starts2[i-1] + nprocs[i-1];
126530dcb7c9SBarry Smith   /*
126630dcb7c9SBarry Smith      Each message is 1 + nprocs[i] long, and consists of
126730dcb7c9SBarry Smith        (0) the number of nodes being sent back
126830dcb7c9SBarry Smith        (1) the local node number,
126930dcb7c9SBarry Smith        (2) the number of processors sharing it,
127030dcb7c9SBarry Smith        (3) the processors sharing it
127130dcb7c9SBarry Smith   */
127230dcb7c9SBarry Smith   for (i=0; i<nsends2; i++) {
127330dcb7c9SBarry Smith     cnt = 1;
127430dcb7c9SBarry Smith     sends2[starts2[i]] = 0;
127530dcb7c9SBarry Smith     for (j=0; j<len[i]; j++) {
127630dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
127730dcb7c9SBarry Smith       if (nownedsenders[node] > 1) {
127830dcb7c9SBarry Smith         sends2[starts2[i]]++;
127930dcb7c9SBarry Smith         sends2[starts2[i]+cnt++] = recvs[2*i*nmax+2*j+1];
128030dcb7c9SBarry Smith         sends2[starts2[i]+cnt++] = nownedsenders[node];
1281580bdb30SBarry Smith         ierr = PetscArraycpy(&sends2[starts2[i]+cnt],&ownedsenders[starts[node]],nownedsenders[node]);CHKERRQ(ierr);
128230dcb7c9SBarry Smith         cnt += nownedsenders[node];
128330dcb7c9SBarry Smith       }
128430dcb7c9SBarry Smith     }
128530dcb7c9SBarry Smith   }
128630dcb7c9SBarry Smith 
128730dcb7c9SBarry Smith   /* receive the message lengths */
128830dcb7c9SBarry Smith   nrecvs2 = nsends;
1289854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&lens2);CHKERRQ(ierr);
1290854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&starts3);CHKERRQ(ierr);
1291854ce69bSBarry Smith   ierr    = PetscMalloc1(nrecvs2+1,&recv_waits);CHKERRQ(ierr);
129230dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
1293ffc4695bSBarry Smith     ierr = MPI_Irecv(&lens2[i],1,MPIU_INT,dest[i],tag2,comm,recv_waits+i);CHKERRMPI(ierr);
129430dcb7c9SBarry Smith   }
1295d44834fbSBarry Smith 
12968a8e0b3aSBarry Smith   /* send the message lengths */
12978a8e0b3aSBarry Smith   for (i=0; i<nsends2; i++) {
1298ffc4695bSBarry Smith     ierr = MPI_Send(&nprocs[i],1,MPIU_INT,source[i],tag2,comm);CHKERRMPI(ierr);
12998a8e0b3aSBarry Smith   }
13008a8e0b3aSBarry Smith 
1301d44834fbSBarry Smith   /* wait on receives of lens */
13020c468ba9SBarry Smith   if (nrecvs2) {
1303785e854fSJed Brown     ierr = PetscMalloc1(nrecvs2,&recv_statuses);CHKERRQ(ierr);
1304ffc4695bSBarry Smith     ierr = MPI_Waitall(nrecvs2,recv_waits,recv_statuses);CHKERRMPI(ierr);
1305d44834fbSBarry Smith     ierr = PetscFree(recv_statuses);CHKERRQ(ierr);
13060c468ba9SBarry Smith   }
1307a2ea699eSBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
1308d44834fbSBarry Smith 
130930dcb7c9SBarry Smith   starts3[0] = 0;
1310d44834fbSBarry Smith   nt         = 0;
131130dcb7c9SBarry Smith   for (i=0; i<nrecvs2-1; i++) {
131230dcb7c9SBarry Smith     starts3[i+1] = starts3[i] + lens2[i];
1313d44834fbSBarry Smith     nt          += lens2[i];
131430dcb7c9SBarry Smith   }
131576466f69SStefano Zampini   if (nrecvs2) nt += lens2[nrecvs2-1];
1316d44834fbSBarry Smith 
1317854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,&recvs2);CHKERRQ(ierr);
1318854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs2+1,&recv_waits);CHKERRQ(ierr);
131952b72c4aSBarry Smith   for (i=0; i<nrecvs2; i++) {
1320ffc4695bSBarry Smith     ierr = MPI_Irecv(recvs2+starts3[i],lens2[i],MPIU_INT,dest[i],tag3,comm,recv_waits+i);CHKERRMPI(ierr);
132130dcb7c9SBarry Smith   }
132230dcb7c9SBarry Smith 
132330dcb7c9SBarry Smith   /* send the messages */
1324854ce69bSBarry Smith   ierr = PetscMalloc1(nsends2+1,&send_waits);CHKERRQ(ierr);
132530dcb7c9SBarry Smith   for (i=0; i<nsends2; i++) {
1326ffc4695bSBarry Smith     ierr = MPI_Isend(sends2+starts2[i],nprocs[i],MPIU_INT,source[i],tag3,comm,send_waits+i);CHKERRMPI(ierr);
132730dcb7c9SBarry Smith   }
132830dcb7c9SBarry Smith 
132930dcb7c9SBarry Smith   /* wait on receives */
13300c468ba9SBarry Smith   if (nrecvs2) {
1331785e854fSJed Brown     ierr = PetscMalloc1(nrecvs2,&recv_statuses);CHKERRQ(ierr);
1332ffc4695bSBarry Smith     ierr = MPI_Waitall(nrecvs2,recv_waits,recv_statuses);CHKERRMPI(ierr);
133330dcb7c9SBarry Smith     ierr = PetscFree(recv_statuses);CHKERRQ(ierr);
13340c468ba9SBarry Smith   }
133530dcb7c9SBarry Smith   ierr = PetscFree(recv_waits);CHKERRQ(ierr);
133630dcb7c9SBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
133730dcb7c9SBarry Smith 
133807b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
133930dcb7c9SBarry Smith     cnt = 0;
134030dcb7c9SBarry Smith     for (i=0; i<nrecvs2; i++) {
134130dcb7c9SBarry Smith       nt = recvs2[cnt++];
134230dcb7c9SBarry Smith       for (j=0; j<nt; j++) {
13437904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"[%d] local node %D number of subdomains %D: ",rank,recvs2[cnt],recvs2[cnt+1]);CHKERRQ(ierr);
134430dcb7c9SBarry Smith         for (k=0; k<recvs2[cnt+1]; k++) {
13457904a332SBarry Smith           ierr = PetscSynchronizedPrintf(comm,"%D ",recvs2[cnt+2+k]);CHKERRQ(ierr);
134630dcb7c9SBarry Smith         }
134730dcb7c9SBarry Smith         cnt += 2 + recvs2[cnt+1];
134830dcb7c9SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
134930dcb7c9SBarry Smith       }
135030dcb7c9SBarry Smith     }
13510ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
135207b52d57SBarry Smith   } /* -----------------------------------  */
135330dcb7c9SBarry Smith 
135430dcb7c9SBarry Smith   /* count number subdomains for each local node */
1355580bdb30SBarry Smith   ierr = PetscCalloc1(size,&nprocs);CHKERRQ(ierr);
135630dcb7c9SBarry Smith   cnt  = 0;
135730dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
135830dcb7c9SBarry Smith     nt = recvs2[cnt++];
135930dcb7c9SBarry Smith     for (j=0; j<nt; j++) {
1360f6e5521dSKarl Rupp       for (k=0; k<recvs2[cnt+1]; k++) nprocs[recvs2[cnt+2+k]]++;
136130dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt+1];
136230dcb7c9SBarry Smith     }
136330dcb7c9SBarry Smith   }
136430dcb7c9SBarry Smith   nt = 0; for (i=0; i<size; i++) nt += (nprocs[i] > 0);
136530dcb7c9SBarry Smith   *nproc    = nt;
1366854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,procs);CHKERRQ(ierr);
1367854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,numprocs);CHKERRQ(ierr);
1368854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,indices);CHKERRQ(ierr);
13690298fd71SBarry Smith   for (i=0;i<nt+1;i++) (*indices)[i]=NULL;
1370785e854fSJed Brown   ierr = PetscMalloc1(size,&bprocs);CHKERRQ(ierr);
137130dcb7c9SBarry Smith   cnt  = 0;
137230dcb7c9SBarry Smith   for (i=0; i<size; i++) {
137330dcb7c9SBarry Smith     if (nprocs[i] > 0) {
137430dcb7c9SBarry Smith       bprocs[i]        = cnt;
137530dcb7c9SBarry Smith       (*procs)[cnt]    = i;
137630dcb7c9SBarry Smith       (*numprocs)[cnt] = nprocs[i];
1377785e854fSJed Brown       ierr             = PetscMalloc1(nprocs[i],&(*indices)[cnt]);CHKERRQ(ierr);
137830dcb7c9SBarry Smith       cnt++;
137930dcb7c9SBarry Smith     }
138030dcb7c9SBarry Smith   }
138130dcb7c9SBarry Smith 
138230dcb7c9SBarry Smith   /* make the list of subdomains for each nontrivial local node */
1383580bdb30SBarry Smith   ierr = PetscArrayzero(*numprocs,nt);CHKERRQ(ierr);
138430dcb7c9SBarry Smith   cnt  = 0;
138530dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
138630dcb7c9SBarry Smith     nt = recvs2[cnt++];
138730dcb7c9SBarry Smith     for (j=0; j<nt; j++) {
1388f6e5521dSKarl Rupp       for (k=0; k<recvs2[cnt+1]; k++) (*indices)[bprocs[recvs2[cnt+2+k]]][(*numprocs)[bprocs[recvs2[cnt+2+k]]]++] = recvs2[cnt];
138930dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt+1];
139030dcb7c9SBarry Smith     }
139130dcb7c9SBarry Smith   }
139230dcb7c9SBarry Smith   ierr = PetscFree(bprocs);CHKERRQ(ierr);
139307b52d57SBarry Smith   ierr = PetscFree(recvs2);CHKERRQ(ierr);
139430dcb7c9SBarry Smith 
139507b52d57SBarry Smith   /* sort the node indexing by their global numbers */
139607b52d57SBarry Smith   nt = *nproc;
139707b52d57SBarry Smith   for (i=0; i<nt; i++) {
1398854ce69bSBarry Smith     ierr = PetscMalloc1((*numprocs)[i],&tmp);CHKERRQ(ierr);
1399f6e5521dSKarl Rupp     for (j=0; j<(*numprocs)[i]; j++) tmp[j] = lindices[(*indices)[i][j]];
140007b52d57SBarry Smith     ierr = PetscSortIntWithArray((*numprocs)[i],tmp,(*indices)[i]);CHKERRQ(ierr);
140107b52d57SBarry Smith     ierr = PetscFree(tmp);CHKERRQ(ierr);
140207b52d57SBarry Smith   }
140307b52d57SBarry Smith 
140407b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
140530dcb7c9SBarry Smith     nt = *nproc;
140630dcb7c9SBarry Smith     for (i=0; i<nt; i++) {
14077904a332SBarry Smith       ierr = PetscSynchronizedPrintf(comm,"[%d] subdomain %D number of indices %D: ",rank,(*procs)[i],(*numprocs)[i]);CHKERRQ(ierr);
140830dcb7c9SBarry Smith       for (j=0; j<(*numprocs)[i]; j++) {
14097904a332SBarry Smith         ierr = PetscSynchronizedPrintf(comm,"%D ",(*indices)[i][j]);CHKERRQ(ierr);
141030dcb7c9SBarry Smith       }
141130dcb7c9SBarry Smith       ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr);
141230dcb7c9SBarry Smith     }
14130ec8b6e3SBarry Smith     ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr);
141407b52d57SBarry Smith   } /* -----------------------------------  */
141530dcb7c9SBarry Smith 
141630dcb7c9SBarry Smith   /* wait on sends */
141730dcb7c9SBarry Smith   if (nsends2) {
1418785e854fSJed Brown     ierr = PetscMalloc1(nsends2,&send_status);CHKERRQ(ierr);
1419ffc4695bSBarry Smith     ierr = MPI_Waitall(nsends2,send_waits,send_status);CHKERRMPI(ierr);
142030dcb7c9SBarry Smith     ierr = PetscFree(send_status);CHKERRQ(ierr);
142130dcb7c9SBarry Smith   }
142230dcb7c9SBarry Smith 
142330dcb7c9SBarry Smith   ierr = PetscFree(starts3);CHKERRQ(ierr);
142430dcb7c9SBarry Smith   ierr = PetscFree(dest);CHKERRQ(ierr);
142530dcb7c9SBarry Smith   ierr = PetscFree(send_waits);CHKERRQ(ierr);
14263677ff5aSBarry Smith 
1427bc8ff85bSBarry Smith   ierr = PetscFree(nownedsenders);CHKERRQ(ierr);
1428bc8ff85bSBarry Smith   ierr = PetscFree(ownedsenders);CHKERRQ(ierr);
1429bc8ff85bSBarry Smith   ierr = PetscFree(starts);CHKERRQ(ierr);
143030dcb7c9SBarry Smith   ierr = PetscFree(starts2);CHKERRQ(ierr);
143130dcb7c9SBarry Smith   ierr = PetscFree(lens2);CHKERRQ(ierr);
143289d82c54SBarry Smith 
143389d82c54SBarry Smith   ierr = PetscFree(source);CHKERRQ(ierr);
143497f1f81fSBarry Smith   ierr = PetscFree(len);CHKERRQ(ierr);
143589d82c54SBarry Smith   ierr = PetscFree(recvs);CHKERRQ(ierr);
14363a96401aSBarry Smith   ierr = PetscFree(nprocs);CHKERRQ(ierr);
143730dcb7c9SBarry Smith   ierr = PetscFree(sends2);CHKERRQ(ierr);
143824cf384cSBarry Smith 
143924cf384cSBarry Smith   /* put the information about myself as the first entry in the list */
144024cf384cSBarry Smith   first_procs    = (*procs)[0];
144124cf384cSBarry Smith   first_numprocs = (*numprocs)[0];
144224cf384cSBarry Smith   first_indices  = (*indices)[0];
144324cf384cSBarry Smith   for (i=0; i<*nproc; i++) {
144424cf384cSBarry Smith     if ((*procs)[i] == rank) {
144524cf384cSBarry Smith       (*procs)[0]    = (*procs)[i];
144624cf384cSBarry Smith       (*numprocs)[0] = (*numprocs)[i];
144724cf384cSBarry Smith       (*indices)[0]  = (*indices)[i];
144824cf384cSBarry Smith       (*procs)[i]    = first_procs;
144924cf384cSBarry Smith       (*numprocs)[i] = first_numprocs;
145024cf384cSBarry Smith       (*indices)[i]  = first_indices;
145124cf384cSBarry Smith       break;
145224cf384cSBarry Smith     }
145324cf384cSBarry Smith   }
1454268a049cSStefano Zampini 
1455268a049cSStefano Zampini   /* save info for reuse */
1456268a049cSStefano Zampini   mapping->info_nproc = *nproc;
1457268a049cSStefano Zampini   mapping->info_procs = *procs;
1458268a049cSStefano Zampini   mapping->info_numprocs = *numprocs;
1459268a049cSStefano Zampini   mapping->info_indices = *indices;
1460268a049cSStefano Zampini   mapping->info_cached = PETSC_TRUE;
146189d82c54SBarry Smith   PetscFunctionReturn(0);
146289d82c54SBarry Smith }
146389d82c54SBarry Smith 
14646a818285SBarry Smith /*@C
14656a818285SBarry Smith     ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by ISLocalToGlobalMappingGetBlockInfo()
14666a818285SBarry Smith 
14676a818285SBarry Smith     Collective on ISLocalToGlobalMapping
14686a818285SBarry Smith 
14696a818285SBarry Smith     Input Parameters:
14706a818285SBarry Smith .   mapping - the mapping from local to global indexing
14716a818285SBarry Smith 
14726a818285SBarry Smith     Output Parameter:
14736a818285SBarry Smith +   nproc - number of processors that are connected to this one
14746a818285SBarry Smith .   proc - neighboring processors
14756a818285SBarry Smith .   numproc - number of indices for each processor
14766a818285SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
14776a818285SBarry Smith 
14786a818285SBarry Smith     Level: advanced
14796a818285SBarry Smith 
14806a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
14816a818285SBarry Smith           ISLocalToGlobalMappingGetInfo()
14826a818285SBarry Smith @*/
14836a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
14846a818285SBarry Smith {
14856a818285SBarry Smith   PetscErrorCode ierr;
14866a818285SBarry Smith 
14876a818285SBarry Smith   PetscFunctionBegin;
1488cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1489268a049cSStefano Zampini   if (mapping->info_free) {
14906a818285SBarry Smith     ierr = PetscFree(*numprocs);CHKERRQ(ierr);
14916a818285SBarry Smith     if (*indices) {
1492268a049cSStefano Zampini       PetscInt i;
1493268a049cSStefano Zampini 
14946a818285SBarry Smith       ierr = PetscFree((*indices)[0]);CHKERRQ(ierr);
14956a818285SBarry Smith       for (i=1; i<*nproc; i++) {
14966a818285SBarry Smith         ierr = PetscFree((*indices)[i]);CHKERRQ(ierr);
14976a818285SBarry Smith       }
14986a818285SBarry Smith       ierr = PetscFree(*indices);CHKERRQ(ierr);
14996a818285SBarry Smith     }
1500268a049cSStefano Zampini   }
1501268a049cSStefano Zampini   *nproc    = 0;
1502268a049cSStefano Zampini   *procs    = NULL;
1503268a049cSStefano Zampini   *numprocs = NULL;
1504268a049cSStefano Zampini   *indices  = NULL;
15056a818285SBarry Smith   PetscFunctionReturn(0);
15066a818285SBarry Smith }
15076a818285SBarry Smith 
15086a818285SBarry Smith /*@C
15096a818285SBarry Smith     ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each processor and
15106a818285SBarry Smith      each index shared by more than one processor
15116a818285SBarry Smith 
15126a818285SBarry Smith     Collective on ISLocalToGlobalMapping
15136a818285SBarry Smith 
15146a818285SBarry Smith     Input Parameters:
15156a818285SBarry Smith .   mapping - the mapping from local to global indexing
15166a818285SBarry Smith 
15176a818285SBarry Smith     Output Parameter:
15186a818285SBarry Smith +   nproc - number of processors that are connected to this one
15196a818285SBarry Smith .   proc - neighboring processors
15206a818285SBarry Smith .   numproc - number of indices for each subdomain (processor)
15216a818285SBarry Smith -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
15226a818285SBarry Smith 
15236a818285SBarry Smith     Level: advanced
15246a818285SBarry Smith 
15251bd0b88eSStefano Zampini     Notes: The user needs to call ISLocalToGlobalMappingRestoreInfo when the data is no longer needed.
15261bd0b88eSStefano Zampini 
15276a818285SBarry Smith     Fortran Usage:
15286a818285SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
15296a818285SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
15306a818285SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
15316a818285SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
15326a818285SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
15336a818285SBarry Smith 
15346a818285SBarry Smith 
15356a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
15366a818285SBarry Smith           ISLocalToGlobalMappingRestoreInfo()
15376a818285SBarry Smith @*/
15386a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
15396a818285SBarry Smith {
15406a818285SBarry Smith   PetscErrorCode ierr;
1541268a049cSStefano Zampini   PetscInt       **bindices = NULL,*bnumprocs = NULL,bs = mapping->bs,i,j,k;
15426a818285SBarry Smith 
15436a818285SBarry Smith   PetscFunctionBegin;
15446a818285SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1545268a049cSStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockInfo(mapping,nproc,procs,&bnumprocs,&bindices);CHKERRQ(ierr);
1546268a049cSStefano Zampini   if (bs > 1) { /* we need to expand the cached info */
1547732f65e3SBarry Smith     ierr = PetscCalloc1(*nproc,&*indices);CHKERRQ(ierr);
1548268a049cSStefano Zampini     ierr = PetscCalloc1(*nproc,&*numprocs);CHKERRQ(ierr);
15496a818285SBarry Smith     for (i=0; i<*nproc; i++) {
1550268a049cSStefano Zampini       ierr = PetscMalloc1(bs*bnumprocs[i],&(*indices)[i]);CHKERRQ(ierr);
1551268a049cSStefano Zampini       for (j=0; j<bnumprocs[i]; j++) {
15526a818285SBarry Smith         for (k=0; k<bs; k++) {
15536a818285SBarry Smith           (*indices)[i][j*bs+k] = bs*bindices[i][j] + k;
15546a818285SBarry Smith         }
15556a818285SBarry Smith       }
1556268a049cSStefano Zampini       (*numprocs)[i] = bnumprocs[i]*bs;
15576a818285SBarry Smith     }
1558268a049cSStefano Zampini     mapping->info_free = PETSC_TRUE;
1559268a049cSStefano Zampini   } else {
1560268a049cSStefano Zampini     *numprocs = bnumprocs;
1561268a049cSStefano Zampini     *indices  = bindices;
15626a818285SBarry Smith   }
15636a818285SBarry Smith   PetscFunctionReturn(0);
15646a818285SBarry Smith }
15656a818285SBarry Smith 
156607b52d57SBarry Smith /*@C
156707b52d57SBarry Smith     ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by ISLocalToGlobalMappingGetInfo()
156889d82c54SBarry Smith 
156907b52d57SBarry Smith     Collective on ISLocalToGlobalMapping
157007b52d57SBarry Smith 
157107b52d57SBarry Smith     Input Parameters:
157207b52d57SBarry Smith .   mapping - the mapping from local to global indexing
157307b52d57SBarry Smith 
157407b52d57SBarry Smith     Output Parameter:
157507b52d57SBarry Smith +   nproc - number of processors that are connected to this one
157607b52d57SBarry Smith .   proc - neighboring processors
157707b52d57SBarry Smith .   numproc - number of indices for each processor
157807b52d57SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
157907b52d57SBarry Smith 
158007b52d57SBarry Smith     Level: advanced
158107b52d57SBarry Smith 
158207b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
158307b52d57SBarry Smith           ISLocalToGlobalMappingGetInfo()
158407b52d57SBarry Smith @*/
15857087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
158607b52d57SBarry Smith {
15876849ba73SBarry Smith   PetscErrorCode ierr;
158807b52d57SBarry Smith 
158907b52d57SBarry Smith   PetscFunctionBegin;
15906a818285SBarry Smith   ierr = ISLocalToGlobalMappingRestoreBlockInfo(mapping,nproc,procs,numprocs,indices);CHKERRQ(ierr);
159107b52d57SBarry Smith   PetscFunctionReturn(0);
159207b52d57SBarry Smith }
159386994e45SJed Brown 
159486994e45SJed Brown /*@C
15951bd0b88eSStefano Zampini     ISLocalToGlobalMappingGetNodeInfo - Gets the neighbor information for each node
15961bd0b88eSStefano Zampini 
15971bd0b88eSStefano Zampini     Collective on ISLocalToGlobalMapping
15981bd0b88eSStefano Zampini 
15991bd0b88eSStefano Zampini     Input Parameters:
16001bd0b88eSStefano Zampini .   mapping - the mapping from local to global indexing
16011bd0b88eSStefano Zampini 
16021bd0b88eSStefano Zampini     Output Parameter:
16031bd0b88eSStefano Zampini +   nnodes - number of local nodes (same ISLocalToGlobalMappingGetSize())
16041bd0b88eSStefano Zampini .   count - number of neighboring processors per node
16051bd0b88eSStefano Zampini -   indices - indices of processes sharing the node (sorted)
16061bd0b88eSStefano Zampini 
16071bd0b88eSStefano Zampini     Level: advanced
16081bd0b88eSStefano Zampini 
16091bd0b88eSStefano Zampini     Notes: The user needs to call ISLocalToGlobalMappingRestoreInfo when the data is no longer needed.
16101bd0b88eSStefano Zampini 
16111bd0b88eSStefano Zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
16121bd0b88eSStefano Zampini           ISLocalToGlobalMappingGetInfo(), ISLocalToGlobalMappingRestoreNodeInfo()
16131bd0b88eSStefano Zampini @*/
16141bd0b88eSStefano Zampini PetscErrorCode  ISLocalToGlobalMappingGetNodeInfo(ISLocalToGlobalMapping mapping,PetscInt *nnodes,PetscInt *count[],PetscInt **indices[])
16151bd0b88eSStefano Zampini {
16161bd0b88eSStefano Zampini   PetscInt       n;
16171bd0b88eSStefano Zampini   PetscErrorCode ierr;
16181bd0b88eSStefano Zampini 
16191bd0b88eSStefano Zampini   PetscFunctionBegin;
16201bd0b88eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
16211bd0b88eSStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(mapping,&n);CHKERRQ(ierr);
16221bd0b88eSStefano Zampini   if (!mapping->info_nodec) {
16231bd0b88eSStefano Zampini     PetscInt i,m,n_neigh,*neigh,*n_shared,**shared;
16241bd0b88eSStefano Zampini 
1625071fcb05SBarry Smith     ierr = PetscMalloc2(n+1,&mapping->info_nodec,n,&mapping->info_nodei);CHKERRQ(ierr);
16261bd0b88eSStefano Zampini     ierr = ISLocalToGlobalMappingGetInfo(mapping,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
1627071fcb05SBarry Smith     for (i=0;i<n;i++) { mapping->info_nodec[i] = 1;}
1628071fcb05SBarry Smith     m = n;
1629071fcb05SBarry Smith     mapping->info_nodec[n] = 0;
16301bd0b88eSStefano Zampini     for (i=1;i<n_neigh;i++) {
16311bd0b88eSStefano Zampini       PetscInt j;
16321bd0b88eSStefano Zampini 
16331bd0b88eSStefano Zampini       m += n_shared[i];
16341bd0b88eSStefano Zampini       for (j=0;j<n_shared[i];j++) mapping->info_nodec[shared[i][j]] += 1;
16351bd0b88eSStefano Zampini     }
16361bd0b88eSStefano Zampini     if (n) { ierr = PetscMalloc1(m,&mapping->info_nodei[0]);CHKERRQ(ierr); }
16371bd0b88eSStefano Zampini     for (i=1;i<n;i++) mapping->info_nodei[i] = mapping->info_nodei[i-1] + mapping->info_nodec[i-1];
1638580bdb30SBarry Smith     ierr = PetscArrayzero(mapping->info_nodec,n);CHKERRQ(ierr);
16391bd0b88eSStefano Zampini     for (i=0;i<n;i++) { mapping->info_nodec[i] = 1; mapping->info_nodei[i][0] = neigh[0]; }
16401bd0b88eSStefano Zampini     for (i=1;i<n_neigh;i++) {
16411bd0b88eSStefano Zampini       PetscInt j;
16421bd0b88eSStefano Zampini 
16431bd0b88eSStefano Zampini       for (j=0;j<n_shared[i];j++) {
16441bd0b88eSStefano Zampini         PetscInt k = shared[i][j];
16451bd0b88eSStefano Zampini 
16461bd0b88eSStefano Zampini         mapping->info_nodei[k][mapping->info_nodec[k]] = neigh[i];
16471bd0b88eSStefano Zampini         mapping->info_nodec[k] += 1;
16481bd0b88eSStefano Zampini       }
16491bd0b88eSStefano Zampini     }
16501bd0b88eSStefano Zampini     for (i=0;i<n;i++) { ierr = PetscSortRemoveDupsInt(&mapping->info_nodec[i],mapping->info_nodei[i]);CHKERRQ(ierr); }
16511bd0b88eSStefano Zampini     ierr = ISLocalToGlobalMappingRestoreInfo(mapping,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
16521bd0b88eSStefano Zampini   }
16531bd0b88eSStefano Zampini   if (nnodes)  *nnodes  = n;
16541bd0b88eSStefano Zampini   if (count)   *count   = mapping->info_nodec;
16551bd0b88eSStefano Zampini   if (indices) *indices = mapping->info_nodei;
16561bd0b88eSStefano Zampini   PetscFunctionReturn(0);
16571bd0b88eSStefano Zampini }
16581bd0b88eSStefano Zampini 
16591bd0b88eSStefano Zampini /*@C
16601bd0b88eSStefano Zampini     ISLocalToGlobalMappingRestoreNodeInfo - Frees the memory allocated by ISLocalToGlobalMappingGetNodeInfo()
16611bd0b88eSStefano Zampini 
16621bd0b88eSStefano Zampini     Collective on ISLocalToGlobalMapping
16631bd0b88eSStefano Zampini 
16641bd0b88eSStefano Zampini     Input Parameters:
16651bd0b88eSStefano Zampini .   mapping - the mapping from local to global indexing
16661bd0b88eSStefano Zampini 
16671bd0b88eSStefano Zampini     Output Parameter:
16681bd0b88eSStefano Zampini +   nnodes - number of local nodes
16691bd0b88eSStefano Zampini .   count - number of neighboring processors per node
16701bd0b88eSStefano Zampini -   indices - indices of processes sharing the node (sorted)
16711bd0b88eSStefano Zampini 
16721bd0b88eSStefano Zampini     Level: advanced
16731bd0b88eSStefano Zampini 
16741bd0b88eSStefano Zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
16751bd0b88eSStefano Zampini           ISLocalToGlobalMappingGetInfo()
16761bd0b88eSStefano Zampini @*/
16771bd0b88eSStefano Zampini PetscErrorCode  ISLocalToGlobalMappingRestoreNodeInfo(ISLocalToGlobalMapping mapping,PetscInt *nnodes,PetscInt *count[],PetscInt **indices[])
16781bd0b88eSStefano Zampini {
16791bd0b88eSStefano Zampini   PetscFunctionBegin;
16801bd0b88eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
16811bd0b88eSStefano Zampini   if (nnodes)  *nnodes  = 0;
16821bd0b88eSStefano Zampini   if (count)   *count   = NULL;
16831bd0b88eSStefano Zampini   if (indices) *indices = NULL;
16841bd0b88eSStefano Zampini   PetscFunctionReturn(0);
16851bd0b88eSStefano Zampini }
16861bd0b88eSStefano Zampini 
16871bd0b88eSStefano Zampini /*@C
1688107e9a97SBarry Smith    ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped
168986994e45SJed Brown 
169086994e45SJed Brown    Not Collective
169186994e45SJed Brown 
169286994e45SJed Brown    Input Arguments:
169386994e45SJed Brown . ltog - local to global mapping
169486994e45SJed Brown 
169586994e45SJed Brown    Output Arguments:
1696565245c5SBarry Smith . array - array of indices, the length of this array may be obtained with ISLocalToGlobalMappingGetSize()
169786994e45SJed Brown 
169886994e45SJed Brown    Level: advanced
169986994e45SJed Brown 
170095452b02SPatrick Sanan    Notes:
170195452b02SPatrick Sanan     ISLocalToGlobalMappingGetSize() returns the length the this array
1702107e9a97SBarry Smith 
1703107e9a97SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreIndices(), ISLocalToGlobalMappingGetBlockIndices(), ISLocalToGlobalMappingRestoreBlockIndices()
170486994e45SJed Brown @*/
17057087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
170686994e45SJed Brown {
170786994e45SJed Brown   PetscFunctionBegin;
170886994e45SJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
170986994e45SJed Brown   PetscValidPointer(array,2);
171045b6f7e9SBarry Smith   if (ltog->bs == 1) {
171186994e45SJed Brown     *array = ltog->indices;
171245b6f7e9SBarry Smith   } else {
171345b6f7e9SBarry Smith     PetscInt       *jj,k,i,j,n = ltog->n, bs = ltog->bs;
171445b6f7e9SBarry Smith     const PetscInt *ii;
171545b6f7e9SBarry Smith     PetscErrorCode ierr;
171645b6f7e9SBarry Smith 
171745b6f7e9SBarry Smith     ierr = PetscMalloc1(bs*n,&jj);CHKERRQ(ierr);
171845b6f7e9SBarry Smith     *array = jj;
171945b6f7e9SBarry Smith     k    = 0;
172045b6f7e9SBarry Smith     ii   = ltog->indices;
172145b6f7e9SBarry Smith     for (i=0; i<n; i++)
172245b6f7e9SBarry Smith       for (j=0; j<bs; j++)
172345b6f7e9SBarry Smith         jj[k++] = bs*ii[i] + j;
172445b6f7e9SBarry Smith   }
172586994e45SJed Brown   PetscFunctionReturn(0);
172686994e45SJed Brown }
172786994e45SJed Brown 
172886994e45SJed Brown /*@C
1729193a2b41SJulian Andrej    ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with ISLocalToGlobalMappingGetIndices()
173086994e45SJed Brown 
173186994e45SJed Brown    Not Collective
173286994e45SJed Brown 
173386994e45SJed Brown    Input Arguments:
173486994e45SJed Brown + ltog - local to global mapping
173586994e45SJed Brown - array - array of indices
173686994e45SJed Brown 
173786994e45SJed Brown    Level: advanced
173886994e45SJed Brown 
173986994e45SJed Brown .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices()
174086994e45SJed Brown @*/
17417087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
174286994e45SJed Brown {
174386994e45SJed Brown   PetscFunctionBegin;
174486994e45SJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
174586994e45SJed Brown   PetscValidPointer(array,2);
174645b6f7e9SBarry Smith   if (ltog->bs == 1 && *array != ltog->indices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
174745b6f7e9SBarry Smith 
174845b6f7e9SBarry Smith   if (ltog->bs > 1) {
174945b6f7e9SBarry Smith     PetscErrorCode ierr;
175045b6f7e9SBarry Smith     ierr = PetscFree(*(void**)array);CHKERRQ(ierr);
175145b6f7e9SBarry Smith   }
175245b6f7e9SBarry Smith   PetscFunctionReturn(0);
175345b6f7e9SBarry Smith }
175445b6f7e9SBarry Smith 
175545b6f7e9SBarry Smith /*@C
175645b6f7e9SBarry Smith    ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block
175745b6f7e9SBarry Smith 
175845b6f7e9SBarry Smith    Not Collective
175945b6f7e9SBarry Smith 
176045b6f7e9SBarry Smith    Input Arguments:
176145b6f7e9SBarry Smith . ltog - local to global mapping
176245b6f7e9SBarry Smith 
176345b6f7e9SBarry Smith    Output Arguments:
176445b6f7e9SBarry Smith . array - array of indices
176545b6f7e9SBarry Smith 
176645b6f7e9SBarry Smith    Level: advanced
176745b6f7e9SBarry Smith 
176845b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreBlockIndices()
176945b6f7e9SBarry Smith @*/
177045b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
177145b6f7e9SBarry Smith {
177245b6f7e9SBarry Smith   PetscFunctionBegin;
177345b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
177445b6f7e9SBarry Smith   PetscValidPointer(array,2);
177545b6f7e9SBarry Smith   *array = ltog->indices;
177645b6f7e9SBarry Smith   PetscFunctionReturn(0);
177745b6f7e9SBarry Smith }
177845b6f7e9SBarry Smith 
177945b6f7e9SBarry Smith /*@C
178045b6f7e9SBarry Smith    ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with ISLocalToGlobalMappingGetBlockIndices()
178145b6f7e9SBarry Smith 
178245b6f7e9SBarry Smith    Not Collective
178345b6f7e9SBarry Smith 
178445b6f7e9SBarry Smith    Input Arguments:
178545b6f7e9SBarry Smith + ltog - local to global mapping
178645b6f7e9SBarry Smith - array - array of indices
178745b6f7e9SBarry Smith 
178845b6f7e9SBarry Smith    Level: advanced
178945b6f7e9SBarry Smith 
179045b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices()
179145b6f7e9SBarry Smith @*/
179245b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
179345b6f7e9SBarry Smith {
179445b6f7e9SBarry Smith   PetscFunctionBegin;
179545b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
179645b6f7e9SBarry Smith   PetscValidPointer(array,2);
179786994e45SJed Brown   if (*array != ltog->indices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
17980298fd71SBarry Smith   *array = NULL;
179986994e45SJed Brown   PetscFunctionReturn(0);
180086994e45SJed Brown }
1801f7efa3c7SJed Brown 
1802f7efa3c7SJed Brown /*@C
1803f7efa3c7SJed Brown    ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings
1804f7efa3c7SJed Brown 
1805f7efa3c7SJed Brown    Not Collective
1806f7efa3c7SJed Brown 
1807f7efa3c7SJed Brown    Input Arguments:
1808f7efa3c7SJed Brown + comm - communicator for the new mapping, must contain the communicator of every mapping to concatenate
1809f7efa3c7SJed Brown . n - number of mappings to concatenate
1810f7efa3c7SJed Brown - ltogs - local to global mappings
1811f7efa3c7SJed Brown 
1812f7efa3c7SJed Brown    Output Arguments:
1813f7efa3c7SJed Brown . ltogcat - new mapping
1814f7efa3c7SJed Brown 
18159d90f715SBarry Smith    Note: this currently always returns a mapping with block size of 1
18169d90f715SBarry Smith 
18179d90f715SBarry Smith    Developer Note: If all the input mapping have the same block size we could easily handle that as a special case
18189d90f715SBarry Smith 
1819f7efa3c7SJed Brown    Level: advanced
1820f7efa3c7SJed Brown 
1821f7efa3c7SJed Brown .seealso: ISLocalToGlobalMappingCreate()
1822f7efa3c7SJed Brown @*/
1823f7efa3c7SJed Brown PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm,PetscInt n,const ISLocalToGlobalMapping ltogs[],ISLocalToGlobalMapping *ltogcat)
1824f7efa3c7SJed Brown {
1825f7efa3c7SJed Brown   PetscInt       i,cnt,m,*idx;
1826f7efa3c7SJed Brown   PetscErrorCode ierr;
1827f7efa3c7SJed Brown 
1828f7efa3c7SJed Brown   PetscFunctionBegin;
1829f7efa3c7SJed Brown   if (n < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have a non-negative number of mappings, given %D",n);
1830f7efa3c7SJed Brown   if (n > 0) PetscValidPointer(ltogs,3);
1831f7efa3c7SJed Brown   for (i=0; i<n; i++) PetscValidHeaderSpecific(ltogs[i],IS_LTOGM_CLASSID,3);
1832f7efa3c7SJed Brown   PetscValidPointer(ltogcat,4);
1833f7efa3c7SJed Brown   for (cnt=0,i=0; i<n; i++) {
1834f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
1835f7efa3c7SJed Brown     cnt += m;
1836f7efa3c7SJed Brown   }
1837785e854fSJed Brown   ierr = PetscMalloc1(cnt,&idx);CHKERRQ(ierr);
1838f7efa3c7SJed Brown   for (cnt=0,i=0; i<n; i++) {
1839f7efa3c7SJed Brown     const PetscInt *subidx;
1840f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr);
1841f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingGetIndices(ltogs[i],&subidx);CHKERRQ(ierr);
1842580bdb30SBarry Smith     ierr = PetscArraycpy(&idx[cnt],subidx,m);CHKERRQ(ierr);
1843f7efa3c7SJed Brown     ierr = ISLocalToGlobalMappingRestoreIndices(ltogs[i],&subidx);CHKERRQ(ierr);
1844f7efa3c7SJed Brown     cnt += m;
1845f7efa3c7SJed Brown   }
1846f0413b6fSBarry Smith   ierr = ISLocalToGlobalMappingCreate(comm,1,cnt,idx,PETSC_OWN_POINTER,ltogcat);CHKERRQ(ierr);
1847f7efa3c7SJed Brown   PetscFunctionReturn(0);
1848f7efa3c7SJed Brown }
184904a59952SBarry Smith 
1850413f72f0SBarry Smith /*MC
1851413f72f0SBarry Smith       ISLOCALTOGLOBALMAPPINGBASIC - basic implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is
1852413f72f0SBarry Smith                                     used this is good for only small and moderate size problems.
1853413f72f0SBarry Smith 
1854413f72f0SBarry Smith    Options Database Keys:
1855a2b725a8SWilliam Gropp .   -islocaltoglobalmapping_type basic - select this method
1856413f72f0SBarry Smith 
1857413f72f0SBarry Smith    Level: beginner
1858413f72f0SBarry Smith 
1859413f72f0SBarry Smith .seealso:  ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType(), ISLOCALTOGLOBALMAPPINGHASH
1860413f72f0SBarry Smith M*/
1861413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Basic(ISLocalToGlobalMapping ltog)
1862413f72f0SBarry Smith {
1863413f72f0SBarry Smith   PetscFunctionBegin;
1864413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Basic;
1865413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Basic;
1866413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Basic;
1867413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Basic;
1868413f72f0SBarry Smith   PetscFunctionReturn(0);
1869413f72f0SBarry Smith }
1870413f72f0SBarry Smith 
1871413f72f0SBarry Smith /*MC
1872413f72f0SBarry Smith       ISLOCALTOGLOBALMAPPINGHASH - hash implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is
1873ed56e8eaSBarry Smith                                     used this is good for large memory problems.
1874413f72f0SBarry Smith 
1875413f72f0SBarry Smith    Options Database Keys:
1876a2b725a8SWilliam Gropp .   -islocaltoglobalmapping_type hash - select this method
1877413f72f0SBarry Smith 
1878ed56e8eaSBarry Smith 
187995452b02SPatrick Sanan    Notes:
188095452b02SPatrick Sanan     This is selected automatically for large problems if the user does not set the type.
1881ed56e8eaSBarry Smith 
1882413f72f0SBarry Smith    Level: beginner
1883413f72f0SBarry Smith 
1884413f72f0SBarry Smith .seealso:  ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType(), ISLOCALTOGLOBALMAPPINGHASH
1885413f72f0SBarry Smith M*/
1886413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Hash(ISLocalToGlobalMapping ltog)
1887413f72f0SBarry Smith {
1888413f72f0SBarry Smith   PetscFunctionBegin;
1889413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Hash;
1890413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Hash;
1891413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Hash;
1892413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Hash;
1893413f72f0SBarry Smith   PetscFunctionReturn(0);
1894413f72f0SBarry Smith }
1895413f72f0SBarry Smith 
1896413f72f0SBarry Smith 
1897413f72f0SBarry Smith /*@C
1898413f72f0SBarry Smith     ISLocalToGlobalMappingRegister -  Adds a method for applying a global to local mapping with an ISLocalToGlobalMapping
1899413f72f0SBarry Smith 
1900413f72f0SBarry Smith    Not Collective
1901413f72f0SBarry Smith 
1902413f72f0SBarry Smith    Input Parameters:
1903413f72f0SBarry Smith +  sname - name of a new method
1904413f72f0SBarry Smith -  routine_create - routine to create method context
1905413f72f0SBarry Smith 
1906413f72f0SBarry Smith    Notes:
1907ed56e8eaSBarry Smith    ISLocalToGlobalMappingRegister() may be called multiple times to add several user-defined mappings.
1908413f72f0SBarry Smith 
1909413f72f0SBarry Smith    Sample usage:
1910413f72f0SBarry Smith .vb
1911413f72f0SBarry Smith    ISLocalToGlobalMappingRegister("my_mapper",MyCreate);
1912413f72f0SBarry Smith .ve
1913413f72f0SBarry Smith 
1914ed56e8eaSBarry Smith    Then, your mapping can be chosen with the procedural interface via
1915413f72f0SBarry Smith $     ISLocalToGlobalMappingSetType(ltog,"my_mapper")
1916413f72f0SBarry Smith    or at runtime via the option
1917ed56e8eaSBarry Smith $     -islocaltoglobalmapping_type my_mapper
1918413f72f0SBarry Smith 
1919413f72f0SBarry Smith    Level: advanced
1920413f72f0SBarry Smith 
1921413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingRegisterAll(), ISLocalToGlobalMappingRegisterDestroy(), ISLOCALTOGLOBALMAPPINGBASIC, ISLOCALTOGLOBALMAPPINGHASH
1922413f72f0SBarry Smith 
1923413f72f0SBarry Smith @*/
1924413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRegister(const char sname[],PetscErrorCode (*function)(ISLocalToGlobalMapping))
1925413f72f0SBarry Smith {
1926413f72f0SBarry Smith   PetscErrorCode ierr;
1927413f72f0SBarry Smith 
1928413f72f0SBarry Smith   PetscFunctionBegin;
19291d36bdfdSBarry Smith   ierr = ISInitializePackage();CHKERRQ(ierr);
1930413f72f0SBarry Smith   ierr = PetscFunctionListAdd(&ISLocalToGlobalMappingList,sname,function);CHKERRQ(ierr);
1931413f72f0SBarry Smith   PetscFunctionReturn(0);
1932413f72f0SBarry Smith }
1933413f72f0SBarry Smith 
1934413f72f0SBarry Smith /*@C
1935ed56e8eaSBarry Smith    ISLocalToGlobalMappingSetType - Builds ISLocalToGlobalMapping for a particular global to local mapping approach.
1936413f72f0SBarry Smith 
1937413f72f0SBarry Smith    Logically Collective on ISLocalToGlobalMapping
1938413f72f0SBarry Smith 
1939413f72f0SBarry Smith    Input Parameters:
1940413f72f0SBarry Smith +  ltog - the ISLocalToGlobalMapping object
1941413f72f0SBarry Smith -  type - a known method
1942413f72f0SBarry Smith 
1943413f72f0SBarry Smith    Options Database Key:
1944ed56e8eaSBarry Smith .  -islocaltoglobalmapping_type  <method> - Sets the method; use -help for a list
1945413f72f0SBarry Smith     of available methods (for instance, basic or hash)
1946413f72f0SBarry Smith 
1947413f72f0SBarry Smith    Notes:
1948413f72f0SBarry Smith    See "petsc/include/petscis.h" for available methods
1949413f72f0SBarry Smith 
1950413f72f0SBarry Smith   Normally, it is best to use the ISLocalToGlobalMappingSetFromOptions() command and
1951413f72f0SBarry Smith   then set the ISLocalToGlobalMapping type from the options database rather than by using
1952413f72f0SBarry Smith   this routine.
1953413f72f0SBarry Smith 
1954413f72f0SBarry Smith   Level: intermediate
1955413f72f0SBarry Smith 
1956413f72f0SBarry Smith   Developer Note: ISLocalToGlobalMappingRegister() is used to add new types to ISLocalToGlobalMappingList from which they
1957413f72f0SBarry Smith   are accessed by ISLocalToGlobalMappingSetType().
1958413f72f0SBarry Smith 
1959413f72f0SBarry Smith .seealso: PCSetType(), ISLocalToGlobalMappingType, ISLocalToGlobalMappingRegister(), ISLocalToGlobalMappingCreate()
1960413f72f0SBarry Smith 
1961413f72f0SBarry Smith @*/
1962413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType type)
1963413f72f0SBarry Smith {
1964413f72f0SBarry Smith   PetscErrorCode ierr,(*r)(ISLocalToGlobalMapping);
1965413f72f0SBarry Smith   PetscBool      match;
1966413f72f0SBarry Smith 
1967413f72f0SBarry Smith   PetscFunctionBegin;
1968413f72f0SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
1969413f72f0SBarry Smith   PetscValidCharPointer(type,2);
1970413f72f0SBarry Smith 
1971413f72f0SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)ltog,type,&match);CHKERRQ(ierr);
1972413f72f0SBarry Smith   if (match) PetscFunctionReturn(0);
1973413f72f0SBarry Smith 
1974413f72f0SBarry Smith   ierr =  PetscFunctionListFind(ISLocalToGlobalMappingList,type,&r);CHKERRQ(ierr);
1975413f72f0SBarry Smith   if (!r) SETERRQ1(PetscObjectComm((PetscObject)ltog),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested ISLocalToGlobalMapping type %s",type);
1976413f72f0SBarry Smith   /* Destroy the previous private LTOG context */
1977413f72f0SBarry Smith   if (ltog->ops->destroy) {
1978413f72f0SBarry Smith     ierr              = (*ltog->ops->destroy)(ltog);CHKERRQ(ierr);
1979413f72f0SBarry Smith     ltog->ops->destroy = NULL;
1980413f72f0SBarry Smith   }
1981413f72f0SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)ltog,type);CHKERRQ(ierr);
1982413f72f0SBarry Smith   ierr = (*r)(ltog);CHKERRQ(ierr);
1983413f72f0SBarry Smith   PetscFunctionReturn(0);
1984413f72f0SBarry Smith }
1985413f72f0SBarry Smith 
1986413f72f0SBarry Smith PetscBool ISLocalToGlobalMappingRegisterAllCalled = PETSC_FALSE;
1987413f72f0SBarry Smith 
1988413f72f0SBarry Smith /*@C
1989413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAll - Registers all of the local to global mapping components in the IS package.
1990413f72f0SBarry Smith 
1991413f72f0SBarry Smith   Not Collective
1992413f72f0SBarry Smith 
1993413f72f0SBarry Smith   Level: advanced
1994413f72f0SBarry Smith 
1995413f72f0SBarry Smith .seealso:  ISRegister(),  ISLocalToGlobalRegister()
1996413f72f0SBarry Smith @*/
1997413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRegisterAll(void)
1998413f72f0SBarry Smith {
1999413f72f0SBarry Smith   PetscErrorCode ierr;
2000413f72f0SBarry Smith 
2001413f72f0SBarry Smith   PetscFunctionBegin;
2002413f72f0SBarry Smith   if (ISLocalToGlobalMappingRegisterAllCalled) PetscFunctionReturn(0);
2003413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAllCalled = PETSC_TRUE;
2004413f72f0SBarry Smith 
2005413f72f0SBarry Smith   ierr = ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGBASIC, ISLocalToGlobalMappingCreate_Basic);CHKERRQ(ierr);
2006413f72f0SBarry Smith   ierr = ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGHASH, ISLocalToGlobalMappingCreate_Hash);CHKERRQ(ierr);
2007413f72f0SBarry Smith   PetscFunctionReturn(0);
2008413f72f0SBarry Smith }
2009