xref: /petsc/src/vec/is/utils/isltog.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
41db781477SPatrick Sanan .seealso: `ISRestorePointRange()`, `ISGetPointSubrange()`, `ISGetIndices()`, `ISCreateStride()`
426528b96dSMatthew G. Knepley @*/
439371c9d4SSatish Balay PetscErrorCode ISGetPointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points) {
449305a4c7SMatthew G. Knepley   PetscInt  numCells, step = 1;
459305a4c7SMatthew G. Knepley   PetscBool isStride;
469305a4c7SMatthew G. Knepley 
479305a4c7SMatthew G. Knepley   PetscFunctionBeginHot;
489305a4c7SMatthew G. Knepley   *pStart = 0;
499305a4c7SMatthew G. Knepley   *points = NULL;
509566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(pointIS, &numCells));
519566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pointIS, ISSTRIDE, &isStride));
529566063dSJacob Faibussowitsch   if (isStride) PetscCall(ISStrideGetInfo(pointIS, pStart, &step));
539305a4c7SMatthew G. Knepley   *pEnd = *pStart + numCells;
549566063dSJacob Faibussowitsch   if (!isStride || step != 1) PetscCall(ISGetIndices(pointIS, points));
559305a4c7SMatthew G. Knepley   PetscFunctionReturn(0);
569305a4c7SMatthew G. Knepley }
579305a4c7SMatthew G. Knepley 
586528b96dSMatthew G. Knepley /*@C
596528b96dSMatthew G. Knepley   ISRestorePointRange - Destroys the traversal description
606528b96dSMatthew G. Knepley 
616528b96dSMatthew G. Knepley   Not collective
626528b96dSMatthew G. Knepley 
636528b96dSMatthew G. Knepley   Input Parameters:
646528b96dSMatthew G. Knepley + pointIS - The IS object
656528b96dSMatthew G. Knepley . pStart  - The first index, from ISGetPointRange()
666528b96dSMatthew G. Knepley . pEnd    - One past the last index, from ISGetPointRange()
676528b96dSMatthew G. Knepley - points  - The indices, from ISGetPointRange()
686528b96dSMatthew G. Knepley 
696528b96dSMatthew G. Knepley   Notes:
706528b96dSMatthew 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
716528b96dSMatthew G. Knepley $ ISGetPointRange(is, &pStart, &pEnd, &points);
726528b96dSMatthew G. Knepley $ for (p = pStart; p < pEnd; ++p) {
736528b96dSMatthew G. Knepley $   const PetscInt point = points ? points[p] : p;
746528b96dSMatthew G. Knepley $ }
756528b96dSMatthew G. Knepley $ ISRestorePointRange(is, &pstart, &pEnd, &points);
766528b96dSMatthew G. Knepley 
776528b96dSMatthew G. Knepley   Level: intermediate
786528b96dSMatthew G. Knepley 
79db781477SPatrick Sanan .seealso: `ISGetPointRange()`, `ISGetPointSubrange()`, `ISGetIndices()`, `ISCreateStride()`
806528b96dSMatthew G. Knepley @*/
819371c9d4SSatish Balay PetscErrorCode ISRestorePointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points) {
829305a4c7SMatthew G. Knepley   PetscInt  step = 1;
839305a4c7SMatthew G. Knepley   PetscBool isStride;
849305a4c7SMatthew G. Knepley 
859305a4c7SMatthew G. Knepley   PetscFunctionBeginHot;
869566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pointIS, ISSTRIDE, &isStride));
879566063dSJacob Faibussowitsch   if (isStride) PetscCall(ISStrideGetInfo(pointIS, pStart, &step));
889566063dSJacob Faibussowitsch   if (!isStride || step != 1) PetscCall(ISGetIndices(pointIS, points));
899305a4c7SMatthew G. Knepley   PetscFunctionReturn(0);
909305a4c7SMatthew G. Knepley }
919305a4c7SMatthew G. Knepley 
926528b96dSMatthew G. Knepley /*@C
936528b96dSMatthew G. Knepley   ISGetPointSubrange - Configures the input IS to be a subrange for the traversal information given
946528b96dSMatthew G. Knepley 
956528b96dSMatthew G. Knepley   Not collective
966528b96dSMatthew G. Knepley 
976528b96dSMatthew G. Knepley   Input Parameters:
986528b96dSMatthew G. Knepley + subpointIS - The IS object to be configured
996528b96dSMatthew G. Knepley . pStar   t  - The first index of the subrange
1006528b96dSMatthew G. Knepley . pEnd       - One past the last index for the subrange
1016528b96dSMatthew G. Knepley - points     - The indices for the entire range, from ISGetPointRange()
1026528b96dSMatthew G. Knepley 
1036528b96dSMatthew G. Knepley   Output Parameters:
1046528b96dSMatthew G. Knepley . subpointIS - The IS object now configured to be a subrange
1056528b96dSMatthew G. Knepley 
1066528b96dSMatthew G. Knepley   Notes:
1076528b96dSMatthew G. Knepley   The input IS will now respond properly to calls to ISGetPointRange() and return the subrange.
1086528b96dSMatthew G. Knepley 
1096528b96dSMatthew G. Knepley   Level: intermediate
1106528b96dSMatthew G. Knepley 
111db781477SPatrick Sanan .seealso: `ISGetPointRange()`, `ISRestorePointRange()`, `ISGetIndices()`, `ISCreateStride()`
1126528b96dSMatthew G. Knepley @*/
1139371c9d4SSatish Balay PetscErrorCode ISGetPointSubrange(IS subpointIS, PetscInt pStart, PetscInt pEnd, const PetscInt *points) {
1149305a4c7SMatthew G. Knepley   PetscFunctionBeginHot;
1159305a4c7SMatthew G. Knepley   if (points) {
1169566063dSJacob Faibussowitsch     PetscCall(ISSetType(subpointIS, ISGENERAL));
1179566063dSJacob Faibussowitsch     PetscCall(ISGeneralSetIndices(subpointIS, pEnd - pStart, &points[pStart], PETSC_USE_POINTER));
1189305a4c7SMatthew G. Knepley   } else {
1199566063dSJacob Faibussowitsch     PetscCall(ISSetType(subpointIS, ISSTRIDE));
1209566063dSJacob Faibussowitsch     PetscCall(ISStrideSetStride(subpointIS, pEnd - pStart, pStart, 1));
1219305a4c7SMatthew G. Knepley   }
1229305a4c7SMatthew G. Knepley   PetscFunctionReturn(0);
1239305a4c7SMatthew G. Knepley }
1249305a4c7SMatthew G. Knepley 
125413f72f0SBarry Smith /* -----------------------------------------------------------------------------------------*/
126413f72f0SBarry Smith 
127413f72f0SBarry Smith /*
128413f72f0SBarry Smith     Creates the global mapping information in the ISLocalToGlobalMapping structure
129413f72f0SBarry Smith 
130413f72f0SBarry Smith     If the user has not selected how to handle the global to local mapping then use HASH for "large" problems
131413f72f0SBarry Smith */
1329371c9d4SSatish Balay static PetscErrorCode ISGlobalToLocalMappingSetUp(ISLocalToGlobalMapping mapping) {
133413f72f0SBarry Smith   PetscInt i, *idx = mapping->indices, n = mapping->n, end, start;
134413f72f0SBarry Smith 
135413f72f0SBarry Smith   PetscFunctionBegin;
136413f72f0SBarry Smith   if (mapping->data) PetscFunctionReturn(0);
137413f72f0SBarry Smith   end   = 0;
138413f72f0SBarry Smith   start = PETSC_MAX_INT;
139413f72f0SBarry Smith 
140413f72f0SBarry Smith   for (i = 0; i < n; i++) {
141413f72f0SBarry Smith     if (idx[i] < 0) continue;
142413f72f0SBarry Smith     if (idx[i] < start) start = idx[i];
143413f72f0SBarry Smith     if (idx[i] > end) end = idx[i];
144413f72f0SBarry Smith   }
1459371c9d4SSatish Balay   if (start > end) {
1469371c9d4SSatish Balay     start = 0;
1479371c9d4SSatish Balay     end   = -1;
1489371c9d4SSatish Balay   }
149413f72f0SBarry Smith   mapping->globalstart = start;
150413f72f0SBarry Smith   mapping->globalend   = end;
151413f72f0SBarry Smith   if (!((PetscObject)mapping)->type_name) {
152413f72f0SBarry Smith     if ((end - start) > PetscMax(4 * n, 1000000)) {
1539566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingSetType(mapping, ISLOCALTOGLOBALMAPPINGHASH));
154413f72f0SBarry Smith     } else {
1559566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingSetType(mapping, ISLOCALTOGLOBALMAPPINGBASIC));
156413f72f0SBarry Smith     }
157413f72f0SBarry Smith   }
158dbbe0bcdSBarry Smith   PetscTryTypeMethod(mapping, globaltolocalmappingsetup);
159413f72f0SBarry Smith   PetscFunctionReturn(0);
160413f72f0SBarry Smith }
161413f72f0SBarry Smith 
1629371c9d4SSatish Balay static PetscErrorCode ISGlobalToLocalMappingSetUp_Basic(ISLocalToGlobalMapping mapping) {
163413f72f0SBarry Smith   PetscInt                      i, *idx = mapping->indices, n = mapping->n, end, start, *globals;
164413f72f0SBarry Smith   ISLocalToGlobalMapping_Basic *map;
165413f72f0SBarry Smith 
166413f72f0SBarry Smith   PetscFunctionBegin;
167413f72f0SBarry Smith   start = mapping->globalstart;
168413f72f0SBarry Smith   end   = mapping->globalend;
1699566063dSJacob Faibussowitsch   PetscCall(PetscNew(&map));
1709566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(end - start + 2, &globals));
171413f72f0SBarry Smith   map->globals = globals;
172413f72f0SBarry Smith   for (i = 0; i < end - start + 1; i++) globals[i] = -1;
173413f72f0SBarry Smith   for (i = 0; i < n; i++) {
174413f72f0SBarry Smith     if (idx[i] < 0) continue;
175413f72f0SBarry Smith     globals[idx[i] - start] = i;
176413f72f0SBarry Smith   }
177413f72f0SBarry Smith   mapping->data = (void *)map;
1789566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)mapping, (end - start + 1) * sizeof(PetscInt)));
179413f72f0SBarry Smith   PetscFunctionReturn(0);
180413f72f0SBarry Smith }
181413f72f0SBarry Smith 
1829371c9d4SSatish Balay static PetscErrorCode ISGlobalToLocalMappingSetUp_Hash(ISLocalToGlobalMapping mapping) {
183413f72f0SBarry Smith   PetscInt                     i, *idx = mapping->indices, n = mapping->n;
184413f72f0SBarry Smith   ISLocalToGlobalMapping_Hash *map;
185413f72f0SBarry Smith 
186413f72f0SBarry Smith   PetscFunctionBegin;
1879566063dSJacob Faibussowitsch   PetscCall(PetscNew(&map));
1889566063dSJacob Faibussowitsch   PetscCall(PetscHMapICreate(&map->globalht));
189413f72f0SBarry Smith   for (i = 0; i < n; i++) {
190413f72f0SBarry Smith     if (idx[i] < 0) continue;
1919566063dSJacob Faibussowitsch     PetscCall(PetscHMapISet(map->globalht, idx[i], i));
192413f72f0SBarry Smith   }
193413f72f0SBarry Smith   mapping->data = (void *)map;
1949566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)mapping, 2 * n * sizeof(PetscInt)));
195413f72f0SBarry Smith   PetscFunctionReturn(0);
196413f72f0SBarry Smith }
197413f72f0SBarry Smith 
1989371c9d4SSatish Balay static PetscErrorCode ISLocalToGlobalMappingDestroy_Basic(ISLocalToGlobalMapping mapping) {
199413f72f0SBarry Smith   ISLocalToGlobalMapping_Basic *map = (ISLocalToGlobalMapping_Basic *)mapping->data;
200413f72f0SBarry Smith 
201413f72f0SBarry Smith   PetscFunctionBegin;
202413f72f0SBarry Smith   if (!map) PetscFunctionReturn(0);
2039566063dSJacob Faibussowitsch   PetscCall(PetscFree(map->globals));
2049566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping->data));
205413f72f0SBarry Smith   PetscFunctionReturn(0);
206413f72f0SBarry Smith }
207413f72f0SBarry Smith 
2089371c9d4SSatish Balay static PetscErrorCode ISLocalToGlobalMappingDestroy_Hash(ISLocalToGlobalMapping mapping) {
209413f72f0SBarry Smith   ISLocalToGlobalMapping_Hash *map = (ISLocalToGlobalMapping_Hash *)mapping->data;
210413f72f0SBarry Smith 
211413f72f0SBarry Smith   PetscFunctionBegin;
212413f72f0SBarry Smith   if (!map) PetscFunctionReturn(0);
2139566063dSJacob Faibussowitsch   PetscCall(PetscHMapIDestroy(&map->globalht));
2149566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping->data));
215413f72f0SBarry Smith   PetscFunctionReturn(0);
216413f72f0SBarry Smith }
217413f72f0SBarry Smith 
218413f72f0SBarry Smith #define GTOLTYPE _Basic
219413f72f0SBarry Smith #define GTOLNAME _Basic
220541bf97eSAdrian Croucher #define GTOLBS   mapping->bs
2219371c9d4SSatish Balay #define GTOL(g, local) \
2229371c9d4SSatish Balay   do { \
223413f72f0SBarry Smith     local = map->globals[g / bs - start]; \
2240040bde1SJunchao Zhang     if (local >= 0) local = bs * local + (g % bs); \
225413f72f0SBarry Smith   } while (0)
226541bf97eSAdrian Croucher 
227413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
228413f72f0SBarry Smith 
229413f72f0SBarry Smith #define GTOLTYPE _Basic
230413f72f0SBarry Smith #define GTOLNAME Block_Basic
231541bf97eSAdrian Croucher #define GTOLBS   1
2329371c9d4SSatish Balay #define GTOL(g, local) \
2339371c9d4SSatish Balay   do { local = map->globals[g - start]; } while (0)
234413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
235413f72f0SBarry Smith 
236413f72f0SBarry Smith #define GTOLTYPE _Hash
237413f72f0SBarry Smith #define GTOLNAME _Hash
238541bf97eSAdrian Croucher #define GTOLBS   mapping->bs
2399371c9d4SSatish Balay #define GTOL(g, local) \
2409371c9d4SSatish Balay   do { \
241e8f14785SLisandro Dalcin     (void)PetscHMapIGet(map->globalht, g / bs, &local); \
2420040bde1SJunchao Zhang     if (local >= 0) local = bs * local + (g % bs); \
243413f72f0SBarry Smith   } while (0)
244413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
245413f72f0SBarry Smith 
246413f72f0SBarry Smith #define GTOLTYPE _Hash
247413f72f0SBarry Smith #define GTOLNAME Block_Hash
248541bf97eSAdrian Croucher #define GTOLBS   1
2499371c9d4SSatish Balay #define GTOL(g, local) \
2509371c9d4SSatish Balay   do { (void)PetscHMapIGet(map->globalht, g, &local); } while (0)
251413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
252413f72f0SBarry Smith 
2536658fb44Sstefano_zampini /*@
2546658fb44Sstefano_zampini     ISLocalToGlobalMappingDuplicate - Duplicates the local to global mapping object
2556658fb44Sstefano_zampini 
2566658fb44Sstefano_zampini     Not Collective
2576658fb44Sstefano_zampini 
2586658fb44Sstefano_zampini     Input Parameter:
2596658fb44Sstefano_zampini .   ltog - local to global mapping
2606658fb44Sstefano_zampini 
2616658fb44Sstefano_zampini     Output Parameter:
2626658fb44Sstefano_zampini .   nltog - the duplicated local to global mapping
2636658fb44Sstefano_zampini 
2646658fb44Sstefano_zampini     Level: advanced
2656658fb44Sstefano_zampini 
266db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`
2676658fb44Sstefano_zampini @*/
2689371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingDuplicate(ISLocalToGlobalMapping ltog, ISLocalToGlobalMapping *nltog) {
269a0d79125SStefano Zampini   ISLocalToGlobalMappingType l2gtype;
2706658fb44Sstefano_zampini 
2716658fb44Sstefano_zampini   PetscFunctionBegin;
2726658fb44Sstefano_zampini   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
2739566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)ltog), ltog->bs, ltog->n, ltog->indices, PETSC_COPY_VALUES, nltog));
2749566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetType(ltog, &l2gtype));
2759566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingSetType(*nltog, l2gtype));
2766658fb44Sstefano_zampini   PetscFunctionReturn(0);
2776658fb44Sstefano_zampini }
2786658fb44Sstefano_zampini 
279565245c5SBarry Smith /*@
280107e9a97SBarry Smith     ISLocalToGlobalMappingGetSize - Gets the local size of a local to global mapping
2813b9aefa3SBarry Smith 
2823b9aefa3SBarry Smith     Not Collective
2833b9aefa3SBarry Smith 
2843b9aefa3SBarry Smith     Input Parameter:
2853b9aefa3SBarry Smith .   ltog - local to global mapping
2863b9aefa3SBarry Smith 
2873b9aefa3SBarry Smith     Output Parameter:
288107e9a97SBarry Smith .   n - the number of entries in the local mapping, ISLocalToGlobalMappingGetIndices() returns an array of this length
2893b9aefa3SBarry Smith 
2903b9aefa3SBarry Smith     Level: advanced
2913b9aefa3SBarry Smith 
292db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`
2933b9aefa3SBarry Smith @*/
2949371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping mapping, PetscInt *n) {
2953b9aefa3SBarry Smith   PetscFunctionBegin;
2960700a824SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
2974482741eSBarry Smith   PetscValidIntPointer(n, 2);
298107e9a97SBarry Smith   *n = mapping->bs * mapping->n;
2993b9aefa3SBarry Smith   PetscFunctionReturn(0);
3003b9aefa3SBarry Smith }
3013b9aefa3SBarry Smith 
3025a5d4f66SBarry Smith /*@C
303fe2efc57SMark    ISLocalToGlobalMappingViewFromOptions - View from Options
304fe2efc57SMark 
305fe2efc57SMark    Collective on ISLocalToGlobalMapping
306fe2efc57SMark 
307fe2efc57SMark    Input Parameters:
308fe2efc57SMark +  A - the local to global mapping object
309736c3998SJose E. Roman .  obj - Optional object
310736c3998SJose E. Roman -  name - command line option
311fe2efc57SMark 
312fe2efc57SMark    Level: intermediate
313db781477SPatrick Sanan .seealso: `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingView`, `PetscObjectViewFromOptions()`, `ISLocalToGlobalMappingCreate()`
314fe2efc57SMark @*/
3159371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingViewFromOptions(ISLocalToGlobalMapping A, PetscObject obj, const char name[]) {
316fe2efc57SMark   PetscFunctionBegin;
317fe2efc57SMark   PetscValidHeaderSpecific(A, IS_LTOGM_CLASSID, 1);
3189566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
319fe2efc57SMark   PetscFunctionReturn(0);
320fe2efc57SMark }
321fe2efc57SMark 
322fe2efc57SMark /*@C
3235a5d4f66SBarry Smith     ISLocalToGlobalMappingView - View a local to global mapping
3245a5d4f66SBarry Smith 
325b9cd556bSLois Curfman McInnes     Not Collective
326b9cd556bSLois Curfman McInnes 
3275a5d4f66SBarry Smith     Input Parameters:
3283b9aefa3SBarry Smith +   ltog - local to global mapping
3293b9aefa3SBarry Smith -   viewer - viewer
3305a5d4f66SBarry Smith 
331a997ad1aSLois Curfman McInnes     Level: advanced
332a997ad1aSLois Curfman McInnes 
333db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`
3345a5d4f66SBarry Smith @*/
3359371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingView(ISLocalToGlobalMapping mapping, PetscViewer viewer) {
33632dcc486SBarry Smith   PetscInt    i;
33732dcc486SBarry Smith   PetscMPIInt rank;
338ace3abfcSBarry Smith   PetscBool   iascii;
3395a5d4f66SBarry Smith 
3405a5d4f66SBarry Smith   PetscFunctionBegin;
3410700a824SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
342*48a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mapping), &viewer));
3430700a824SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
3445a5d4f66SBarry Smith 
3459566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mapping), &rank));
3469566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
34732077d6dSBarry Smith   if (iascii) {
3489566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)mapping, viewer));
3499566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
350*48a46eb9SPierre Jolivet     for (i = 0; i < mapping->n; i++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %" PetscInt_FMT " %" PetscInt_FMT "\n", rank, i, mapping->indices[i]));
3519566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
3529566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
3531575c14dSBarry Smith   }
3545a5d4f66SBarry Smith   PetscFunctionReturn(0);
3555a5d4f66SBarry Smith }
3565a5d4f66SBarry Smith 
3571f428162SBarry Smith /*@
3582bdab257SBarry Smith     ISLocalToGlobalMappingCreateIS - Creates a mapping between a local (0 to n)
3592bdab257SBarry Smith     ordering and a global parallel ordering.
3602bdab257SBarry Smith 
3610f5bd95cSBarry Smith     Not collective
362b9cd556bSLois Curfman McInnes 
363a997ad1aSLois Curfman McInnes     Input Parameter:
3648c03b21aSDmitry Karpeev .   is - index set containing the global numbers for each local number
3652bdab257SBarry Smith 
366a997ad1aSLois Curfman McInnes     Output Parameter:
3672bdab257SBarry Smith .   mapping - new mapping data structure
3682bdab257SBarry Smith 
36995452b02SPatrick Sanan     Notes:
37095452b02SPatrick Sanan     the block size of the IS determines the block size of the mapping
371a997ad1aSLois Curfman McInnes     Level: advanced
372a997ad1aSLois Curfman McInnes 
373db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetFromOptions()`
3742bdab257SBarry Smith @*/
3759371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingCreateIS(IS is, ISLocalToGlobalMapping *mapping) {
3763bbf0e92SBarry Smith   PetscInt        n, bs;
3775d0c19d7SBarry Smith   const PetscInt *indices;
3782bdab257SBarry Smith   MPI_Comm        comm;
3793bbf0e92SBarry Smith   PetscBool       isblock;
3803a40ed3dSBarry Smith 
3813a40ed3dSBarry Smith   PetscFunctionBegin;
3820700a824SBarry Smith   PetscValidHeaderSpecific(is, IS_CLASSID, 1);
3834482741eSBarry Smith   PetscValidPointer(mapping, 2);
3842bdab257SBarry Smith 
3859566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)is, &comm));
3869566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is, &n));
3879566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)is, ISBLOCK, &isblock));
3886006e8d2SBarry Smith   if (!isblock) {
3899566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(is, &indices));
3909566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(comm, 1, n, indices, PETSC_COPY_VALUES, mapping));
3919566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(is, &indices));
3926006e8d2SBarry Smith   } else {
3939566063dSJacob Faibussowitsch     PetscCall(ISGetBlockSize(is, &bs));
3949566063dSJacob Faibussowitsch     PetscCall(ISBlockGetIndices(is, &indices));
3959566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(comm, bs, n / bs, indices, PETSC_COPY_VALUES, mapping));
3969566063dSJacob Faibussowitsch     PetscCall(ISBlockRestoreIndices(is, &indices));
3976006e8d2SBarry Smith   }
3983a40ed3dSBarry Smith   PetscFunctionReturn(0);
3992bdab257SBarry Smith }
4005a5d4f66SBarry Smith 
401a4d96a55SJed Brown /*@C
402a4d96a55SJed Brown     ISLocalToGlobalMappingCreateSF - Creates a mapping between a local (0 to n)
403a4d96a55SJed Brown     ordering and a global parallel ordering.
404a4d96a55SJed Brown 
405a4d96a55SJed Brown     Collective
406a4d96a55SJed Brown 
407d8d19677SJose E. Roman     Input Parameters:
408a4d96a55SJed Brown +   sf - star forest mapping contiguous local indices to (rank, offset)
4099a535bafSVaclav Hapla -   start - first global index on this process, or PETSC_DECIDE to compute contiguous global numbering automatically
410a4d96a55SJed Brown 
411a4d96a55SJed Brown     Output Parameter:
412a4d96a55SJed Brown .   mapping - new mapping data structure
413a4d96a55SJed Brown 
414a4d96a55SJed Brown     Level: advanced
415a4d96a55SJed Brown 
4169a535bafSVaclav Hapla     Notes:
4179a535bafSVaclav Hapla     If any processor calls this with start = PETSC_DECIDE then all processors must, otherwise the program will hang.
4189a535bafSVaclav Hapla 
419db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingSetFromOptions()`
420a4d96a55SJed Brown @*/
4219371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF sf, PetscInt start, ISLocalToGlobalMapping *mapping) {
422a4d96a55SJed Brown   PetscInt i, maxlocal, nroots, nleaves, *globals, *ltog;
423a4d96a55SJed Brown   MPI_Comm comm;
424a4d96a55SJed Brown 
425a4d96a55SJed Brown   PetscFunctionBegin;
426a4d96a55SJed Brown   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
427a4d96a55SJed Brown   PetscValidPointer(mapping, 3);
4289566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
42941f4c31fSVaclav Hapla   PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL));
4309a535bafSVaclav Hapla   if (start == PETSC_DECIDE) {
4319a535bafSVaclav Hapla     start = 0;
4329566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Exscan(&nroots, &start, 1, MPIU_INT, MPI_SUM, comm));
43341f4c31fSVaclav Hapla   } else PetscCheck(start >= 0, comm, PETSC_ERR_ARG_OUTOFRANGE, "start must be nonnegative or PETSC_DECIDE");
43441f4c31fSVaclav Hapla   PetscCall(PetscSFGetLeafRange(sf, NULL, &maxlocal));
43541f4c31fSVaclav Hapla   ++maxlocal;
4369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nroots, &globals));
4379566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxlocal, &ltog));
438a4d96a55SJed Brown   for (i = 0; i < nroots; i++) globals[i] = start + i;
439a4d96a55SJed Brown   for (i = 0; i < maxlocal; i++) ltog[i] = -1;
4409566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sf, MPIU_INT, globals, ltog, MPI_REPLACE));
4419566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sf, MPIU_INT, globals, ltog, MPI_REPLACE));
4429566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(comm, 1, maxlocal, ltog, PETSC_OWN_POINTER, mapping));
4439566063dSJacob Faibussowitsch   PetscCall(PetscFree(globals));
444a4d96a55SJed Brown   PetscFunctionReturn(0);
445a4d96a55SJed Brown }
446b46b645bSBarry Smith 
44763fa5c83Sstefano_zampini /*@
44863fa5c83Sstefano_zampini     ISLocalToGlobalMappingSetBlockSize - Sets the blocksize of the mapping
44963fa5c83Sstefano_zampini 
45063fa5c83Sstefano_zampini     Not collective
45163fa5c83Sstefano_zampini 
45263fa5c83Sstefano_zampini     Input Parameters:
453a2b725a8SWilliam Gropp +   mapping - mapping data structure
454a2b725a8SWilliam Gropp -   bs - the blocksize
45563fa5c83Sstefano_zampini 
45663fa5c83Sstefano_zampini     Level: advanced
45763fa5c83Sstefano_zampini 
458db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`
45963fa5c83Sstefano_zampini @*/
4609371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingSetBlockSize(ISLocalToGlobalMapping mapping, PetscInt bs) {
461a59f3c4dSstefano_zampini   PetscInt       *nid;
462a59f3c4dSstefano_zampini   const PetscInt *oid;
463a59f3c4dSstefano_zampini   PetscInt        i, cn, on, obs, nn;
46463fa5c83Sstefano_zampini 
46563fa5c83Sstefano_zampini   PetscFunctionBegin;
46663fa5c83Sstefano_zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
46708401ef6SPierre Jolivet   PetscCheck(bs >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid block size %" PetscInt_FMT, bs);
46863fa5c83Sstefano_zampini   if (bs == mapping->bs) PetscFunctionReturn(0);
46963fa5c83Sstefano_zampini   on  = mapping->n;
47063fa5c83Sstefano_zampini   obs = mapping->bs;
47163fa5c83Sstefano_zampini   oid = mapping->indices;
47263fa5c83Sstefano_zampini   nn  = (on * obs) / bs;
47308401ef6SPierre Jolivet   PetscCheck((on * obs) % bs == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Block size %" PetscInt_FMT " is inconsistent with block size %" PetscInt_FMT " and number of block indices %" PetscInt_FMT, bs, obs, on);
474a59f3c4dSstefano_zampini 
4759566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nn, &nid));
4769566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetIndices(mapping, &oid));
477a59f3c4dSstefano_zampini   for (i = 0; i < nn; i++) {
478a59f3c4dSstefano_zampini     PetscInt j;
479a59f3c4dSstefano_zampini     for (j = 0, cn = 0; j < bs - 1; j++) {
4809371c9d4SSatish Balay       if (oid[i * bs + j] < 0) {
4819371c9d4SSatish Balay         cn++;
4829371c9d4SSatish Balay         continue;
4839371c9d4SSatish Balay       }
48408401ef6SPierre Jolivet       PetscCheck(oid[i * bs + j] == oid[i * bs + j + 1] - 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Block sizes %" PetscInt_FMT " and %" PetscInt_FMT " are incompatible with the block indices: non consecutive indices %" PetscInt_FMT " %" PetscInt_FMT, bs, obs, oid[i * bs + j], oid[i * bs + j + 1]);
485a59f3c4dSstefano_zampini     }
486a59f3c4dSstefano_zampini     if (oid[i * bs + j] < 0) cn++;
4878b7cb0e6Sstefano_zampini     if (cn) {
48808401ef6SPierre Jolivet       PetscCheck(cn == bs, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Block sizes %" PetscInt_FMT " and %" PetscInt_FMT " are incompatible with the block indices: invalid number of negative entries in block %" PetscInt_FMT, bs, obs, cn);
489a59f3c4dSstefano_zampini       nid[i] = -1;
4908b7cb0e6Sstefano_zampini     } else {
491a59f3c4dSstefano_zampini       nid[i] = oid[i * bs] / bs;
49263fa5c83Sstefano_zampini     }
49363fa5c83Sstefano_zampini   }
4949566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreIndices(mapping, &oid));
495a59f3c4dSstefano_zampini 
49663fa5c83Sstefano_zampini   mapping->n  = nn;
49763fa5c83Sstefano_zampini   mapping->bs = bs;
4989566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping->indices));
49963fa5c83Sstefano_zampini   mapping->indices     = nid;
500c9345713Sstefano_zampini   mapping->globalstart = 0;
501c9345713Sstefano_zampini   mapping->globalend   = 0;
5021bd0b88eSStefano Zampini 
5031bd0b88eSStefano Zampini   /* reset the cached information */
5049566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping->info_procs));
5059566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping->info_numprocs));
5061bd0b88eSStefano Zampini   if (mapping->info_indices) {
5071bd0b88eSStefano Zampini     PetscInt i;
5081bd0b88eSStefano Zampini 
5099566063dSJacob Faibussowitsch     PetscCall(PetscFree((mapping->info_indices)[0]));
510*48a46eb9SPierre Jolivet     for (i = 1; i < mapping->info_nproc; i++) PetscCall(PetscFree(mapping->info_indices[i]));
5119566063dSJacob Faibussowitsch     PetscCall(PetscFree(mapping->info_indices));
5121bd0b88eSStefano Zampini   }
5131bd0b88eSStefano Zampini   mapping->info_cached = PETSC_FALSE;
5141bd0b88eSStefano Zampini 
515dbbe0bcdSBarry Smith   PetscTryTypeMethod(mapping, destroy);
51663fa5c83Sstefano_zampini   PetscFunctionReturn(0);
51763fa5c83Sstefano_zampini }
51863fa5c83Sstefano_zampini 
51945b6f7e9SBarry Smith /*@
52045b6f7e9SBarry Smith     ISLocalToGlobalMappingGetBlockSize - Gets the blocksize of the mapping
52145b6f7e9SBarry Smith     ordering and a global parallel ordering.
52245b6f7e9SBarry Smith 
52345b6f7e9SBarry Smith     Not Collective
52445b6f7e9SBarry Smith 
52545b6f7e9SBarry Smith     Input Parameters:
52645b6f7e9SBarry Smith .   mapping - mapping data structure
52745b6f7e9SBarry Smith 
52845b6f7e9SBarry Smith     Output Parameter:
52945b6f7e9SBarry Smith .   bs - the blocksize
53045b6f7e9SBarry Smith 
53145b6f7e9SBarry Smith     Level: advanced
53245b6f7e9SBarry Smith 
533db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`
53445b6f7e9SBarry Smith @*/
5359371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping mapping, PetscInt *bs) {
53645b6f7e9SBarry Smith   PetscFunctionBegin;
537cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
53845b6f7e9SBarry Smith   *bs = mapping->bs;
53945b6f7e9SBarry Smith   PetscFunctionReturn(0);
54045b6f7e9SBarry Smith }
54145b6f7e9SBarry Smith 
542ba5bb76aSSatish Balay /*@
54390f02eecSBarry Smith     ISLocalToGlobalMappingCreate - Creates a mapping between a local (0 to n)
54490f02eecSBarry Smith     ordering and a global parallel ordering.
5452362add9SBarry Smith 
54689d82c54SBarry Smith     Not Collective, but communicator may have more than one process
547b9cd556bSLois Curfman McInnes 
5482362add9SBarry Smith     Input Parameters:
54989d82c54SBarry Smith +   comm - MPI communicator
550f0413b6fSBarry Smith .   bs - the block size
55128bc9809SBarry Smith .   n - the number of local elements divided by the block size, or equivalently the number of block indices
55228bc9809SBarry 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
553d5ad8652SBarry Smith -   mode - see PetscCopyMode
5542362add9SBarry Smith 
555a997ad1aSLois Curfman McInnes     Output Parameter:
55690f02eecSBarry Smith .   mapping - new mapping data structure
5572362add9SBarry Smith 
55895452b02SPatrick Sanan     Notes:
55995452b02SPatrick Sanan     There is one integer value in indices per block and it represents the actual indices bs*idx + j, where j=0,..,bs-1
560413f72f0SBarry Smith 
5619a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
562413f72f0SBarry 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.
563413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
564413f72f0SBarry Smith 
565a997ad1aSLois Curfman McInnes     Level: advanced
566a997ad1aSLois Curfman McInnes 
567db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingSetFromOptions()`, `ISLOCALTOGLOBALMAPPINGBASIC`, `ISLOCALTOGLOBALMAPPINGHASH`
568db781477SPatrick Sanan           `ISLocalToGlobalMappingSetType()`, `ISLocalToGlobalMappingType`
5692362add9SBarry Smith @*/
5709371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingCreate(MPI_Comm comm, PetscInt bs, PetscInt n, const PetscInt indices[], PetscCopyMode mode, ISLocalToGlobalMapping *mapping) {
57132dcc486SBarry Smith   PetscInt *in;
572b46b645bSBarry Smith 
573b46b645bSBarry Smith   PetscFunctionBegin;
574064a246eSJacob Faibussowitsch   if (n) PetscValidIntPointer(indices, 4);
575064a246eSJacob Faibussowitsch   PetscValidPointer(mapping, 6);
576b46b645bSBarry Smith 
5770298fd71SBarry Smith   *mapping = NULL;
5789566063dSJacob Faibussowitsch   PetscCall(ISInitializePackage());
5792362add9SBarry Smith 
5809566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(*mapping, IS_LTOGM_CLASSID, "ISLocalToGlobalMapping", "Local to global mapping", "IS", comm, ISLocalToGlobalMappingDestroy, ISLocalToGlobalMappingView));
581d4bb536fSBarry Smith   (*mapping)->n  = n;
582f0413b6fSBarry Smith   (*mapping)->bs = bs;
583d5ad8652SBarry Smith   if (mode == PETSC_COPY_VALUES) {
5849566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &in));
5859566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(in, indices, n));
586d5ad8652SBarry Smith     (*mapping)->indices         = in;
58771910c26SVaclav Hapla     (*mapping)->dealloc_indices = PETSC_TRUE;
5889566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectMemory((PetscObject)*mapping, n * sizeof(PetscInt)));
5896389a1a1SBarry Smith   } else if (mode == PETSC_OWN_POINTER) {
5906389a1a1SBarry Smith     (*mapping)->indices         = (PetscInt *)indices;
59171910c26SVaclav Hapla     (*mapping)->dealloc_indices = PETSC_TRUE;
5929566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectMemory((PetscObject)*mapping, n * sizeof(PetscInt)));
59371910c26SVaclav Hapla   } else if (mode == PETSC_USE_POINTER) {
59471910c26SVaclav Hapla     (*mapping)->indices = (PetscInt *)indices;
5959371c9d4SSatish Balay   } else SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid mode %d", mode);
5963a40ed3dSBarry Smith   PetscFunctionReturn(0);
5972362add9SBarry Smith }
5982362add9SBarry Smith 
599413f72f0SBarry Smith PetscFunctionList ISLocalToGlobalMappingList = NULL;
600413f72f0SBarry Smith 
60190f02eecSBarry Smith /*@
6027e99dc12SLawrence Mitchell    ISLocalToGlobalMappingSetFromOptions - Set mapping options from the options database.
6037e99dc12SLawrence Mitchell 
6047e99dc12SLawrence Mitchell    Not collective
6057e99dc12SLawrence Mitchell 
6067e99dc12SLawrence Mitchell    Input Parameters:
6077e99dc12SLawrence Mitchell .  mapping - mapping data structure
6087e99dc12SLawrence Mitchell 
6097e99dc12SLawrence Mitchell    Level: advanced
6107e99dc12SLawrence Mitchell 
6117e99dc12SLawrence Mitchell @*/
6129371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingSetFromOptions(ISLocalToGlobalMapping mapping) {
613413f72f0SBarry Smith   char                       type[256];
614413f72f0SBarry Smith   ISLocalToGlobalMappingType defaulttype = "Not set";
6157e99dc12SLawrence Mitchell   PetscBool                  flg;
6167e99dc12SLawrence Mitchell 
6177e99dc12SLawrence Mitchell   PetscFunctionBegin;
6187e99dc12SLawrence Mitchell   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
6199566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRegisterAll());
620d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)mapping);
6219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-islocaltoglobalmapping_type", "ISLocalToGlobalMapping method", "ISLocalToGlobalMappingSetType", ISLocalToGlobalMappingList, (char *)(((PetscObject)mapping)->type_name) ? ((PetscObject)mapping)->type_name : defaulttype, type, 256, &flg));
6221baa6e33SBarry Smith   if (flg) PetscCall(ISLocalToGlobalMappingSetType(mapping, type));
623d0609cedSBarry Smith   PetscOptionsEnd();
6247e99dc12SLawrence Mitchell   PetscFunctionReturn(0);
6257e99dc12SLawrence Mitchell }
6267e99dc12SLawrence Mitchell 
6277e99dc12SLawrence Mitchell /*@
62890f02eecSBarry Smith    ISLocalToGlobalMappingDestroy - Destroys a mapping between a local (0 to n)
62990f02eecSBarry Smith    ordering and a global parallel ordering.
63090f02eecSBarry Smith 
6310f5bd95cSBarry Smith    Note Collective
632b9cd556bSLois Curfman McInnes 
63390f02eecSBarry Smith    Input Parameters:
63490f02eecSBarry Smith .  mapping - mapping data structure
63590f02eecSBarry Smith 
636a997ad1aSLois Curfman McInnes    Level: advanced
637a997ad1aSLois Curfman McInnes 
638db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`
63990f02eecSBarry Smith @*/
6409371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping *mapping) {
6413a40ed3dSBarry Smith   PetscFunctionBegin;
6426bf464f9SBarry Smith   if (!*mapping) PetscFunctionReturn(0);
6436bf464f9SBarry Smith   PetscValidHeaderSpecific((*mapping), IS_LTOGM_CLASSID, 1);
6449371c9d4SSatish Balay   if (--((PetscObject)(*mapping))->refct > 0) {
6459371c9d4SSatish Balay     *mapping = NULL;
6469371c9d4SSatish Balay     PetscFunctionReturn(0);
64771910c26SVaclav Hapla   }
648*48a46eb9SPierre Jolivet   if ((*mapping)->dealloc_indices) PetscCall(PetscFree((*mapping)->indices));
6499566063dSJacob Faibussowitsch   PetscCall(PetscFree((*mapping)->info_procs));
6509566063dSJacob Faibussowitsch   PetscCall(PetscFree((*mapping)->info_numprocs));
651268a049cSStefano Zampini   if ((*mapping)->info_indices) {
652268a049cSStefano Zampini     PetscInt i;
653268a049cSStefano Zampini 
6549566063dSJacob Faibussowitsch     PetscCall(PetscFree(((*mapping)->info_indices)[0]));
655*48a46eb9SPierre Jolivet     for (i = 1; i < (*mapping)->info_nproc; i++) PetscCall(PetscFree(((*mapping)->info_indices)[i]));
6569566063dSJacob Faibussowitsch     PetscCall(PetscFree((*mapping)->info_indices));
657268a049cSStefano Zampini   }
658*48a46eb9SPierre Jolivet   if ((*mapping)->info_nodei) PetscCall(PetscFree(((*mapping)->info_nodei)[0]));
6599566063dSJacob Faibussowitsch   PetscCall(PetscFree2((*mapping)->info_nodec, (*mapping)->info_nodei));
660*48a46eb9SPierre Jolivet   if ((*mapping)->ops->destroy) PetscCall((*(*mapping)->ops->destroy)(*mapping));
6619566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(mapping));
6624c8fdceaSLisandro Dalcin   *mapping = NULL;
6633a40ed3dSBarry Smith   PetscFunctionReturn(0);
66490f02eecSBarry Smith }
66590f02eecSBarry Smith 
66690f02eecSBarry Smith /*@
6673acfe500SLois Curfman McInnes     ISLocalToGlobalMappingApplyIS - Creates from an IS in the local numbering
6683acfe500SLois Curfman McInnes     a new index set using the global numbering defined in an ISLocalToGlobalMapping
6693acfe500SLois Curfman McInnes     context.
67090f02eecSBarry Smith 
6714cb36875SStefano Zampini     Collective on is
672b9cd556bSLois Curfman McInnes 
67390f02eecSBarry Smith     Input Parameters:
674b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
675b9cd556bSLois Curfman McInnes -   is - index set in local numbering
67690f02eecSBarry Smith 
67790f02eecSBarry Smith     Output Parameters:
67890f02eecSBarry Smith .   newis - index set in global numbering
67990f02eecSBarry Smith 
6804cb36875SStefano Zampini     Notes:
6814cb36875SStefano Zampini     The output IS will have the same communicator of the input IS.
6824cb36875SStefano Zampini 
683a997ad1aSLois Curfman McInnes     Level: advanced
684a997ad1aSLois Curfman McInnes 
685db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingCreate()`,
686db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`, `ISGlobalToLocalMappingApply()`
68790f02eecSBarry Smith @*/
6889371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping, IS is, IS *newis) {
689e24637baSBarry Smith   PetscInt        n, *idxout;
6905d0c19d7SBarry Smith   const PetscInt *idxin;
6913a40ed3dSBarry Smith 
6923a40ed3dSBarry Smith   PetscFunctionBegin;
6930700a824SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
6940700a824SBarry Smith   PetscValidHeaderSpecific(is, IS_CLASSID, 2);
6954482741eSBarry Smith   PetscValidPointer(newis, 3);
69690f02eecSBarry Smith 
6979566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is, &n));
6989566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is, &idxin));
6999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &idxout));
7009566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(mapping, n, idxin, idxout));
7019566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is, &idxin));
7029566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)is), n, idxout, PETSC_OWN_POINTER, newis));
7033a40ed3dSBarry Smith   PetscFunctionReturn(0);
70490f02eecSBarry Smith }
70590f02eecSBarry Smith 
706b89cb25eSSatish Balay /*@
7073acfe500SLois Curfman McInnes    ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering
7083acfe500SLois Curfman McInnes    and converts them to the global numbering.
70990f02eecSBarry Smith 
710b9cd556bSLois Curfman McInnes    Not collective
711b9cd556bSLois Curfman McInnes 
712bb25748dSBarry Smith    Input Parameters:
713b9cd556bSLois Curfman McInnes +  mapping - the local to global mapping context
714bb25748dSBarry Smith .  N - number of integers
715b9cd556bSLois Curfman McInnes -  in - input indices in local numbering
716bb25748dSBarry Smith 
717bb25748dSBarry Smith    Output Parameter:
718bb25748dSBarry Smith .  out - indices in global numbering
719bb25748dSBarry Smith 
720b9cd556bSLois Curfman McInnes    Notes:
721b9cd556bSLois Curfman McInnes    The in and out array parameters may be identical.
722d4bb536fSBarry Smith 
723a997ad1aSLois Curfman McInnes    Level: advanced
724a997ad1aSLois Curfman McInnes 
725c2e3fba1SPatrick Sanan .seealso: `ISLocalToGlobalMappingApplyBlock()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingDestroy()`,
726c2e3fba1SPatrick Sanan           `ISLocalToGlobalMappingApplyIS()`, `AOCreateBasic()`, `AOApplicationToPetsc()`,
727db781477SPatrick Sanan           `AOPetscToApplication()`, `ISGlobalToLocalMappingApply()`
728bb25748dSBarry Smith 
729afcb2eb5SJed Brown @*/
7309371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping, PetscInt N, const PetscInt in[], PetscInt out[]) {
731cbc1caf0SMatthew G. Knepley   PetscInt i, bs, Nmax;
73245b6f7e9SBarry Smith 
73345b6f7e9SBarry Smith   PetscFunctionBegin;
734cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
735cbc1caf0SMatthew G. Knepley   bs   = mapping->bs;
736cbc1caf0SMatthew G. Knepley   Nmax = bs * mapping->n;
73745b6f7e9SBarry Smith   if (bs == 1) {
738cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
73945b6f7e9SBarry Smith     for (i = 0; i < N; i++) {
74045b6f7e9SBarry Smith       if (in[i] < 0) {
74145b6f7e9SBarry Smith         out[i] = in[i];
74245b6f7e9SBarry Smith         continue;
74345b6f7e9SBarry Smith       }
74408401ef6SPierre Jolivet       PetscCheck(in[i] < Nmax, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local index %" PetscInt_FMT " too large %" PetscInt_FMT " (max) at %" PetscInt_FMT, in[i], Nmax - 1, i);
74545b6f7e9SBarry Smith       out[i] = idx[in[i]];
74645b6f7e9SBarry Smith     }
74745b6f7e9SBarry Smith   } else {
748cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
74945b6f7e9SBarry Smith     for (i = 0; i < N; i++) {
75045b6f7e9SBarry Smith       if (in[i] < 0) {
75145b6f7e9SBarry Smith         out[i] = in[i];
75245b6f7e9SBarry Smith         continue;
75345b6f7e9SBarry Smith       }
75408401ef6SPierre Jolivet       PetscCheck(in[i] < Nmax, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local index %" PetscInt_FMT " too large %" PetscInt_FMT " (max) at %" PetscInt_FMT, in[i], Nmax - 1, i);
75545b6f7e9SBarry Smith       out[i] = idx[in[i] / bs] * bs + (in[i] % bs);
75645b6f7e9SBarry Smith     }
75745b6f7e9SBarry Smith   }
75845b6f7e9SBarry Smith   PetscFunctionReturn(0);
75945b6f7e9SBarry Smith }
76045b6f7e9SBarry Smith 
76145b6f7e9SBarry Smith /*@
7626006e8d2SBarry Smith    ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local block numbering and converts them to the global block numbering
76345b6f7e9SBarry Smith 
76445b6f7e9SBarry Smith    Not collective
76545b6f7e9SBarry Smith 
76645b6f7e9SBarry Smith    Input Parameters:
76745b6f7e9SBarry Smith +  mapping - the local to global mapping context
76845b6f7e9SBarry Smith .  N - number of integers
7696006e8d2SBarry Smith -  in - input indices in local block numbering
77045b6f7e9SBarry Smith 
77145b6f7e9SBarry Smith    Output Parameter:
7726006e8d2SBarry Smith .  out - indices in global block numbering
77345b6f7e9SBarry Smith 
77445b6f7e9SBarry Smith    Notes:
77545b6f7e9SBarry Smith    The in and out array parameters may be identical.
77645b6f7e9SBarry Smith 
7776006e8d2SBarry Smith    Example:
7786006e8d2SBarry 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
7796006e8d2SBarry Smith      (the first block) would produce 0 and the mapping applied to 1 (the second block) would produce 3.
7806006e8d2SBarry Smith 
78145b6f7e9SBarry Smith    Level: advanced
78245b6f7e9SBarry Smith 
783c2e3fba1SPatrick Sanan .seealso: `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingDestroy()`,
784c2e3fba1SPatrick Sanan           `ISLocalToGlobalMappingApplyIS()`, `AOCreateBasic()`, `AOApplicationToPetsc()`,
785db781477SPatrick Sanan           `AOPetscToApplication()`, `ISGlobalToLocalMappingApply()`
78645b6f7e9SBarry Smith 
78745b6f7e9SBarry Smith @*/
7889371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping, PetscInt N, const PetscInt in[], PetscInt out[]) {
7898a1f772fSStefano Zampini   PetscInt        i, Nmax;
7908a1f772fSStefano Zampini   const PetscInt *idx;
791d4bb536fSBarry Smith 
792a0d79125SStefano Zampini   PetscFunctionBegin;
793a0d79125SStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
7948a1f772fSStefano Zampini   Nmax = mapping->n;
7958a1f772fSStefano Zampini   idx  = mapping->indices;
796afcb2eb5SJed Brown   for (i = 0; i < N; i++) {
797afcb2eb5SJed Brown     if (in[i] < 0) {
798afcb2eb5SJed Brown       out[i] = in[i];
799afcb2eb5SJed Brown       continue;
800afcb2eb5SJed Brown     }
80108401ef6SPierre Jolivet     PetscCheck(in[i] < Nmax, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local block index %" PetscInt_FMT " too large %" PetscInt_FMT " (max) at %" PetscInt_FMT, in[i], Nmax - 1, i);
802afcb2eb5SJed Brown     out[i] = idx[in[i]];
803afcb2eb5SJed Brown   }
804afcb2eb5SJed Brown   PetscFunctionReturn(0);
805afcb2eb5SJed Brown }
806d4bb536fSBarry Smith 
8077e99dc12SLawrence Mitchell /*@
808a997ad1aSLois Curfman McInnes     ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers
809a997ad1aSLois Curfman McInnes     specified with a global numbering.
810d4bb536fSBarry Smith 
811b9cd556bSLois Curfman McInnes     Not collective
812b9cd556bSLois Curfman McInnes 
813d4bb536fSBarry Smith     Input Parameters:
814b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
8150040bde1SJunchao Zhang .   type - IS_GTOLM_MASK - maps global indices with no local value to -1 in the output list (i.e., mask them)
816d4bb536fSBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
817d4bb536fSBarry Smith .   n - number of global indices to map
818b9cd556bSLois Curfman McInnes -   idx - global indices to map
819d4bb536fSBarry Smith 
820d4bb536fSBarry Smith     Output Parameters:
821b9cd556bSLois Curfman McInnes +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
822b9cd556bSLois Curfman McInnes -   idxout - local index of each global index, one must pass in an array long enough
823e182c471SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApply() with
8240298fd71SBarry Smith              idxout == NULL to determine the required length (returned in nout)
825e182c471SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApply()
826e182c471SBarry Smith              a second time to set the values.
827d4bb536fSBarry Smith 
828b9cd556bSLois Curfman McInnes     Notes:
8290298fd71SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
830d4bb536fSBarry Smith 
8319a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
832413f72f0SBarry 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.
833413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
8340f5bd95cSBarry Smith 
835a997ad1aSLois Curfman McInnes     Level: advanced
836a997ad1aSLois Curfman McInnes 
83732fd6b96SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
83832fd6b96SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
83932fd6b96SBarry Smith 
840db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingApply()`, `ISGlobalToLocalMappingApplyBlock()`, `ISLocalToGlobalMappingCreate()`,
841db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`
842d4bb536fSBarry Smith @*/
8439371c9d4SSatish Balay PetscErrorCode ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, PetscInt n, const PetscInt idx[], PetscInt *nout, PetscInt idxout[]) {
8449d90f715SBarry Smith   PetscFunctionBegin;
8459d90f715SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
846*48a46eb9SPierre Jolivet   if (!mapping->data) PetscCall(ISGlobalToLocalMappingSetUp(mapping));
847dbbe0bcdSBarry Smith   PetscUseTypeMethod(mapping, globaltolocalmappingapply, type, n, idx, nout, idxout);
8489d90f715SBarry Smith   PetscFunctionReturn(0);
8499d90f715SBarry Smith }
8509d90f715SBarry Smith 
851d4fe737eSStefano Zampini /*@
852d4fe737eSStefano Zampini     ISGlobalToLocalMappingApplyIS - Creates from an IS in the global numbering
853d4fe737eSStefano Zampini     a new index set using the local numbering defined in an ISLocalToGlobalMapping
854d4fe737eSStefano Zampini     context.
855d4fe737eSStefano Zampini 
856d4fe737eSStefano Zampini     Not collective
857d4fe737eSStefano Zampini 
858d4fe737eSStefano Zampini     Input Parameters:
859d4fe737eSStefano Zampini +   mapping - mapping between local and global numbering
8600040bde1SJunchao Zhang .   type - IS_GTOLM_MASK - maps global indices with no local value to -1 in the output list (i.e., mask them)
8612785b321SStefano Zampini            IS_GTOLM_DROP - drops the indices with no local value from the output list
862d4fe737eSStefano Zampini -   is - index set in global numbering
863d4fe737eSStefano Zampini 
864d4fe737eSStefano Zampini     Output Parameters:
865d4fe737eSStefano Zampini .   newis - index set in local numbering
866d4fe737eSStefano Zampini 
8674cb36875SStefano Zampini     Notes:
8684cb36875SStefano Zampini     The output IS will be sequential, as it encodes a purely local operation
8694cb36875SStefano Zampini 
870d4fe737eSStefano Zampini     Level: advanced
871d4fe737eSStefano Zampini 
872db781477SPatrick Sanan .seealso: `ISGlobalToLocalMappingApply()`, `ISLocalToGlobalMappingCreate()`,
873db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`
874d4fe737eSStefano Zampini @*/
8759371c9d4SSatish Balay PetscErrorCode ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, IS is, IS *newis) {
876d4fe737eSStefano Zampini   PetscInt        n, nout, *idxout;
877d4fe737eSStefano Zampini   const PetscInt *idxin;
878d4fe737eSStefano Zampini 
879d4fe737eSStefano Zampini   PetscFunctionBegin;
880d4fe737eSStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
881d4fe737eSStefano Zampini   PetscValidHeaderSpecific(is, IS_CLASSID, 3);
882d4fe737eSStefano Zampini   PetscValidPointer(newis, 4);
883d4fe737eSStefano Zampini 
8849566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is, &n));
8859566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is, &idxin));
886d4fe737eSStefano Zampini   if (type == IS_GTOLM_MASK) {
8879566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &idxout));
888d4fe737eSStefano Zampini   } else {
8899566063dSJacob Faibussowitsch     PetscCall(ISGlobalToLocalMappingApply(mapping, type, n, idxin, &nout, NULL));
8909566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nout, &idxout));
891d4fe737eSStefano Zampini   }
8929566063dSJacob Faibussowitsch   PetscCall(ISGlobalToLocalMappingApply(mapping, type, n, idxin, &nout, idxout));
8939566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is, &idxin));
8949566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nout, idxout, PETSC_OWN_POINTER, newis));
895d4fe737eSStefano Zampini   PetscFunctionReturn(0);
896d4fe737eSStefano Zampini }
897d4fe737eSStefano Zampini 
8989d90f715SBarry Smith /*@
8999d90f715SBarry Smith     ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers
9009d90f715SBarry Smith     specified with a block global numbering.
9019d90f715SBarry Smith 
9029d90f715SBarry Smith     Not collective
9039d90f715SBarry Smith 
9049d90f715SBarry Smith     Input Parameters:
9059d90f715SBarry Smith +   mapping - mapping between local and global numbering
9060040bde1SJunchao Zhang .   type - IS_GTOLM_MASK - maps global indices with no local value to -1 in the output list (i.e., mask them)
9079d90f715SBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
9089d90f715SBarry Smith .   n - number of global indices to map
9099d90f715SBarry Smith -   idx - global indices to map
9109d90f715SBarry Smith 
9119d90f715SBarry Smith     Output Parameters:
9129d90f715SBarry Smith +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
9139d90f715SBarry Smith -   idxout - local index of each global index, one must pass in an array long enough
9149d90f715SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApplyBlock() with
9159d90f715SBarry Smith              idxout == NULL to determine the required length (returned in nout)
9169d90f715SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApplyBlock()
9179d90f715SBarry Smith              a second time to set the values.
9189d90f715SBarry Smith 
9199d90f715SBarry Smith     Notes:
9209d90f715SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
9219d90f715SBarry Smith 
9229a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
923413f72f0SBarry 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.
924413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
9259d90f715SBarry Smith 
9269d90f715SBarry Smith     Level: advanced
9279d90f715SBarry Smith 
9289d90f715SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
9299d90f715SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
9309d90f715SBarry Smith 
931db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingApply()`, `ISGlobalToLocalMappingApply()`, `ISLocalToGlobalMappingCreate()`,
932db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`
9339d90f715SBarry Smith @*/
9349371c9d4SSatish Balay PetscErrorCode ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, PetscInt n, const PetscInt idx[], PetscInt *nout, PetscInt idxout[]) {
9353a40ed3dSBarry Smith   PetscFunctionBegin;
9360700a824SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
937*48a46eb9SPierre Jolivet   if (!mapping->data) PetscCall(ISGlobalToLocalMappingSetUp(mapping));
938dbbe0bcdSBarry Smith   PetscUseTypeMethod(mapping, globaltolocalmappingapplyblock, type, n, idx, nout, idxout);
9393a40ed3dSBarry Smith   PetscFunctionReturn(0);
940d4bb536fSBarry Smith }
94190f02eecSBarry Smith 
94289d82c54SBarry Smith /*@C
9436a818285SBarry Smith     ISLocalToGlobalMappingGetBlockInfo - Gets the neighbor information for each processor and
94489d82c54SBarry Smith      each index shared by more than one processor
94589d82c54SBarry Smith 
94689d82c54SBarry Smith     Collective on ISLocalToGlobalMapping
94789d82c54SBarry Smith 
948f899ff85SJose E. Roman     Input Parameter:
94989d82c54SBarry Smith .   mapping - the mapping from local to global indexing
95089d82c54SBarry Smith 
951d8d19677SJose E. Roman     Output Parameters:
95289d82c54SBarry Smith +   nproc - number of processors that are connected to this one
95389d82c54SBarry Smith .   proc - neighboring processors
95407b52d57SBarry Smith .   numproc - number of indices for each subdomain (processor)
9553463a7baSJed Brown -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
95689d82c54SBarry Smith 
95789d82c54SBarry Smith     Level: advanced
95889d82c54SBarry Smith 
9592cfcea29SBarry Smith     Fortran Usage:
9602cfcea29SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
9612cfcea29SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
9622cfcea29SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
9632cfcea29SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
9642cfcea29SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
9652cfcea29SBarry Smith 
966db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
967db781477SPatrick Sanan           `ISLocalToGlobalMappingRestoreInfo()`
96889d82c54SBarry Smith @*/
9699371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[]) {
970268a049cSStefano Zampini   PetscFunctionBegin;
971268a049cSStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
972268a049cSStefano Zampini   if (mapping->info_cached) {
973268a049cSStefano Zampini     *nproc    = mapping->info_nproc;
974268a049cSStefano Zampini     *procs    = mapping->info_procs;
975268a049cSStefano Zampini     *numprocs = mapping->info_numprocs;
976268a049cSStefano Zampini     *indices  = mapping->info_indices;
977268a049cSStefano Zampini   } else {
9789566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetBlockInfo_Private(mapping, nproc, procs, numprocs, indices));
979268a049cSStefano Zampini   }
980268a049cSStefano Zampini   PetscFunctionReturn(0);
981268a049cSStefano Zampini }
982268a049cSStefano Zampini 
9839371c9d4SSatish Balay static PetscErrorCode ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[]) {
98497f1f81fSBarry Smith   PetscMPIInt  size, rank, tag1, tag2, tag3, *len, *source, imdex;
98532dcc486SBarry Smith   PetscInt     i, n = mapping->n, Ng, ng, max = 0, *lindices = mapping->indices;
98632dcc486SBarry Smith   PetscInt    *nprocs, *owner, nsends, *sends, j, *starts, nmax, nrecvs, *recvs, proc;
987c599c493SJunchao Zhang   PetscInt     cnt, scale, *ownedsenders, *nownedsenders, rstart;
98832dcc486SBarry Smith   PetscInt     node, nownedm, nt, *sends2, nsends2, *starts2, *lens2, *dest, nrecvs2, *starts3, *recvs2, k, *bprocs, *tmp;
98932dcc486SBarry Smith   PetscInt     first_procs, first_numprocs, *first_indices;
99089d82c54SBarry Smith   MPI_Request *recv_waits, *send_waits;
99130dcb7c9SBarry Smith   MPI_Status   recv_status, *send_status, *recv_statuses;
992ce94432eSBarry Smith   MPI_Comm     comm;
993ace3abfcSBarry Smith   PetscBool    debug = PETSC_FALSE;
99489d82c54SBarry Smith 
99589d82c54SBarry Smith   PetscFunctionBegin;
9969566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mapping, &comm));
9979566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
9989566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
99924cf384cSBarry Smith   if (size == 1) {
100024cf384cSBarry Smith     *nproc = 0;
10010298fd71SBarry Smith     *procs = NULL;
10029566063dSJacob Faibussowitsch     PetscCall(PetscNew(numprocs));
10031e2105dcSBarry Smith     (*numprocs)[0] = 0;
10049566063dSJacob Faibussowitsch     PetscCall(PetscNew(indices));
10050298fd71SBarry Smith     (*indices)[0]          = NULL;
1006268a049cSStefano Zampini     /* save info for reuse */
1007268a049cSStefano Zampini     mapping->info_nproc    = *nproc;
1008268a049cSStefano Zampini     mapping->info_procs    = *procs;
1009268a049cSStefano Zampini     mapping->info_numprocs = *numprocs;
1010268a049cSStefano Zampini     mapping->info_indices  = *indices;
1011268a049cSStefano Zampini     mapping->info_cached   = PETSC_TRUE;
101224cf384cSBarry Smith     PetscFunctionReturn(0);
101324cf384cSBarry Smith   }
101424cf384cSBarry Smith 
10159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)mapping)->options, NULL, "-islocaltoglobalmappinggetinfo_debug", &debug, NULL));
101607b52d57SBarry Smith 
10173677ff5aSBarry Smith   /*
10186a818285SBarry Smith     Notes on ISLocalToGlobalMappingGetBlockInfo
10193677ff5aSBarry Smith 
10203677ff5aSBarry Smith     globally owned node - the nodes that have been assigned to this processor in global
10213677ff5aSBarry Smith            numbering, just for this routine.
10223677ff5aSBarry Smith 
10233677ff5aSBarry Smith     nontrivial globally owned node - node assigned to this processor that is on a subdomain
10243677ff5aSBarry Smith            boundary (i.e. is has more than one local owner)
10253677ff5aSBarry Smith 
10263677ff5aSBarry Smith     locally owned node - node that exists on this processors subdomain
10273677ff5aSBarry Smith 
10283677ff5aSBarry Smith     nontrivial locally owned node - node that is not in the interior (i.e. has more than one
10293677ff5aSBarry Smith            local subdomain
10303677ff5aSBarry Smith   */
10319566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetNewTag((PetscObject)mapping, &tag1));
10329566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetNewTag((PetscObject)mapping, &tag2));
10339566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetNewTag((PetscObject)mapping, &tag3));
103489d82c54SBarry Smith 
103589d82c54SBarry Smith   for (i = 0; i < n; i++) {
103689d82c54SBarry Smith     if (lindices[i] > max) max = lindices[i];
103789d82c54SBarry Smith   }
10381c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&max, &Ng, 1, MPIU_INT, MPI_MAX, comm));
103978058e43SBarry Smith   Ng++;
10409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
10419566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1042bc8ff85bSBarry Smith   scale = Ng / size + 1;
10439371c9d4SSatish Balay   ng    = scale;
10449371c9d4SSatish Balay   if (rank == size - 1) ng = Ng - scale * (size - 1);
10459371c9d4SSatish Balay   ng     = PetscMax(1, ng);
1046caba0dd0SBarry Smith   rstart = scale * rank;
104789d82c54SBarry Smith 
104889d82c54SBarry Smith   /* determine ownership ranges of global indices */
10499566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * size, &nprocs));
10509566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(nprocs, 2 * size));
105189d82c54SBarry Smith 
105289d82c54SBarry Smith   /* determine owners of each local node  */
10539566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &owner));
105489d82c54SBarry Smith   for (i = 0; i < n; i++) {
10553677ff5aSBarry Smith     proc                 = lindices[i] / scale; /* processor that globally owns this index */
105627c402fcSBarry Smith     nprocs[2 * proc + 1] = 1;                   /* processor globally owns at least one of ours */
10573677ff5aSBarry Smith     owner[i]             = proc;
105827c402fcSBarry Smith     nprocs[2 * proc]++; /* count of how many that processor globally owns of ours */
105989d82c54SBarry Smith   }
10609371c9d4SSatish Balay   nsends = 0;
10619371c9d4SSatish Balay   for (i = 0; i < size; i++) nsends += nprocs[2 * i + 1];
10629566063dSJacob Faibussowitsch   PetscCall(PetscInfo(mapping, "Number of global owners for my local data %" PetscInt_FMT "\n", nsends));
106389d82c54SBarry Smith 
106489d82c54SBarry Smith   /* inform other processors of number of messages and max length*/
10659566063dSJacob Faibussowitsch   PetscCall(PetscMaxSum(comm, nprocs, &nmax, &nrecvs));
10669566063dSJacob Faibussowitsch   PetscCall(PetscInfo(mapping, "Number of local owners for my global data %" PetscInt_FMT "\n", nrecvs));
106789d82c54SBarry Smith 
106889d82c54SBarry Smith   /* post receives for owned rows */
10699566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((2 * nrecvs + 1) * (nmax + 1), &recvs));
10709566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs + 1, &recv_waits));
1071*48a46eb9SPierre Jolivet   for (i = 0; i < nrecvs; i++) PetscCallMPI(MPI_Irecv(recvs + 2 * nmax * i, 2 * nmax, MPIU_INT, MPI_ANY_SOURCE, tag1, comm, recv_waits + i));
107289d82c54SBarry Smith 
107389d82c54SBarry Smith   /* pack messages containing lists of local nodes to owners */
10749566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * n + 1, &sends));
10759566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size + 1, &starts));
107689d82c54SBarry Smith   starts[0] = 0;
1077f6e5521dSKarl Rupp   for (i = 1; i < size; i++) starts[i] = starts[i - 1] + 2 * nprocs[2 * i - 2];
107889d82c54SBarry Smith   for (i = 0; i < n; i++) {
107989d82c54SBarry Smith     sends[starts[owner[i]]++] = lindices[i];
108030dcb7c9SBarry Smith     sends[starts[owner[i]]++] = i;
108189d82c54SBarry Smith   }
10829566063dSJacob Faibussowitsch   PetscCall(PetscFree(owner));
108389d82c54SBarry Smith   starts[0] = 0;
1084f6e5521dSKarl Rupp   for (i = 1; i < size; i++) starts[i] = starts[i - 1] + 2 * nprocs[2 * i - 2];
108589d82c54SBarry Smith 
108689d82c54SBarry Smith   /* send the messages */
10879566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends + 1, &send_waits));
10889566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends + 1, &dest));
108989d82c54SBarry Smith   cnt = 0;
109089d82c54SBarry Smith   for (i = 0; i < size; i++) {
109127c402fcSBarry Smith     if (nprocs[2 * i]) {
10929566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Isend(sends + starts[i], 2 * nprocs[2 * i], MPIU_INT, i, tag1, comm, send_waits + cnt));
109330dcb7c9SBarry Smith       dest[cnt] = i;
109489d82c54SBarry Smith       cnt++;
109589d82c54SBarry Smith     }
109689d82c54SBarry Smith   }
10979566063dSJacob Faibussowitsch   PetscCall(PetscFree(starts));
109889d82c54SBarry Smith 
109989d82c54SBarry Smith   /* wait on receives */
11009566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs + 1, &source));
11019566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs + 1, &len));
110289d82c54SBarry Smith   cnt = nrecvs;
11039566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(ng + 1, &nownedsenders));
110489d82c54SBarry Smith   while (cnt) {
11059566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitany(nrecvs, recv_waits, &imdex, &recv_status));
110689d82c54SBarry Smith     /* unpack receives into our local space */
11079566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Get_count(&recv_status, MPIU_INT, &len[imdex]));
110889d82c54SBarry Smith     source[imdex] = recv_status.MPI_SOURCE;
110930dcb7c9SBarry Smith     len[imdex]    = len[imdex] / 2;
1110caba0dd0SBarry Smith     /* count how many local owners for each of my global owned indices */
111130dcb7c9SBarry Smith     for (i = 0; i < len[imdex]; i++) nownedsenders[recvs[2 * imdex * nmax + 2 * i] - rstart]++;
111289d82c54SBarry Smith     cnt--;
111389d82c54SBarry Smith   }
11149566063dSJacob Faibussowitsch   PetscCall(PetscFree(recv_waits));
111589d82c54SBarry Smith 
111630dcb7c9SBarry Smith   /* count how many globally owned indices are on an edge multiplied by how many processors own them. */
1117bc8ff85bSBarry Smith   nownedm = 0;
1118bc8ff85bSBarry Smith   for (i = 0; i < ng; i++) {
1119c599c493SJunchao Zhang     if (nownedsenders[i] > 1) nownedm += nownedsenders[i];
1120bc8ff85bSBarry Smith   }
1121bc8ff85bSBarry Smith 
1122bc8ff85bSBarry Smith   /* create single array to contain rank of all local owners of each globally owned index */
11239566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nownedm + 1, &ownedsenders));
11249566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ng + 1, &starts));
1125bc8ff85bSBarry Smith   starts[0] = 0;
1126bc8ff85bSBarry Smith   for (i = 1; i < ng; i++) {
1127bc8ff85bSBarry Smith     if (nownedsenders[i - 1] > 1) starts[i] = starts[i - 1] + nownedsenders[i - 1];
1128bc8ff85bSBarry Smith     else starts[i] = starts[i - 1];
1129bc8ff85bSBarry Smith   }
1130bc8ff85bSBarry Smith 
11316aad120cSJose E. Roman   /* for each nontrivial globally owned node list all arriving processors */
1132bc8ff85bSBarry Smith   for (i = 0; i < nrecvs; i++) {
1133bc8ff85bSBarry Smith     for (j = 0; j < len[i]; j++) {
113430dcb7c9SBarry Smith       node = recvs[2 * i * nmax + 2 * j] - rstart;
1135f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) ownedsenders[starts[node]++] = source[i];
1136bc8ff85bSBarry Smith     }
1137bc8ff85bSBarry Smith   }
1138bc8ff85bSBarry Smith 
113907b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
114030dcb7c9SBarry Smith     starts[0] = 0;
114130dcb7c9SBarry Smith     for (i = 1; i < ng; i++) {
114230dcb7c9SBarry Smith       if (nownedsenders[i - 1] > 1) starts[i] = starts[i - 1] + nownedsenders[i - 1];
114330dcb7c9SBarry Smith       else starts[i] = starts[i - 1];
114430dcb7c9SBarry Smith     }
114530dcb7c9SBarry Smith     for (i = 0; i < ng; i++) {
114630dcb7c9SBarry Smith       if (nownedsenders[i] > 1) {
11479566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(comm, "[%d] global node %" PetscInt_FMT " local owner processors: ", rank, i + rstart));
1148*48a46eb9SPierre Jolivet         for (j = 0; j < nownedsenders[i]; j++) PetscCall(PetscSynchronizedPrintf(comm, "%" PetscInt_FMT " ", ownedsenders[starts[i] + j]));
11499566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(comm, "\n"));
115030dcb7c9SBarry Smith       }
115130dcb7c9SBarry Smith     }
11529566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
115307b52d57SBarry Smith   } /* -----------------------------------  */
115430dcb7c9SBarry Smith 
11553677ff5aSBarry Smith   /* wait on original sends */
11563a96401aSBarry Smith   if (nsends) {
11579566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nsends, &send_status));
11589566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nsends, send_waits, send_status));
11599566063dSJacob Faibussowitsch     PetscCall(PetscFree(send_status));
11603a96401aSBarry Smith   }
11619566063dSJacob Faibussowitsch   PetscCall(PetscFree(send_waits));
11629566063dSJacob Faibussowitsch   PetscCall(PetscFree(sends));
11639566063dSJacob Faibussowitsch   PetscCall(PetscFree(nprocs));
11643677ff5aSBarry Smith 
11653677ff5aSBarry Smith   /* pack messages to send back to local owners */
116630dcb7c9SBarry Smith   starts[0] = 0;
116730dcb7c9SBarry Smith   for (i = 1; i < ng; i++) {
116830dcb7c9SBarry Smith     if (nownedsenders[i - 1] > 1) starts[i] = starts[i - 1] + nownedsenders[i - 1];
116930dcb7c9SBarry Smith     else starts[i] = starts[i - 1];
117030dcb7c9SBarry Smith   }
117130dcb7c9SBarry Smith   nsends2 = nrecvs;
11729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends2 + 1, &nprocs)); /* length of each message */
117330dcb7c9SBarry Smith   for (i = 0; i < nrecvs; i++) {
117430dcb7c9SBarry Smith     nprocs[i] = 1;
117530dcb7c9SBarry Smith     for (j = 0; j < len[i]; j++) {
117630dcb7c9SBarry Smith       node = recvs[2 * i * nmax + 2 * j] - rstart;
1177f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) nprocs[i] += 2 + nownedsenders[node];
117830dcb7c9SBarry Smith     }
117930dcb7c9SBarry Smith   }
1180f6e5521dSKarl Rupp   nt = 0;
1181f6e5521dSKarl Rupp   for (i = 0; i < nsends2; i++) nt += nprocs[i];
1182f6e5521dSKarl Rupp 
11839566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt + 1, &sends2));
11849566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends2 + 1, &starts2));
1185f6e5521dSKarl Rupp 
1186f6e5521dSKarl Rupp   starts2[0] = 0;
1187f6e5521dSKarl Rupp   for (i = 1; i < nsends2; i++) starts2[i] = starts2[i - 1] + nprocs[i - 1];
118830dcb7c9SBarry Smith   /*
118930dcb7c9SBarry Smith      Each message is 1 + nprocs[i] long, and consists of
119030dcb7c9SBarry Smith        (0) the number of nodes being sent back
119130dcb7c9SBarry Smith        (1) the local node number,
119230dcb7c9SBarry Smith        (2) the number of processors sharing it,
119330dcb7c9SBarry Smith        (3) the processors sharing it
119430dcb7c9SBarry Smith   */
119530dcb7c9SBarry Smith   for (i = 0; i < nsends2; i++) {
119630dcb7c9SBarry Smith     cnt                = 1;
119730dcb7c9SBarry Smith     sends2[starts2[i]] = 0;
119830dcb7c9SBarry Smith     for (j = 0; j < len[i]; j++) {
119930dcb7c9SBarry Smith       node = recvs[2 * i * nmax + 2 * j] - rstart;
120030dcb7c9SBarry Smith       if (nownedsenders[node] > 1) {
120130dcb7c9SBarry Smith         sends2[starts2[i]]++;
120230dcb7c9SBarry Smith         sends2[starts2[i] + cnt++] = recvs[2 * i * nmax + 2 * j + 1];
120330dcb7c9SBarry Smith         sends2[starts2[i] + cnt++] = nownedsenders[node];
12049566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(&sends2[starts2[i] + cnt], &ownedsenders[starts[node]], nownedsenders[node]));
120530dcb7c9SBarry Smith         cnt += nownedsenders[node];
120630dcb7c9SBarry Smith       }
120730dcb7c9SBarry Smith     }
120830dcb7c9SBarry Smith   }
120930dcb7c9SBarry Smith 
121030dcb7c9SBarry Smith   /* receive the message lengths */
121130dcb7c9SBarry Smith   nrecvs2 = nsends;
12129566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs2 + 1, &lens2));
12139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs2 + 1, &starts3));
12149566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs2 + 1, &recv_waits));
1215*48a46eb9SPierre Jolivet   for (i = 0; i < nrecvs2; i++) PetscCallMPI(MPI_Irecv(&lens2[i], 1, MPIU_INT, dest[i], tag2, comm, recv_waits + i));
1216d44834fbSBarry Smith 
12178a8e0b3aSBarry Smith   /* send the message lengths */
1218*48a46eb9SPierre Jolivet   for (i = 0; i < nsends2; i++) PetscCallMPI(MPI_Send(&nprocs[i], 1, MPIU_INT, source[i], tag2, comm));
12198a8e0b3aSBarry Smith 
1220d44834fbSBarry Smith   /* wait on receives of lens */
12210c468ba9SBarry Smith   if (nrecvs2) {
12229566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrecvs2, &recv_statuses));
12239566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrecvs2, recv_waits, recv_statuses));
12249566063dSJacob Faibussowitsch     PetscCall(PetscFree(recv_statuses));
12250c468ba9SBarry Smith   }
12269566063dSJacob Faibussowitsch   PetscCall(PetscFree(recv_waits));
1227d44834fbSBarry Smith 
122830dcb7c9SBarry Smith   starts3[0] = 0;
1229d44834fbSBarry Smith   nt         = 0;
123030dcb7c9SBarry Smith   for (i = 0; i < nrecvs2 - 1; i++) {
123130dcb7c9SBarry Smith     starts3[i + 1] = starts3[i] + lens2[i];
1232d44834fbSBarry Smith     nt += lens2[i];
123330dcb7c9SBarry Smith   }
123476466f69SStefano Zampini   if (nrecvs2) nt += lens2[nrecvs2 - 1];
1235d44834fbSBarry Smith 
12369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt + 1, &recvs2));
12379566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs2 + 1, &recv_waits));
1238*48a46eb9SPierre Jolivet   for (i = 0; i < nrecvs2; i++) PetscCallMPI(MPI_Irecv(recvs2 + starts3[i], lens2[i], MPIU_INT, dest[i], tag3, comm, recv_waits + i));
123930dcb7c9SBarry Smith 
124030dcb7c9SBarry Smith   /* send the messages */
12419566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends2 + 1, &send_waits));
1242*48a46eb9SPierre Jolivet   for (i = 0; i < nsends2; i++) PetscCallMPI(MPI_Isend(sends2 + starts2[i], nprocs[i], MPIU_INT, source[i], tag3, comm, send_waits + i));
124330dcb7c9SBarry Smith 
124430dcb7c9SBarry Smith   /* wait on receives */
12450c468ba9SBarry Smith   if (nrecvs2) {
12469566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrecvs2, &recv_statuses));
12479566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrecvs2, recv_waits, recv_statuses));
12489566063dSJacob Faibussowitsch     PetscCall(PetscFree(recv_statuses));
12490c468ba9SBarry Smith   }
12509566063dSJacob Faibussowitsch   PetscCall(PetscFree(recv_waits));
12519566063dSJacob Faibussowitsch   PetscCall(PetscFree(nprocs));
125230dcb7c9SBarry Smith 
125307b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
125430dcb7c9SBarry Smith     cnt = 0;
125530dcb7c9SBarry Smith     for (i = 0; i < nrecvs2; i++) {
125630dcb7c9SBarry Smith       nt = recvs2[cnt++];
125730dcb7c9SBarry Smith       for (j = 0; j < nt; j++) {
12589566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(comm, "[%d] local node %" PetscInt_FMT " number of subdomains %" PetscInt_FMT ": ", rank, recvs2[cnt], recvs2[cnt + 1]));
1259*48a46eb9SPierre Jolivet         for (k = 0; k < recvs2[cnt + 1]; k++) PetscCall(PetscSynchronizedPrintf(comm, "%" PetscInt_FMT " ", recvs2[cnt + 2 + k]));
126030dcb7c9SBarry Smith         cnt += 2 + recvs2[cnt + 1];
12619566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(comm, "\n"));
126230dcb7c9SBarry Smith       }
126330dcb7c9SBarry Smith     }
12649566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
126507b52d57SBarry Smith   } /* -----------------------------------  */
126630dcb7c9SBarry Smith 
126730dcb7c9SBarry Smith   /* count number subdomains for each local node */
12689566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(size, &nprocs));
126930dcb7c9SBarry Smith   cnt = 0;
127030dcb7c9SBarry Smith   for (i = 0; i < nrecvs2; i++) {
127130dcb7c9SBarry Smith     nt = recvs2[cnt++];
127230dcb7c9SBarry Smith     for (j = 0; j < nt; j++) {
1273f6e5521dSKarl Rupp       for (k = 0; k < recvs2[cnt + 1]; k++) nprocs[recvs2[cnt + 2 + k]]++;
127430dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt + 1];
127530dcb7c9SBarry Smith     }
127630dcb7c9SBarry Smith   }
12779371c9d4SSatish Balay   nt = 0;
12789371c9d4SSatish Balay   for (i = 0; i < size; i++) nt += (nprocs[i] > 0);
127930dcb7c9SBarry Smith   *nproc = nt;
12809566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt + 1, procs));
12819566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt + 1, numprocs));
12829566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt + 1, indices));
12830298fd71SBarry Smith   for (i = 0; i < nt + 1; i++) (*indices)[i] = NULL;
12849566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size, &bprocs));
128530dcb7c9SBarry Smith   cnt = 0;
128630dcb7c9SBarry Smith   for (i = 0; i < size; i++) {
128730dcb7c9SBarry Smith     if (nprocs[i] > 0) {
128830dcb7c9SBarry Smith       bprocs[i]        = cnt;
128930dcb7c9SBarry Smith       (*procs)[cnt]    = i;
129030dcb7c9SBarry Smith       (*numprocs)[cnt] = nprocs[i];
12919566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nprocs[i], &(*indices)[cnt]));
129230dcb7c9SBarry Smith       cnt++;
129330dcb7c9SBarry Smith     }
129430dcb7c9SBarry Smith   }
129530dcb7c9SBarry Smith 
129630dcb7c9SBarry Smith   /* make the list of subdomains for each nontrivial local node */
12979566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(*numprocs, nt));
129830dcb7c9SBarry Smith   cnt = 0;
129930dcb7c9SBarry Smith   for (i = 0; i < nrecvs2; i++) {
130030dcb7c9SBarry Smith     nt = recvs2[cnt++];
130130dcb7c9SBarry Smith     for (j = 0; j < nt; j++) {
1302f6e5521dSKarl Rupp       for (k = 0; k < recvs2[cnt + 1]; k++) (*indices)[bprocs[recvs2[cnt + 2 + k]]][(*numprocs)[bprocs[recvs2[cnt + 2 + k]]]++] = recvs2[cnt];
130330dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt + 1];
130430dcb7c9SBarry Smith     }
130530dcb7c9SBarry Smith   }
13069566063dSJacob Faibussowitsch   PetscCall(PetscFree(bprocs));
13079566063dSJacob Faibussowitsch   PetscCall(PetscFree(recvs2));
130830dcb7c9SBarry Smith 
130907b52d57SBarry Smith   /* sort the node indexing by their global numbers */
131007b52d57SBarry Smith   nt = *nproc;
131107b52d57SBarry Smith   for (i = 0; i < nt; i++) {
13129566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1((*numprocs)[i], &tmp));
1313f6e5521dSKarl Rupp     for (j = 0; j < (*numprocs)[i]; j++) tmp[j] = lindices[(*indices)[i][j]];
13149566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithArray((*numprocs)[i], tmp, (*indices)[i]));
13159566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmp));
131607b52d57SBarry Smith   }
131707b52d57SBarry Smith 
131807b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
131930dcb7c9SBarry Smith     nt = *nproc;
132030dcb7c9SBarry Smith     for (i = 0; i < nt; i++) {
13219566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedPrintf(comm, "[%d] subdomain %" PetscInt_FMT " number of indices %" PetscInt_FMT ": ", rank, (*procs)[i], (*numprocs)[i]));
1322*48a46eb9SPierre Jolivet       for (j = 0; j < (*numprocs)[i]; j++) PetscCall(PetscSynchronizedPrintf(comm, "%" PetscInt_FMT " ", (*indices)[i][j]));
13239566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedPrintf(comm, "\n"));
132430dcb7c9SBarry Smith     }
13259566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
132607b52d57SBarry Smith   } /* -----------------------------------  */
132730dcb7c9SBarry Smith 
132830dcb7c9SBarry Smith   /* wait on sends */
132930dcb7c9SBarry Smith   if (nsends2) {
13309566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nsends2, &send_status));
13319566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nsends2, send_waits, send_status));
13329566063dSJacob Faibussowitsch     PetscCall(PetscFree(send_status));
133330dcb7c9SBarry Smith   }
133430dcb7c9SBarry Smith 
13359566063dSJacob Faibussowitsch   PetscCall(PetscFree(starts3));
13369566063dSJacob Faibussowitsch   PetscCall(PetscFree(dest));
13379566063dSJacob Faibussowitsch   PetscCall(PetscFree(send_waits));
13383677ff5aSBarry Smith 
13399566063dSJacob Faibussowitsch   PetscCall(PetscFree(nownedsenders));
13409566063dSJacob Faibussowitsch   PetscCall(PetscFree(ownedsenders));
13419566063dSJacob Faibussowitsch   PetscCall(PetscFree(starts));
13429566063dSJacob Faibussowitsch   PetscCall(PetscFree(starts2));
13439566063dSJacob Faibussowitsch   PetscCall(PetscFree(lens2));
134489d82c54SBarry Smith 
13459566063dSJacob Faibussowitsch   PetscCall(PetscFree(source));
13469566063dSJacob Faibussowitsch   PetscCall(PetscFree(len));
13479566063dSJacob Faibussowitsch   PetscCall(PetscFree(recvs));
13489566063dSJacob Faibussowitsch   PetscCall(PetscFree(nprocs));
13499566063dSJacob Faibussowitsch   PetscCall(PetscFree(sends2));
135024cf384cSBarry Smith 
135124cf384cSBarry Smith   /* put the information about myself as the first entry in the list */
135224cf384cSBarry Smith   first_procs    = (*procs)[0];
135324cf384cSBarry Smith   first_numprocs = (*numprocs)[0];
135424cf384cSBarry Smith   first_indices  = (*indices)[0];
135524cf384cSBarry Smith   for (i = 0; i < *nproc; i++) {
135624cf384cSBarry Smith     if ((*procs)[i] == rank) {
135724cf384cSBarry Smith       (*procs)[0]    = (*procs)[i];
135824cf384cSBarry Smith       (*numprocs)[0] = (*numprocs)[i];
135924cf384cSBarry Smith       (*indices)[0]  = (*indices)[i];
136024cf384cSBarry Smith       (*procs)[i]    = first_procs;
136124cf384cSBarry Smith       (*numprocs)[i] = first_numprocs;
136224cf384cSBarry Smith       (*indices)[i]  = first_indices;
136324cf384cSBarry Smith       break;
136424cf384cSBarry Smith     }
136524cf384cSBarry Smith   }
1366268a049cSStefano Zampini 
1367268a049cSStefano Zampini   /* save info for reuse */
1368268a049cSStefano Zampini   mapping->info_nproc    = *nproc;
1369268a049cSStefano Zampini   mapping->info_procs    = *procs;
1370268a049cSStefano Zampini   mapping->info_numprocs = *numprocs;
1371268a049cSStefano Zampini   mapping->info_indices  = *indices;
1372268a049cSStefano Zampini   mapping->info_cached   = PETSC_TRUE;
137389d82c54SBarry Smith   PetscFunctionReturn(0);
137489d82c54SBarry Smith }
137589d82c54SBarry Smith 
13766a818285SBarry Smith /*@C
13776a818285SBarry Smith     ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by ISLocalToGlobalMappingGetBlockInfo()
13786a818285SBarry Smith 
13796a818285SBarry Smith     Collective on ISLocalToGlobalMapping
13806a818285SBarry Smith 
1381f899ff85SJose E. Roman     Input Parameter:
13826a818285SBarry Smith .   mapping - the mapping from local to global indexing
13836a818285SBarry Smith 
1384d8d19677SJose E. Roman     Output Parameters:
13856a818285SBarry Smith +   nproc - number of processors that are connected to this one
13866a818285SBarry Smith .   proc - neighboring processors
13876a818285SBarry Smith .   numproc - number of indices for each processor
13886a818285SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
13896a818285SBarry Smith 
13906a818285SBarry Smith     Level: advanced
13916a818285SBarry Smith 
1392db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1393db781477SPatrick Sanan           `ISLocalToGlobalMappingGetInfo()`
13946a818285SBarry Smith @*/
13959371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[]) {
13966a818285SBarry Smith   PetscFunctionBegin;
1397cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
1398268a049cSStefano Zampini   if (mapping->info_free) {
13999566063dSJacob Faibussowitsch     PetscCall(PetscFree(*numprocs));
14006a818285SBarry Smith     if (*indices) {
1401268a049cSStefano Zampini       PetscInt i;
1402268a049cSStefano Zampini 
14039566063dSJacob Faibussowitsch       PetscCall(PetscFree((*indices)[0]));
1404*48a46eb9SPierre Jolivet       for (i = 1; i < *nproc; i++) PetscCall(PetscFree((*indices)[i]));
14059566063dSJacob Faibussowitsch       PetscCall(PetscFree(*indices));
14066a818285SBarry Smith     }
1407268a049cSStefano Zampini   }
1408268a049cSStefano Zampini   *nproc    = 0;
1409268a049cSStefano Zampini   *procs    = NULL;
1410268a049cSStefano Zampini   *numprocs = NULL;
1411268a049cSStefano Zampini   *indices  = NULL;
14126a818285SBarry Smith   PetscFunctionReturn(0);
14136a818285SBarry Smith }
14146a818285SBarry Smith 
14156a818285SBarry Smith /*@C
14166a818285SBarry Smith     ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each processor and
14176a818285SBarry Smith      each index shared by more than one processor
14186a818285SBarry Smith 
14196a818285SBarry Smith     Collective on ISLocalToGlobalMapping
14206a818285SBarry Smith 
1421f899ff85SJose E. Roman     Input Parameter:
14226a818285SBarry Smith .   mapping - the mapping from local to global indexing
14236a818285SBarry Smith 
1424d8d19677SJose E. Roman     Output Parameters:
14256a818285SBarry Smith +   nproc - number of processors that are connected to this one
14266a818285SBarry Smith .   proc - neighboring processors
14276a818285SBarry Smith .   numproc - number of indices for each subdomain (processor)
14286a818285SBarry Smith -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
14296a818285SBarry Smith 
14306a818285SBarry Smith     Level: advanced
14316a818285SBarry Smith 
14321bd0b88eSStefano Zampini     Notes: The user needs to call ISLocalToGlobalMappingRestoreInfo when the data is no longer needed.
14331bd0b88eSStefano Zampini 
14346a818285SBarry Smith     Fortran Usage:
14356a818285SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
14366a818285SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
14376a818285SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
14386a818285SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
14396a818285SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
14406a818285SBarry Smith 
1441db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1442db781477SPatrick Sanan           `ISLocalToGlobalMappingRestoreInfo()`
14436a818285SBarry Smith @*/
14449371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[]) {
14458a1f772fSStefano Zampini   PetscInt **bindices = NULL, *bnumprocs = NULL, bs, i, j, k;
14466a818285SBarry Smith 
14476a818285SBarry Smith   PetscFunctionBegin;
14486a818285SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
14498a1f772fSStefano Zampini   bs = mapping->bs;
14509566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockInfo(mapping, nproc, procs, &bnumprocs, &bindices));
1451268a049cSStefano Zampini   if (bs > 1) { /* we need to expand the cached info */
14529566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(*nproc, &*indices));
14539566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(*nproc, &*numprocs));
14546a818285SBarry Smith     for (i = 0; i < *nproc; i++) {
14559566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(bs * bnumprocs[i], &(*indices)[i]));
1456268a049cSStefano Zampini       for (j = 0; j < bnumprocs[i]; j++) {
14579371c9d4SSatish Balay         for (k = 0; k < bs; k++) { (*indices)[i][j * bs + k] = bs * bindices[i][j] + k; }
14586a818285SBarry Smith       }
1459268a049cSStefano Zampini       (*numprocs)[i] = bnumprocs[i] * bs;
14606a818285SBarry Smith     }
1461268a049cSStefano Zampini     mapping->info_free = PETSC_TRUE;
1462268a049cSStefano Zampini   } else {
1463268a049cSStefano Zampini     *numprocs = bnumprocs;
1464268a049cSStefano Zampini     *indices  = bindices;
14656a818285SBarry Smith   }
14666a818285SBarry Smith   PetscFunctionReturn(0);
14676a818285SBarry Smith }
14686a818285SBarry Smith 
146907b52d57SBarry Smith /*@C
147007b52d57SBarry Smith     ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by ISLocalToGlobalMappingGetInfo()
147189d82c54SBarry Smith 
147207b52d57SBarry Smith     Collective on ISLocalToGlobalMapping
147307b52d57SBarry Smith 
1474f899ff85SJose E. Roman     Input Parameter:
147507b52d57SBarry Smith .   mapping - the mapping from local to global indexing
147607b52d57SBarry Smith 
1477d8d19677SJose E. Roman     Output Parameters:
147807b52d57SBarry Smith +   nproc - number of processors that are connected to this one
147907b52d57SBarry Smith .   proc - neighboring processors
148007b52d57SBarry Smith .   numproc - number of indices for each processor
148107b52d57SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
148207b52d57SBarry Smith 
148307b52d57SBarry Smith     Level: advanced
148407b52d57SBarry Smith 
1485db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1486db781477SPatrick Sanan           `ISLocalToGlobalMappingGetInfo()`
148707b52d57SBarry Smith @*/
14889371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[]) {
148907b52d57SBarry Smith   PetscFunctionBegin;
14909566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockInfo(mapping, nproc, procs, numprocs, indices));
149107b52d57SBarry Smith   PetscFunctionReturn(0);
149207b52d57SBarry Smith }
149386994e45SJed Brown 
149486994e45SJed Brown /*@C
14951bd0b88eSStefano Zampini     ISLocalToGlobalMappingGetNodeInfo - Gets the neighbor information for each node
14961bd0b88eSStefano Zampini 
14971bd0b88eSStefano Zampini     Collective on ISLocalToGlobalMapping
14981bd0b88eSStefano Zampini 
1499f899ff85SJose E. Roman     Input Parameter:
15001bd0b88eSStefano Zampini .   mapping - the mapping from local to global indexing
15011bd0b88eSStefano Zampini 
1502d8d19677SJose E. Roman     Output Parameters:
15031bd0b88eSStefano Zampini +   nnodes - number of local nodes (same ISLocalToGlobalMappingGetSize())
15041bd0b88eSStefano Zampini .   count - number of neighboring processors per node
15051bd0b88eSStefano Zampini -   indices - indices of processes sharing the node (sorted)
15061bd0b88eSStefano Zampini 
15071bd0b88eSStefano Zampini     Level: advanced
15081bd0b88eSStefano Zampini 
15091bd0b88eSStefano Zampini     Notes: The user needs to call ISLocalToGlobalMappingRestoreInfo when the data is no longer needed.
15101bd0b88eSStefano Zampini 
1511db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1512db781477SPatrick Sanan           `ISLocalToGlobalMappingGetInfo()`, `ISLocalToGlobalMappingRestoreNodeInfo()`
15131bd0b88eSStefano Zampini @*/
15149371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingGetNodeInfo(ISLocalToGlobalMapping mapping, PetscInt *nnodes, PetscInt *count[], PetscInt **indices[]) {
15151bd0b88eSStefano Zampini   PetscInt n;
15161bd0b88eSStefano Zampini 
15171bd0b88eSStefano Zampini   PetscFunctionBegin;
15181bd0b88eSStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
15199566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetSize(mapping, &n));
15201bd0b88eSStefano Zampini   if (!mapping->info_nodec) {
15211bd0b88eSStefano Zampini     PetscInt i, m, n_neigh, *neigh, *n_shared, **shared;
15221bd0b88eSStefano Zampini 
15239566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(n + 1, &mapping->info_nodec, n, &mapping->info_nodei));
15249566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetInfo(mapping, &n_neigh, &neigh, &n_shared, &shared));
1525071fcb05SBarry Smith     for (i = 0; i < n; i++) { mapping->info_nodec[i] = 1; }
1526071fcb05SBarry Smith     m                      = n;
1527071fcb05SBarry Smith     mapping->info_nodec[n] = 0;
15281bd0b88eSStefano Zampini     for (i = 1; i < n_neigh; i++) {
15291bd0b88eSStefano Zampini       PetscInt j;
15301bd0b88eSStefano Zampini 
15311bd0b88eSStefano Zampini       m += n_shared[i];
15321bd0b88eSStefano Zampini       for (j = 0; j < n_shared[i]; j++) mapping->info_nodec[shared[i][j]] += 1;
15331bd0b88eSStefano Zampini     }
15349566063dSJacob Faibussowitsch     if (n) PetscCall(PetscMalloc1(m, &mapping->info_nodei[0]));
15351bd0b88eSStefano Zampini     for (i = 1; i < n; i++) mapping->info_nodei[i] = mapping->info_nodei[i - 1] + mapping->info_nodec[i - 1];
15369566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(mapping->info_nodec, n));
15379371c9d4SSatish Balay     for (i = 0; i < n; i++) {
15389371c9d4SSatish Balay       mapping->info_nodec[i]    = 1;
15399371c9d4SSatish Balay       mapping->info_nodei[i][0] = neigh[0];
15409371c9d4SSatish Balay     }
15411bd0b88eSStefano Zampini     for (i = 1; i < n_neigh; i++) {
15421bd0b88eSStefano Zampini       PetscInt j;
15431bd0b88eSStefano Zampini 
15441bd0b88eSStefano Zampini       for (j = 0; j < n_shared[i]; j++) {
15451bd0b88eSStefano Zampini         PetscInt k = shared[i][j];
15461bd0b88eSStefano Zampini 
15471bd0b88eSStefano Zampini         mapping->info_nodei[k][mapping->info_nodec[k]] = neigh[i];
15481bd0b88eSStefano Zampini         mapping->info_nodec[k] += 1;
15491bd0b88eSStefano Zampini       }
15501bd0b88eSStefano Zampini     }
15519566063dSJacob Faibussowitsch     for (i = 0; i < n; i++) PetscCall(PetscSortRemoveDupsInt(&mapping->info_nodec[i], mapping->info_nodei[i]));
15529566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreInfo(mapping, &n_neigh, &neigh, &n_shared, &shared));
15531bd0b88eSStefano Zampini   }
15541bd0b88eSStefano Zampini   if (nnodes) *nnodes = n;
15551bd0b88eSStefano Zampini   if (count) *count = mapping->info_nodec;
15561bd0b88eSStefano Zampini   if (indices) *indices = mapping->info_nodei;
15571bd0b88eSStefano Zampini   PetscFunctionReturn(0);
15581bd0b88eSStefano Zampini }
15591bd0b88eSStefano Zampini 
15601bd0b88eSStefano Zampini /*@C
15611bd0b88eSStefano Zampini     ISLocalToGlobalMappingRestoreNodeInfo - Frees the memory allocated by ISLocalToGlobalMappingGetNodeInfo()
15621bd0b88eSStefano Zampini 
15631bd0b88eSStefano Zampini     Collective on ISLocalToGlobalMapping
15641bd0b88eSStefano Zampini 
1565f899ff85SJose E. Roman     Input Parameter:
15661bd0b88eSStefano Zampini .   mapping - the mapping from local to global indexing
15671bd0b88eSStefano Zampini 
1568d8d19677SJose E. Roman     Output Parameters:
15691bd0b88eSStefano Zampini +   nnodes - number of local nodes
15701bd0b88eSStefano Zampini .   count - number of neighboring processors per node
15711bd0b88eSStefano Zampini -   indices - indices of processes sharing the node (sorted)
15721bd0b88eSStefano Zampini 
15731bd0b88eSStefano Zampini     Level: advanced
15741bd0b88eSStefano Zampini 
1575db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1576db781477SPatrick Sanan           `ISLocalToGlobalMappingGetInfo()`
15771bd0b88eSStefano Zampini @*/
15789371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingRestoreNodeInfo(ISLocalToGlobalMapping mapping, PetscInt *nnodes, PetscInt *count[], PetscInt **indices[]) {
15791bd0b88eSStefano Zampini   PetscFunctionBegin;
15801bd0b88eSStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
15811bd0b88eSStefano Zampini   if (nnodes) *nnodes = 0;
15821bd0b88eSStefano Zampini   if (count) *count = NULL;
15831bd0b88eSStefano Zampini   if (indices) *indices = NULL;
15841bd0b88eSStefano Zampini   PetscFunctionReturn(0);
15851bd0b88eSStefano Zampini }
15861bd0b88eSStefano Zampini 
15871bd0b88eSStefano Zampini /*@C
1588107e9a97SBarry Smith    ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped
158986994e45SJed Brown 
159086994e45SJed Brown    Not Collective
159186994e45SJed Brown 
15924165533cSJose E. Roman    Input Parameter:
159386994e45SJed Brown . ltog - local to global mapping
159486994e45SJed Brown 
15954165533cSJose E. Roman    Output Parameter:
1596565245c5SBarry Smith . array - array of indices, the length of this array may be obtained with ISLocalToGlobalMappingGetSize()
159786994e45SJed Brown 
159886994e45SJed Brown    Level: advanced
159986994e45SJed Brown 
160095452b02SPatrick Sanan    Notes:
160195452b02SPatrick Sanan     ISLocalToGlobalMappingGetSize() returns the length the this array
1602107e9a97SBarry Smith 
1603db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingRestoreIndices()`, `ISLocalToGlobalMappingGetBlockIndices()`, `ISLocalToGlobalMappingRestoreBlockIndices()`
160486994e45SJed Brown @*/
16059371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog, const PetscInt **array) {
160686994e45SJed Brown   PetscFunctionBegin;
160786994e45SJed Brown   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
160886994e45SJed Brown   PetscValidPointer(array, 2);
160945b6f7e9SBarry Smith   if (ltog->bs == 1) {
161086994e45SJed Brown     *array = ltog->indices;
161145b6f7e9SBarry Smith   } else {
161245b6f7e9SBarry Smith     PetscInt       *jj, k, i, j, n = ltog->n, bs = ltog->bs;
161345b6f7e9SBarry Smith     const PetscInt *ii;
161445b6f7e9SBarry Smith 
16159566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(bs * n, &jj));
161645b6f7e9SBarry Smith     *array = jj;
161745b6f7e9SBarry Smith     k      = 0;
161845b6f7e9SBarry Smith     ii     = ltog->indices;
161945b6f7e9SBarry Smith     for (i = 0; i < n; i++)
16209371c9d4SSatish Balay       for (j = 0; j < bs; j++) jj[k++] = bs * ii[i] + j;
162145b6f7e9SBarry Smith   }
162286994e45SJed Brown   PetscFunctionReturn(0);
162386994e45SJed Brown }
162486994e45SJed Brown 
162586994e45SJed Brown /*@C
1626193a2b41SJulian Andrej    ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with ISLocalToGlobalMappingGetIndices()
162786994e45SJed Brown 
162886994e45SJed Brown    Not Collective
162986994e45SJed Brown 
16304165533cSJose E. Roman    Input Parameters:
163186994e45SJed Brown + ltog - local to global mapping
163286994e45SJed Brown - array - array of indices
163386994e45SJed Brown 
163486994e45SJed Brown    Level: advanced
163586994e45SJed Brown 
1636db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingGetIndices()`
163786994e45SJed Brown @*/
16389371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog, const PetscInt **array) {
163986994e45SJed Brown   PetscFunctionBegin;
164086994e45SJed Brown   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
164186994e45SJed Brown   PetscValidPointer(array, 2);
1642c9cc58a2SBarry Smith   PetscCheck(ltog->bs != 1 || *array == ltog->indices, PETSC_COMM_SELF, PETSC_ERR_ARG_BADPTR, "Trying to return mismatched pointer");
164345b6f7e9SBarry Smith 
1644*48a46eb9SPierre Jolivet   if (ltog->bs > 1) PetscCall(PetscFree(*(void **)array));
164545b6f7e9SBarry Smith   PetscFunctionReturn(0);
164645b6f7e9SBarry Smith }
164745b6f7e9SBarry Smith 
164845b6f7e9SBarry Smith /*@C
164945b6f7e9SBarry Smith    ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block
165045b6f7e9SBarry Smith 
165145b6f7e9SBarry Smith    Not Collective
165245b6f7e9SBarry Smith 
16534165533cSJose E. Roman    Input Parameter:
165445b6f7e9SBarry Smith . ltog - local to global mapping
165545b6f7e9SBarry Smith 
16564165533cSJose E. Roman    Output Parameter:
165745b6f7e9SBarry Smith . array - array of indices
165845b6f7e9SBarry Smith 
165945b6f7e9SBarry Smith    Level: advanced
166045b6f7e9SBarry Smith 
1661db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingRestoreBlockIndices()`
166245b6f7e9SBarry Smith @*/
16639371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog, const PetscInt **array) {
166445b6f7e9SBarry Smith   PetscFunctionBegin;
166545b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
166645b6f7e9SBarry Smith   PetscValidPointer(array, 2);
166745b6f7e9SBarry Smith   *array = ltog->indices;
166845b6f7e9SBarry Smith   PetscFunctionReturn(0);
166945b6f7e9SBarry Smith }
167045b6f7e9SBarry Smith 
167145b6f7e9SBarry Smith /*@C
167245b6f7e9SBarry Smith    ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with ISLocalToGlobalMappingGetBlockIndices()
167345b6f7e9SBarry Smith 
167445b6f7e9SBarry Smith    Not Collective
167545b6f7e9SBarry Smith 
16764165533cSJose E. Roman    Input Parameters:
167745b6f7e9SBarry Smith + ltog - local to global mapping
167845b6f7e9SBarry Smith - array - array of indices
167945b6f7e9SBarry Smith 
168045b6f7e9SBarry Smith    Level: advanced
168145b6f7e9SBarry Smith 
1682db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingGetIndices()`
168345b6f7e9SBarry Smith @*/
16849371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog, const PetscInt **array) {
168545b6f7e9SBarry Smith   PetscFunctionBegin;
168645b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
168745b6f7e9SBarry Smith   PetscValidPointer(array, 2);
168808401ef6SPierre Jolivet   PetscCheck(*array == ltog->indices, PETSC_COMM_SELF, PETSC_ERR_ARG_BADPTR, "Trying to return mismatched pointer");
16890298fd71SBarry Smith   *array = NULL;
169086994e45SJed Brown   PetscFunctionReturn(0);
169186994e45SJed Brown }
1692f7efa3c7SJed Brown 
1693f7efa3c7SJed Brown /*@C
1694f7efa3c7SJed Brown    ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings
1695f7efa3c7SJed Brown 
1696f7efa3c7SJed Brown    Not Collective
1697f7efa3c7SJed Brown 
16984165533cSJose E. Roman    Input Parameters:
1699f7efa3c7SJed Brown + comm - communicator for the new mapping, must contain the communicator of every mapping to concatenate
1700f7efa3c7SJed Brown . n - number of mappings to concatenate
1701f7efa3c7SJed Brown - ltogs - local to global mappings
1702f7efa3c7SJed Brown 
17034165533cSJose E. Roman    Output Parameter:
1704f7efa3c7SJed Brown . ltogcat - new mapping
1705f7efa3c7SJed Brown 
17069d90f715SBarry Smith    Note: this currently always returns a mapping with block size of 1
17079d90f715SBarry Smith 
17089d90f715SBarry Smith    Developer Note: If all the input mapping have the same block size we could easily handle that as a special case
17099d90f715SBarry Smith 
1710f7efa3c7SJed Brown    Level: advanced
1711f7efa3c7SJed Brown 
1712db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`
1713f7efa3c7SJed Brown @*/
17149371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm, PetscInt n, const ISLocalToGlobalMapping ltogs[], ISLocalToGlobalMapping *ltogcat) {
1715f7efa3c7SJed Brown   PetscInt i, cnt, m, *idx;
1716f7efa3c7SJed Brown 
1717f7efa3c7SJed Brown   PetscFunctionBegin;
171808401ef6SPierre Jolivet   PetscCheck(n >= 0, comm, PETSC_ERR_ARG_OUTOFRANGE, "Must have a non-negative number of mappings, given %" PetscInt_FMT, n);
1719f7efa3c7SJed Brown   if (n > 0) PetscValidPointer(ltogs, 3);
1720f7efa3c7SJed Brown   for (i = 0; i < n; i++) PetscValidHeaderSpecific(ltogs[i], IS_LTOGM_CLASSID, 3);
1721f7efa3c7SJed Brown   PetscValidPointer(ltogcat, 4);
1722f7efa3c7SJed Brown   for (cnt = 0, i = 0; i < n; i++) {
17239566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(ltogs[i], &m));
1724f7efa3c7SJed Brown     cnt += m;
1725f7efa3c7SJed Brown   }
17269566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(cnt, &idx));
1727f7efa3c7SJed Brown   for (cnt = 0, i = 0; i < n; i++) {
1728f7efa3c7SJed Brown     const PetscInt *subidx;
17299566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(ltogs[i], &m));
17309566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetIndices(ltogs[i], &subidx));
17319566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(&idx[cnt], subidx, m));
17329566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreIndices(ltogs[i], &subidx));
1733f7efa3c7SJed Brown     cnt += m;
1734f7efa3c7SJed Brown   }
17359566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(comm, 1, cnt, idx, PETSC_OWN_POINTER, ltogcat));
1736f7efa3c7SJed Brown   PetscFunctionReturn(0);
1737f7efa3c7SJed Brown }
173804a59952SBarry Smith 
1739413f72f0SBarry Smith /*MC
1740413f72f0SBarry Smith       ISLOCALTOGLOBALMAPPINGBASIC - basic implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is
1741413f72f0SBarry Smith                                     used this is good for only small and moderate size problems.
1742413f72f0SBarry Smith 
1743413f72f0SBarry Smith    Options Database Keys:
1744a2b725a8SWilliam Gropp .   -islocaltoglobalmapping_type basic - select this method
1745413f72f0SBarry Smith 
1746413f72f0SBarry Smith    Level: beginner
1747413f72f0SBarry Smith 
1748db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`, `ISLOCALTOGLOBALMAPPINGHASH`
1749413f72f0SBarry Smith M*/
17509371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Basic(ISLocalToGlobalMapping ltog) {
1751413f72f0SBarry Smith   PetscFunctionBegin;
1752413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Basic;
1753413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Basic;
1754413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Basic;
1755413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Basic;
1756413f72f0SBarry Smith   PetscFunctionReturn(0);
1757413f72f0SBarry Smith }
1758413f72f0SBarry Smith 
1759413f72f0SBarry Smith /*MC
1760413f72f0SBarry Smith       ISLOCALTOGLOBALMAPPINGHASH - hash implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is
1761ed56e8eaSBarry Smith                                     used this is good for large memory problems.
1762413f72f0SBarry Smith 
1763413f72f0SBarry Smith    Options Database Keys:
1764a2b725a8SWilliam Gropp .   -islocaltoglobalmapping_type hash - select this method
1765413f72f0SBarry Smith 
176695452b02SPatrick Sanan    Notes:
176795452b02SPatrick Sanan     This is selected automatically for large problems if the user does not set the type.
1768ed56e8eaSBarry Smith 
1769413f72f0SBarry Smith    Level: beginner
1770413f72f0SBarry Smith 
1771db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`, `ISLOCALTOGLOBALMAPPINGHASH`
1772413f72f0SBarry Smith M*/
17739371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Hash(ISLocalToGlobalMapping ltog) {
1774413f72f0SBarry Smith   PetscFunctionBegin;
1775413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Hash;
1776413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Hash;
1777413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Hash;
1778413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Hash;
1779413f72f0SBarry Smith   PetscFunctionReturn(0);
1780413f72f0SBarry Smith }
1781413f72f0SBarry Smith 
1782413f72f0SBarry Smith /*@C
1783413f72f0SBarry Smith     ISLocalToGlobalMappingRegister -  Adds a method for applying a global to local mapping with an ISLocalToGlobalMapping
1784413f72f0SBarry Smith 
1785413f72f0SBarry Smith    Not Collective
1786413f72f0SBarry Smith 
1787413f72f0SBarry Smith    Input Parameters:
1788413f72f0SBarry Smith +  sname - name of a new method
1789413f72f0SBarry Smith -  routine_create - routine to create method context
1790413f72f0SBarry Smith 
1791413f72f0SBarry Smith    Notes:
1792ed56e8eaSBarry Smith    ISLocalToGlobalMappingRegister() may be called multiple times to add several user-defined mappings.
1793413f72f0SBarry Smith 
1794413f72f0SBarry Smith    Sample usage:
1795413f72f0SBarry Smith .vb
1796413f72f0SBarry Smith    ISLocalToGlobalMappingRegister("my_mapper",MyCreate);
1797413f72f0SBarry Smith .ve
1798413f72f0SBarry Smith 
1799ed56e8eaSBarry Smith    Then, your mapping can be chosen with the procedural interface via
1800413f72f0SBarry Smith $     ISLocalToGlobalMappingSetType(ltog,"my_mapper")
1801413f72f0SBarry Smith    or at runtime via the option
1802ed56e8eaSBarry Smith $     -islocaltoglobalmapping_type my_mapper
1803413f72f0SBarry Smith 
1804413f72f0SBarry Smith    Level: advanced
1805413f72f0SBarry Smith 
1806db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingRegisterAll()`, `ISLocalToGlobalMappingRegisterDestroy()`, `ISLOCALTOGLOBALMAPPINGBASIC`, `ISLOCALTOGLOBALMAPPINGHASH`
1807413f72f0SBarry Smith 
1808413f72f0SBarry Smith @*/
18099371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingRegister(const char sname[], PetscErrorCode (*function)(ISLocalToGlobalMapping)) {
1810413f72f0SBarry Smith   PetscFunctionBegin;
18119566063dSJacob Faibussowitsch   PetscCall(ISInitializePackage());
18129566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&ISLocalToGlobalMappingList, sname, function));
1813413f72f0SBarry Smith   PetscFunctionReturn(0);
1814413f72f0SBarry Smith }
1815413f72f0SBarry Smith 
1816413f72f0SBarry Smith /*@C
1817ed56e8eaSBarry Smith    ISLocalToGlobalMappingSetType - Builds ISLocalToGlobalMapping for a particular global to local mapping approach.
1818413f72f0SBarry Smith 
1819413f72f0SBarry Smith    Logically Collective on ISLocalToGlobalMapping
1820413f72f0SBarry Smith 
1821413f72f0SBarry Smith    Input Parameters:
1822413f72f0SBarry Smith +  ltog - the ISLocalToGlobalMapping object
1823413f72f0SBarry Smith -  type - a known method
1824413f72f0SBarry Smith 
1825413f72f0SBarry Smith    Options Database Key:
1826ed56e8eaSBarry Smith .  -islocaltoglobalmapping_type  <method> - Sets the method; use -help for a list
1827413f72f0SBarry Smith     of available methods (for instance, basic or hash)
1828413f72f0SBarry Smith 
1829413f72f0SBarry Smith    Notes:
1830413f72f0SBarry Smith    See "petsc/include/petscis.h" for available methods
1831413f72f0SBarry Smith 
1832413f72f0SBarry Smith   Normally, it is best to use the ISLocalToGlobalMappingSetFromOptions() command and
1833413f72f0SBarry Smith   then set the ISLocalToGlobalMapping type from the options database rather than by using
1834413f72f0SBarry Smith   this routine.
1835413f72f0SBarry Smith 
1836413f72f0SBarry Smith   Level: intermediate
1837413f72f0SBarry Smith 
1838413f72f0SBarry Smith   Developer Note: ISLocalToGlobalMappingRegister() is used to add new types to ISLocalToGlobalMappingList from which they
1839413f72f0SBarry Smith   are accessed by ISLocalToGlobalMappingSetType().
1840413f72f0SBarry Smith 
1841db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingRegister()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingGetType()`
1842413f72f0SBarry Smith @*/
18439371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType type) {
1844413f72f0SBarry Smith   PetscBool match;
18455f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(ISLocalToGlobalMapping) = NULL;
1846413f72f0SBarry Smith 
1847413f72f0SBarry Smith   PetscFunctionBegin;
1848413f72f0SBarry Smith   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
1849a0d79125SStefano Zampini   if (type) PetscValidCharPointer(type, 2);
1850413f72f0SBarry Smith 
18519566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ltog, type, &match));
1852413f72f0SBarry Smith   if (match) PetscFunctionReturn(0);
1853413f72f0SBarry Smith 
1854a0d79125SStefano Zampini   /* L2G maps defer type setup at globaltolocal calls, allow passing NULL here */
1855a0d79125SStefano Zampini   if (type) {
18569566063dSJacob Faibussowitsch     PetscCall(PetscFunctionListFind(ISLocalToGlobalMappingList, type, &r));
1857a0d79125SStefano Zampini     PetscCheck(r, PetscObjectComm((PetscObject)ltog), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested ISLocalToGlobalMapping type %s", type);
1858a0d79125SStefano Zampini   }
1859413f72f0SBarry Smith   /* Destroy the previous private LTOG context */
1860dbbe0bcdSBarry Smith   PetscTryTypeMethod(ltog, destroy);
1861413f72f0SBarry Smith   ltog->ops->destroy = NULL;
1862dbbe0bcdSBarry Smith 
18639566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)ltog, type));
18649566063dSJacob Faibussowitsch   if (r) PetscCall((*r)(ltog));
1865a0d79125SStefano Zampini   PetscFunctionReturn(0);
1866a0d79125SStefano Zampini }
1867a0d79125SStefano Zampini 
1868a0d79125SStefano Zampini /*@C
1869a0d79125SStefano Zampini    ISLocalToGlobalMappingGetType - Get the type of the l2g map
1870a0d79125SStefano Zampini 
1871a0d79125SStefano Zampini    Not Collective
1872a0d79125SStefano Zampini 
1873a0d79125SStefano Zampini    Input Parameter:
1874a0d79125SStefano Zampini .  ltog - the ISLocalToGlobalMapping object
1875a0d79125SStefano Zampini 
1876a0d79125SStefano Zampini    Output Parameter:
1877a0d79125SStefano Zampini .  type - the type
1878a0d79125SStefano Zampini 
187949762cbcSSatish Balay    Level: intermediate
188049762cbcSSatish Balay 
1881db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingRegister()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`
1882a0d79125SStefano Zampini @*/
18839371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingGetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType *type) {
1884a0d79125SStefano Zampini   PetscFunctionBegin;
1885a0d79125SStefano Zampini   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
1886a0d79125SStefano Zampini   PetscValidPointer(type, 2);
1887a0d79125SStefano Zampini   *type = ((PetscObject)ltog)->type_name;
1888413f72f0SBarry Smith   PetscFunctionReturn(0);
1889413f72f0SBarry Smith }
1890413f72f0SBarry Smith 
1891413f72f0SBarry Smith PetscBool ISLocalToGlobalMappingRegisterAllCalled = PETSC_FALSE;
1892413f72f0SBarry Smith 
1893413f72f0SBarry Smith /*@C
1894413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAll - Registers all of the local to global mapping components in the IS package.
1895413f72f0SBarry Smith 
1896413f72f0SBarry Smith   Not Collective
1897413f72f0SBarry Smith 
1898413f72f0SBarry Smith   Level: advanced
1899413f72f0SBarry Smith 
1900db781477SPatrick Sanan .seealso: `ISRegister()`, `ISLocalToGlobalRegister()`
1901413f72f0SBarry Smith @*/
19029371c9d4SSatish Balay PetscErrorCode ISLocalToGlobalMappingRegisterAll(void) {
1903413f72f0SBarry Smith   PetscFunctionBegin;
1904413f72f0SBarry Smith   if (ISLocalToGlobalMappingRegisterAllCalled) PetscFunctionReturn(0);
1905413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAllCalled = PETSC_TRUE;
19069566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGBASIC, ISLocalToGlobalMappingCreate_Basic));
19079566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGHASH, ISLocalToGlobalMappingCreate_Hash));
1908413f72f0SBarry Smith   PetscFunctionReturn(0);
1909413f72f0SBarry Smith }
1910