xref: /petsc/src/vec/is/utils/isltog.c (revision 2fe279fdf3e687a416e4eadb7d3c7a82d60442c6)
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
19cab54364SBarry Smith   ISGetPointRange - Returns a description of the points in an `IS` suitable for traversal
20413f72f0SBarry Smith 
2120662ed9SBarry Smith   Not Collective
226528b96dSMatthew G. Knepley 
236528b96dSMatthew G. Knepley   Input Parameter:
24cab54364SBarry Smith . 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   Level: intermediate
326528b96dSMatthew G. Knepley 
33cab54364SBarry Smith   Notes:
3420662ed9SBarry Smith   If the `IS` contains contiguous indices in an `ISSTRIDE`, then the indices are contained in [pStart, pEnd) and points = `NULL`.
3520662ed9SBarry Smith   Otherwise, `pStart = 0`, `pEnd = numIndices`, and points is an array of the indices. This supports the following pattern
36cab54364SBarry Smith .vb
37cab54364SBarry Smith   ISGetPointRange(is, &pStart, &pEnd, &points);
38cab54364SBarry Smith   for (p = pStart; p < pEnd; ++p) {
39cab54364SBarry Smith     const PetscInt point = points ? points[p] : p;
4020662ed9SBarry Smith     // use point
41cab54364SBarry Smith   }
42cab54364SBarry Smith   ISRestorePointRange(is, &pstart, &pEnd, &points);
43cab54364SBarry Smith .ve
4420662ed9SBarry Smith   Hence the same code can be written for `pointIS` being a `ISSTRIDE` or `ISGENERAL`
45cab54364SBarry Smith 
46cab54364SBarry Smith .seealso: [](sec_scatter), `IS`, `ISRestorePointRange()`, `ISGetPointSubrange()`, `ISGetIndices()`, `ISCreateStride()`
476528b96dSMatthew G. Knepley @*/
48d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGetPointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points)
49d71ae5a4SJacob Faibussowitsch {
509305a4c7SMatthew G. Knepley   PetscInt  numCells, step = 1;
519305a4c7SMatthew G. Knepley   PetscBool isStride;
529305a4c7SMatthew G. Knepley 
539305a4c7SMatthew G. Knepley   PetscFunctionBeginHot;
549305a4c7SMatthew G. Knepley   *pStart = 0;
559305a4c7SMatthew G. Knepley   *points = NULL;
569566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(pointIS, &numCells));
579566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pointIS, ISSTRIDE, &isStride));
589566063dSJacob Faibussowitsch   if (isStride) PetscCall(ISStrideGetInfo(pointIS, pStart, &step));
599305a4c7SMatthew G. Knepley   *pEnd = *pStart + numCells;
609566063dSJacob Faibussowitsch   if (!isStride || step != 1) PetscCall(ISGetIndices(pointIS, points));
613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
629305a4c7SMatthew G. Knepley }
639305a4c7SMatthew G. Knepley 
646528b96dSMatthew G. Knepley /*@C
65cab54364SBarry Smith   ISRestorePointRange - Destroys the traversal description created with `ISGetPointRange()`
666528b96dSMatthew G. Knepley 
6720662ed9SBarry Smith   Not Collective
686528b96dSMatthew G. Knepley 
696528b96dSMatthew G. Knepley   Input Parameters:
70cab54364SBarry Smith + pointIS - The `IS` object
71cab54364SBarry Smith . pStart  - The first index, from `ISGetPointRange()`
72cab54364SBarry Smith . pEnd    - One past the last index, from `ISGetPointRange()`
73cab54364SBarry Smith - points  - The indices, from `ISGetPointRange()`
746528b96dSMatthew G. Knepley 
756528b96dSMatthew G. Knepley   Level: intermediate
766528b96dSMatthew G. Knepley 
77cab54364SBarry Smith .seealso: [](sec_scatter), `IS`, `ISGetPointRange()`, `ISGetPointSubrange()`, `ISGetIndices()`, `ISCreateStride()`
786528b96dSMatthew G. Knepley @*/
79d71ae5a4SJacob Faibussowitsch PetscErrorCode ISRestorePointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points)
80d71ae5a4SJacob Faibussowitsch {
819305a4c7SMatthew G. Knepley   PetscInt  step = 1;
829305a4c7SMatthew G. Knepley   PetscBool isStride;
839305a4c7SMatthew G. Knepley 
849305a4c7SMatthew G. Knepley   PetscFunctionBeginHot;
859566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pointIS, ISSTRIDE, &isStride));
869566063dSJacob Faibussowitsch   if (isStride) PetscCall(ISStrideGetInfo(pointIS, pStart, &step));
879566063dSJacob Faibussowitsch   if (!isStride || step != 1) PetscCall(ISGetIndices(pointIS, points));
883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
899305a4c7SMatthew G. Knepley }
909305a4c7SMatthew G. Knepley 
916528b96dSMatthew G. Knepley /*@C
92cab54364SBarry Smith   ISGetPointSubrange - Configures the input `IS` to be a subrange for the traversal information given
936528b96dSMatthew G. Knepley 
9420662ed9SBarry Smith   Not Collective
956528b96dSMatthew G. Knepley 
966528b96dSMatthew G. Knepley   Input Parameters:
97cab54364SBarry Smith + subpointIS - The `IS` object to be configured
986528b96dSMatthew G. Knepley . pStart     - The first index of the subrange
996528b96dSMatthew G. Knepley . pEnd       - One past the last index for the subrange
100cab54364SBarry Smith - points     - The indices for the entire range, from `ISGetPointRange()`
1016528b96dSMatthew G. Knepley 
1026528b96dSMatthew G. Knepley   Output Parameters:
103cab54364SBarry Smith . subpointIS - The `IS` object now configured to be a subrange
1046528b96dSMatthew G. Knepley 
1056528b96dSMatthew G. Knepley   Level: intermediate
1066528b96dSMatthew G. Knepley 
107cab54364SBarry Smith   Note:
108cab54364SBarry Smith   The input `IS` will now respond properly to calls to `ISGetPointRange()` and return the subrange.
109cab54364SBarry Smith 
110cab54364SBarry Smith .seealso: [](sec_scatter), `IS`, `ISGetPointRange()`, `ISRestorePointRange()`, `ISGetIndices()`, `ISCreateStride()`
1116528b96dSMatthew G. Knepley @*/
112d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGetPointSubrange(IS subpointIS, PetscInt pStart, PetscInt pEnd, const PetscInt *points)
113d71ae5a4SJacob Faibussowitsch {
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   }
1223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 */
132d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISGlobalToLocalMappingSetUp(ISLocalToGlobalMapping mapping)
133d71ae5a4SJacob Faibussowitsch {
134413f72f0SBarry Smith   PetscInt i, *idx = mapping->indices, n = mapping->n, end, start;
135413f72f0SBarry Smith 
136413f72f0SBarry Smith   PetscFunctionBegin;
1373ba16761SJacob Faibussowitsch   if (mapping->data) PetscFunctionReturn(PETSC_SUCCESS);
138413f72f0SBarry Smith   end   = 0;
139413f72f0SBarry Smith   start = PETSC_MAX_INT;
140413f72f0SBarry Smith 
141413f72f0SBarry Smith   for (i = 0; i < n; i++) {
142413f72f0SBarry Smith     if (idx[i] < 0) continue;
143413f72f0SBarry Smith     if (idx[i] < start) start = idx[i];
144413f72f0SBarry Smith     if (idx[i] > end) end = idx[i];
145413f72f0SBarry Smith   }
1469371c9d4SSatish Balay   if (start > end) {
1479371c9d4SSatish Balay     start = 0;
1489371c9d4SSatish Balay     end   = -1;
1499371c9d4SSatish Balay   }
150413f72f0SBarry Smith   mapping->globalstart = start;
151413f72f0SBarry Smith   mapping->globalend   = end;
152413f72f0SBarry Smith   if (!((PetscObject)mapping)->type_name) {
153413f72f0SBarry Smith     if ((end - start) > PetscMax(4 * n, 1000000)) {
1549566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingSetType(mapping, ISLOCALTOGLOBALMAPPINGHASH));
155413f72f0SBarry Smith     } else {
1569566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingSetType(mapping, ISLOCALTOGLOBALMAPPINGBASIC));
157413f72f0SBarry Smith     }
158413f72f0SBarry Smith   }
159dbbe0bcdSBarry Smith   PetscTryTypeMethod(mapping, globaltolocalmappingsetup);
1603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
161413f72f0SBarry Smith }
162413f72f0SBarry Smith 
163d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISGlobalToLocalMappingSetUp_Basic(ISLocalToGlobalMapping mapping)
164d71ae5a4SJacob Faibussowitsch {
165413f72f0SBarry Smith   PetscInt                      i, *idx = mapping->indices, n = mapping->n, end, start, *globals;
166413f72f0SBarry Smith   ISLocalToGlobalMapping_Basic *map;
167413f72f0SBarry Smith 
168413f72f0SBarry Smith   PetscFunctionBegin;
169413f72f0SBarry Smith   start = mapping->globalstart;
170413f72f0SBarry Smith   end   = mapping->globalend;
1719566063dSJacob Faibussowitsch   PetscCall(PetscNew(&map));
1729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(end - start + 2, &globals));
173413f72f0SBarry Smith   map->globals = globals;
174413f72f0SBarry Smith   for (i = 0; i < end - start + 1; i++) globals[i] = -1;
175413f72f0SBarry Smith   for (i = 0; i < n; i++) {
176413f72f0SBarry Smith     if (idx[i] < 0) continue;
177413f72f0SBarry Smith     globals[idx[i] - start] = i;
178413f72f0SBarry Smith   }
179413f72f0SBarry Smith   mapping->data = (void *)map;
1803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
181413f72f0SBarry Smith }
182413f72f0SBarry Smith 
183d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISGlobalToLocalMappingSetUp_Hash(ISLocalToGlobalMapping mapping)
184d71ae5a4SJacob Faibussowitsch {
185413f72f0SBarry Smith   PetscInt                     i, *idx = mapping->indices, n = mapping->n;
186413f72f0SBarry Smith   ISLocalToGlobalMapping_Hash *map;
187413f72f0SBarry Smith 
188413f72f0SBarry Smith   PetscFunctionBegin;
1899566063dSJacob Faibussowitsch   PetscCall(PetscNew(&map));
1909566063dSJacob Faibussowitsch   PetscCall(PetscHMapICreate(&map->globalht));
191413f72f0SBarry Smith   for (i = 0; i < n; i++) {
192413f72f0SBarry Smith     if (idx[i] < 0) continue;
1939566063dSJacob Faibussowitsch     PetscCall(PetscHMapISet(map->globalht, idx[i], i));
194413f72f0SBarry Smith   }
195413f72f0SBarry Smith   mapping->data = (void *)map;
1963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
197413f72f0SBarry Smith }
198413f72f0SBarry Smith 
199d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISLocalToGlobalMappingDestroy_Basic(ISLocalToGlobalMapping mapping)
200d71ae5a4SJacob Faibussowitsch {
201413f72f0SBarry Smith   ISLocalToGlobalMapping_Basic *map = (ISLocalToGlobalMapping_Basic *)mapping->data;
202413f72f0SBarry Smith 
203413f72f0SBarry Smith   PetscFunctionBegin;
2043ba16761SJacob Faibussowitsch   if (!map) PetscFunctionReturn(PETSC_SUCCESS);
2059566063dSJacob Faibussowitsch   PetscCall(PetscFree(map->globals));
2069566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping->data));
2073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
208413f72f0SBarry Smith }
209413f72f0SBarry Smith 
210d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISLocalToGlobalMappingDestroy_Hash(ISLocalToGlobalMapping mapping)
211d71ae5a4SJacob Faibussowitsch {
212413f72f0SBarry Smith   ISLocalToGlobalMapping_Hash *map = (ISLocalToGlobalMapping_Hash *)mapping->data;
213413f72f0SBarry Smith 
214413f72f0SBarry Smith   PetscFunctionBegin;
2153ba16761SJacob Faibussowitsch   if (!map) PetscFunctionReturn(PETSC_SUCCESS);
2169566063dSJacob Faibussowitsch   PetscCall(PetscHMapIDestroy(&map->globalht));
2179566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping->data));
2183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
219413f72f0SBarry Smith }
220413f72f0SBarry Smith 
221413f72f0SBarry Smith #define GTOLTYPE _Basic
222413f72f0SBarry Smith #define GTOLNAME _Basic
223541bf97eSAdrian Croucher #define GTOLBS   mapping->bs
2249371c9d4SSatish Balay #define GTOL(g, local) \
2259371c9d4SSatish Balay   do { \
226413f72f0SBarry Smith     local = map->globals[g / bs - start]; \
2270040bde1SJunchao Zhang     if (local >= 0) local = bs * local + (g % bs); \
228413f72f0SBarry Smith   } while (0)
229541bf97eSAdrian Croucher 
230413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
231413f72f0SBarry Smith 
232413f72f0SBarry Smith #define GTOLTYPE _Basic
233413f72f0SBarry Smith #define GTOLNAME Block_Basic
234541bf97eSAdrian Croucher #define GTOLBS   1
2359371c9d4SSatish Balay #define GTOL(g, local) \
236d71ae5a4SJacob Faibussowitsch   do { \
237d71ae5a4SJacob Faibussowitsch     local = map->globals[g - start]; \
238d71ae5a4SJacob Faibussowitsch   } while (0)
239413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
240413f72f0SBarry Smith 
241413f72f0SBarry Smith #define GTOLTYPE _Hash
242413f72f0SBarry Smith #define GTOLNAME _Hash
243541bf97eSAdrian Croucher #define GTOLBS   mapping->bs
2449371c9d4SSatish Balay #define GTOL(g, local) \
2459371c9d4SSatish Balay   do { \
246e8f14785SLisandro Dalcin     (void)PetscHMapIGet(map->globalht, g / bs, &local); \
2470040bde1SJunchao Zhang     if (local >= 0) local = bs * local + (g % bs); \
248413f72f0SBarry Smith   } while (0)
249413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
250413f72f0SBarry Smith 
251413f72f0SBarry Smith #define GTOLTYPE _Hash
252413f72f0SBarry Smith #define GTOLNAME Block_Hash
253541bf97eSAdrian Croucher #define GTOLBS   1
2549371c9d4SSatish Balay #define GTOL(g, local) \
255d71ae5a4SJacob Faibussowitsch   do { \
256d71ae5a4SJacob Faibussowitsch     (void)PetscHMapIGet(map->globalht, g, &local); \
257d71ae5a4SJacob Faibussowitsch   } while (0)
258413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
259413f72f0SBarry Smith 
2606658fb44Sstefano_zampini /*@
2616658fb44Sstefano_zampini     ISLocalToGlobalMappingDuplicate - Duplicates the local to global mapping object
2626658fb44Sstefano_zampini 
2636658fb44Sstefano_zampini     Not Collective
2646658fb44Sstefano_zampini 
2656658fb44Sstefano_zampini     Input Parameter:
2666658fb44Sstefano_zampini .   ltog - local to global mapping
2676658fb44Sstefano_zampini 
2686658fb44Sstefano_zampini     Output Parameter:
2696658fb44Sstefano_zampini .   nltog - the duplicated local to global mapping
2706658fb44Sstefano_zampini 
2716658fb44Sstefano_zampini     Level: advanced
2726658fb44Sstefano_zampini 
273cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`
2746658fb44Sstefano_zampini @*/
275d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingDuplicate(ISLocalToGlobalMapping ltog, ISLocalToGlobalMapping *nltog)
276d71ae5a4SJacob Faibussowitsch {
277a0d79125SStefano Zampini   ISLocalToGlobalMappingType l2gtype;
2786658fb44Sstefano_zampini 
2796658fb44Sstefano_zampini   PetscFunctionBegin;
2806658fb44Sstefano_zampini   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
2819566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)ltog), ltog->bs, ltog->n, ltog->indices, PETSC_COPY_VALUES, nltog));
2829566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetType(ltog, &l2gtype));
2839566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingSetType(*nltog, l2gtype));
2843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2856658fb44Sstefano_zampini }
2866658fb44Sstefano_zampini 
287565245c5SBarry Smith /*@
288107e9a97SBarry Smith     ISLocalToGlobalMappingGetSize - Gets the local size of a local to global mapping
2893b9aefa3SBarry Smith 
2903b9aefa3SBarry Smith     Not Collective
2913b9aefa3SBarry Smith 
2923b9aefa3SBarry Smith     Input Parameter:
2933b9aefa3SBarry Smith .   ltog - local to global mapping
2943b9aefa3SBarry Smith 
2953b9aefa3SBarry Smith     Output Parameter:
296cab54364SBarry Smith .   n - the number of entries in the local mapping, `ISLocalToGlobalMappingGetIndices()` returns an array of this length
2973b9aefa3SBarry Smith 
2983b9aefa3SBarry Smith     Level: advanced
2993b9aefa3SBarry Smith 
300cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`
3013b9aefa3SBarry Smith @*/
302d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping mapping, PetscInt *n)
303d71ae5a4SJacob Faibussowitsch {
3043b9aefa3SBarry Smith   PetscFunctionBegin;
3050700a824SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
3064482741eSBarry Smith   PetscValidIntPointer(n, 2);
307107e9a97SBarry Smith   *n = mapping->bs * mapping->n;
3083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3093b9aefa3SBarry Smith }
3103b9aefa3SBarry Smith 
3115a5d4f66SBarry Smith /*@C
312cab54364SBarry Smith    ISLocalToGlobalMappingViewFromOptions - View an `ISLocalToGlobalMapping` based on values in the options database
313fe2efc57SMark 
314c3339decSBarry Smith    Collective
315fe2efc57SMark 
316fe2efc57SMark    Input Parameters:
317fe2efc57SMark +  A - the local to global mapping object
31820662ed9SBarry Smith .  obj - Optional object that provides the options prefix used for the options database query
319736c3998SJose E. Roman -  name - command line option
320fe2efc57SMark 
321fe2efc57SMark    Level: intermediate
322cab54364SBarry Smith 
32320662ed9SBarry Smith    Note:
32420662ed9SBarry Smith    See `PetscObjectViewFromOptions()` for the available `PetscViewer` and `PetscViewerFormat`
32520662ed9SBarry Smith 
32620662ed9SBarry Smith .seealso: [](sec_scatter), `PetscViewer`, ``ISLocalToGlobalMapping`, `ISLocalToGlobalMappingView`, `PetscObjectViewFromOptions()`, `ISLocalToGlobalMappingCreate()`
327fe2efc57SMark @*/
328d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingViewFromOptions(ISLocalToGlobalMapping A, PetscObject obj, const char name[])
329d71ae5a4SJacob Faibussowitsch {
330fe2efc57SMark   PetscFunctionBegin;
331fe2efc57SMark   PetscValidHeaderSpecific(A, IS_LTOGM_CLASSID, 1);
3329566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
3333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
334fe2efc57SMark }
335fe2efc57SMark 
336fe2efc57SMark /*@C
3375a5d4f66SBarry Smith     ISLocalToGlobalMappingView - View a local to global mapping
3385a5d4f66SBarry Smith 
339b9cd556bSLois Curfman McInnes     Not Collective
340b9cd556bSLois Curfman McInnes 
3415a5d4f66SBarry Smith     Input Parameters:
3423b9aefa3SBarry Smith +   ltog - local to global mapping
3433b9aefa3SBarry Smith -   viewer - viewer
3445a5d4f66SBarry Smith 
345a997ad1aSLois Curfman McInnes     Level: advanced
346a997ad1aSLois Curfman McInnes 
34720662ed9SBarry Smith .seealso: [](sec_scatter), `PetscViewer`, `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`
3485a5d4f66SBarry Smith @*/
349d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingView(ISLocalToGlobalMapping mapping, PetscViewer viewer)
350d71ae5a4SJacob Faibussowitsch {
35132dcc486SBarry Smith   PetscInt    i;
35232dcc486SBarry Smith   PetscMPIInt rank;
353ace3abfcSBarry Smith   PetscBool   iascii;
3545a5d4f66SBarry Smith 
3555a5d4f66SBarry Smith   PetscFunctionBegin;
3560700a824SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
35748a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mapping), &viewer));
3580700a824SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
3595a5d4f66SBarry Smith 
3609566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mapping), &rank));
3619566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
36232077d6dSBarry Smith   if (iascii) {
3639566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)mapping, viewer));
3649566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
365f2c6b1a2SJed Brown     for (i = 0; i < mapping->n; i++) {
366f2c6b1a2SJed Brown       PetscInt bs = mapping->bs, g = mapping->indices[i];
367f2c6b1a2SJed Brown       if (bs == 1) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %" PetscInt_FMT " %" PetscInt_FMT "\n", rank, i, g));
368f2c6b1a2SJed Brown       else PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %" PetscInt_FMT ":%" PetscInt_FMT " %" PetscInt_FMT ":%" PetscInt_FMT "\n", rank, i * bs, (i + 1) * bs, g * bs, (g + 1) * bs));
369f2c6b1a2SJed Brown     }
3709566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
3719566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
3721575c14dSBarry Smith   }
3733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3745a5d4f66SBarry Smith }
3755a5d4f66SBarry Smith 
3761f428162SBarry Smith /*@
3772bdab257SBarry Smith     ISLocalToGlobalMappingCreateIS - Creates a mapping between a local (0 to n)
3782bdab257SBarry Smith     ordering and a global parallel ordering.
3792bdab257SBarry Smith 
38020662ed9SBarry Smith     Not Collective
381b9cd556bSLois Curfman McInnes 
382a997ad1aSLois Curfman McInnes     Input Parameter:
3838c03b21aSDmitry Karpeev .   is - index set containing the global numbers for each local number
3842bdab257SBarry Smith 
385a997ad1aSLois Curfman McInnes     Output Parameter:
3862bdab257SBarry Smith .   mapping - new mapping data structure
3872bdab257SBarry Smith 
388a997ad1aSLois Curfman McInnes     Level: advanced
389a997ad1aSLois Curfman McInnes 
390cab54364SBarry Smith     Note:
391cab54364SBarry Smith     the block size of the `IS` determines the block size of the mapping
392cab54364SBarry Smith 
393cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetFromOptions()`
3942bdab257SBarry Smith @*/
395d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingCreateIS(IS is, ISLocalToGlobalMapping *mapping)
396d71ae5a4SJacob Faibussowitsch {
3973bbf0e92SBarry Smith   PetscInt        n, bs;
3985d0c19d7SBarry Smith   const PetscInt *indices;
3992bdab257SBarry Smith   MPI_Comm        comm;
4003bbf0e92SBarry Smith   PetscBool       isblock;
4013a40ed3dSBarry Smith 
4023a40ed3dSBarry Smith   PetscFunctionBegin;
4030700a824SBarry Smith   PetscValidHeaderSpecific(is, IS_CLASSID, 1);
4044482741eSBarry Smith   PetscValidPointer(mapping, 2);
4052bdab257SBarry Smith 
4069566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)is, &comm));
4079566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is, &n));
4089566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)is, ISBLOCK, &isblock));
4096006e8d2SBarry Smith   if (!isblock) {
4109566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(is, &indices));
4119566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(comm, 1, n, indices, PETSC_COPY_VALUES, mapping));
4129566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(is, &indices));
4136006e8d2SBarry Smith   } else {
4149566063dSJacob Faibussowitsch     PetscCall(ISGetBlockSize(is, &bs));
4159566063dSJacob Faibussowitsch     PetscCall(ISBlockGetIndices(is, &indices));
4169566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(comm, bs, n / bs, indices, PETSC_COPY_VALUES, mapping));
4179566063dSJacob Faibussowitsch     PetscCall(ISBlockRestoreIndices(is, &indices));
4186006e8d2SBarry Smith   }
4193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4202bdab257SBarry Smith }
4215a5d4f66SBarry Smith 
422a4d96a55SJed Brown /*@C
423a4d96a55SJed Brown     ISLocalToGlobalMappingCreateSF - Creates a mapping between a local (0 to n)
424a4d96a55SJed Brown     ordering and a global parallel ordering.
425a4d96a55SJed Brown 
426a4d96a55SJed Brown     Collective
427a4d96a55SJed Brown 
428d8d19677SJose E. Roman     Input Parameters:
429a4d96a55SJed Brown +   sf - star forest mapping contiguous local indices to (rank, offset)
430cab54364SBarry Smith -   start - first global index on this process, or `PETSC_DECIDE` to compute contiguous global numbering automatically
431a4d96a55SJed Brown 
432a4d96a55SJed Brown     Output Parameter:
433a4d96a55SJed Brown .   mapping - new mapping data structure
434a4d96a55SJed Brown 
435a4d96a55SJed Brown     Level: advanced
436a4d96a55SJed Brown 
43720662ed9SBarry Smith     Note:
43820662ed9SBarry Smith     If any processor calls this with `start` = `PETSC_DECIDE` then all processors must, otherwise the program will hang.
4399a535bafSVaclav Hapla 
440cab54364SBarry Smith .seealso: [](sec_scatter), `PetscSF`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingSetFromOptions()`
441a4d96a55SJed Brown @*/
442d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF sf, PetscInt start, ISLocalToGlobalMapping *mapping)
443d71ae5a4SJacob Faibussowitsch {
444a4d96a55SJed Brown   PetscInt i, maxlocal, nroots, nleaves, *globals, *ltog;
445a4d96a55SJed Brown   MPI_Comm comm;
446a4d96a55SJed Brown 
447a4d96a55SJed Brown   PetscFunctionBegin;
448a4d96a55SJed Brown   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
449a4d96a55SJed Brown   PetscValidPointer(mapping, 3);
4509566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
45141f4c31fSVaclav Hapla   PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL));
4529a535bafSVaclav Hapla   if (start == PETSC_DECIDE) {
4539a535bafSVaclav Hapla     start = 0;
4549566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Exscan(&nroots, &start, 1, MPIU_INT, MPI_SUM, comm));
45541f4c31fSVaclav Hapla   } else PetscCheck(start >= 0, comm, PETSC_ERR_ARG_OUTOFRANGE, "start must be nonnegative or PETSC_DECIDE");
45641f4c31fSVaclav Hapla   PetscCall(PetscSFGetLeafRange(sf, NULL, &maxlocal));
45741f4c31fSVaclav Hapla   ++maxlocal;
4589566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nroots, &globals));
4599566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxlocal, &ltog));
460a4d96a55SJed Brown   for (i = 0; i < nroots; i++) globals[i] = start + i;
461a4d96a55SJed Brown   for (i = 0; i < maxlocal; i++) ltog[i] = -1;
4629566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sf, MPIU_INT, globals, ltog, MPI_REPLACE));
4639566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sf, MPIU_INT, globals, ltog, MPI_REPLACE));
4649566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(comm, 1, maxlocal, ltog, PETSC_OWN_POINTER, mapping));
4659566063dSJacob Faibussowitsch   PetscCall(PetscFree(globals));
4663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
467a4d96a55SJed Brown }
468b46b645bSBarry Smith 
46963fa5c83Sstefano_zampini /*@
47063fa5c83Sstefano_zampini     ISLocalToGlobalMappingSetBlockSize - Sets the blocksize of the mapping
47163fa5c83Sstefano_zampini 
47220662ed9SBarry Smith     Not Collective
47363fa5c83Sstefano_zampini 
47463fa5c83Sstefano_zampini     Input Parameters:
475a2b725a8SWilliam Gropp +   mapping - mapping data structure
476a2b725a8SWilliam Gropp -   bs - the blocksize
47763fa5c83Sstefano_zampini 
47863fa5c83Sstefano_zampini     Level: advanced
47963fa5c83Sstefano_zampini 
480cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`
48163fa5c83Sstefano_zampini @*/
482d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingSetBlockSize(ISLocalToGlobalMapping mapping, PetscInt bs)
483d71ae5a4SJacob Faibussowitsch {
484a59f3c4dSstefano_zampini   PetscInt       *nid;
485a59f3c4dSstefano_zampini   const PetscInt *oid;
486a59f3c4dSstefano_zampini   PetscInt        i, cn, on, obs, nn;
48763fa5c83Sstefano_zampini 
48863fa5c83Sstefano_zampini   PetscFunctionBegin;
48963fa5c83Sstefano_zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
49008401ef6SPierre Jolivet   PetscCheck(bs >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid block size %" PetscInt_FMT, bs);
4913ba16761SJacob Faibussowitsch   if (bs == mapping->bs) PetscFunctionReturn(PETSC_SUCCESS);
49263fa5c83Sstefano_zampini   on  = mapping->n;
49363fa5c83Sstefano_zampini   obs = mapping->bs;
49463fa5c83Sstefano_zampini   oid = mapping->indices;
49563fa5c83Sstefano_zampini   nn  = (on * obs) / bs;
49608401ef6SPierre 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);
497a59f3c4dSstefano_zampini 
4989566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nn, &nid));
4999566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetIndices(mapping, &oid));
500a59f3c4dSstefano_zampini   for (i = 0; i < nn; i++) {
501a59f3c4dSstefano_zampini     PetscInt j;
502a59f3c4dSstefano_zampini     for (j = 0, cn = 0; j < bs - 1; j++) {
5039371c9d4SSatish Balay       if (oid[i * bs + j] < 0) {
5049371c9d4SSatish Balay         cn++;
5059371c9d4SSatish Balay         continue;
5069371c9d4SSatish Balay       }
50708401ef6SPierre 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]);
508a59f3c4dSstefano_zampini     }
509a59f3c4dSstefano_zampini     if (oid[i * bs + j] < 0) cn++;
5108b7cb0e6Sstefano_zampini     if (cn) {
51108401ef6SPierre 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);
512a59f3c4dSstefano_zampini       nid[i] = -1;
5138b7cb0e6Sstefano_zampini     } else {
514a59f3c4dSstefano_zampini       nid[i] = oid[i * bs] / bs;
51563fa5c83Sstefano_zampini     }
51663fa5c83Sstefano_zampini   }
5179566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreIndices(mapping, &oid));
518a59f3c4dSstefano_zampini 
51963fa5c83Sstefano_zampini   mapping->n  = nn;
52063fa5c83Sstefano_zampini   mapping->bs = bs;
5219566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping->indices));
52263fa5c83Sstefano_zampini   mapping->indices     = nid;
523c9345713Sstefano_zampini   mapping->globalstart = 0;
524c9345713Sstefano_zampini   mapping->globalend   = 0;
5251bd0b88eSStefano Zampini 
5261bd0b88eSStefano Zampini   /* reset the cached information */
5279566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping->info_procs));
5289566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping->info_numprocs));
5291bd0b88eSStefano Zampini   if (mapping->info_indices) {
5301bd0b88eSStefano Zampini     PetscInt i;
5311bd0b88eSStefano Zampini 
5329566063dSJacob Faibussowitsch     PetscCall(PetscFree((mapping->info_indices)[0]));
53348a46eb9SPierre Jolivet     for (i = 1; i < mapping->info_nproc; i++) PetscCall(PetscFree(mapping->info_indices[i]));
5349566063dSJacob Faibussowitsch     PetscCall(PetscFree(mapping->info_indices));
5351bd0b88eSStefano Zampini   }
5361bd0b88eSStefano Zampini   mapping->info_cached = PETSC_FALSE;
5371bd0b88eSStefano Zampini 
538dbbe0bcdSBarry Smith   PetscTryTypeMethod(mapping, destroy);
5393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
54063fa5c83Sstefano_zampini }
54163fa5c83Sstefano_zampini 
54245b6f7e9SBarry Smith /*@
54345b6f7e9SBarry Smith     ISLocalToGlobalMappingGetBlockSize - Gets the blocksize of the mapping
54445b6f7e9SBarry Smith     ordering and a global parallel ordering.
54545b6f7e9SBarry Smith 
54645b6f7e9SBarry Smith     Not Collective
54745b6f7e9SBarry Smith 
548*2fe279fdSBarry Smith     Input Parameter:
54945b6f7e9SBarry Smith .   mapping - mapping data structure
55045b6f7e9SBarry Smith 
55145b6f7e9SBarry Smith     Output Parameter:
55245b6f7e9SBarry Smith .   bs - the blocksize
55345b6f7e9SBarry Smith 
55445b6f7e9SBarry Smith     Level: advanced
55545b6f7e9SBarry Smith 
556cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`
55745b6f7e9SBarry Smith @*/
558d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping mapping, PetscInt *bs)
559d71ae5a4SJacob Faibussowitsch {
56045b6f7e9SBarry Smith   PetscFunctionBegin;
561cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
56245b6f7e9SBarry Smith   *bs = mapping->bs;
5633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
56445b6f7e9SBarry Smith }
56545b6f7e9SBarry Smith 
566ba5bb76aSSatish Balay /*@
56790f02eecSBarry Smith     ISLocalToGlobalMappingCreate - Creates a mapping between a local (0 to n)
56890f02eecSBarry Smith     ordering and a global parallel ordering.
5692362add9SBarry Smith 
57089d82c54SBarry Smith     Not Collective, but communicator may have more than one process
571b9cd556bSLois Curfman McInnes 
5722362add9SBarry Smith     Input Parameters:
57389d82c54SBarry Smith +   comm - MPI communicator
574f0413b6fSBarry Smith .   bs - the block size
57528bc9809SBarry Smith .   n - the number of local elements divided by the block size, or equivalently the number of block indices
57628bc9809SBarry 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
577d5ad8652SBarry Smith -   mode - see PetscCopyMode
5782362add9SBarry Smith 
579a997ad1aSLois Curfman McInnes     Output Parameter:
58090f02eecSBarry Smith .   mapping - new mapping data structure
5812362add9SBarry Smith 
582cab54364SBarry Smith     Level: advanced
583cab54364SBarry Smith 
58495452b02SPatrick Sanan     Notes:
58595452b02SPatrick Sanan     There is one integer value in indices per block and it represents the actual indices bs*idx + j, where j=0,..,bs-1
586413f72f0SBarry Smith 
587cab54364SBarry Smith     For "small" problems when using `ISGlobalToLocalMappingApply()` and `ISGlobalToLocalMappingApplyBlock()`, the `ISLocalToGlobalMappingType`
588cab54364SBarry Smith     of `ISLOCALTOGLOBALMAPPINGBASIC` will be used; this uses more memory but is faster; this approach is not scalable for extremely large mappings.
589413f72f0SBarry Smith 
590cab54364SBarry Smith     For large problems `ISLOCALTOGLOBALMAPPINGHASH` is used, this is scalable.
59120662ed9SBarry Smith     Use `ISLocalToGlobalMappingSetType()` or call `ISLocalToGlobalMappingSetFromOptions()` with the option
59220662ed9SBarry Smith     `-islocaltoglobalmapping_type` <`basic`,`hash`> to control which is used.
593a997ad1aSLois Curfman McInnes 
59420662ed9SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingSetFromOptions()`,
59520662ed9SBarry Smith           `ISLOCALTOGLOBALMAPPINGBASIC`, `ISLOCALTOGLOBALMAPPINGHASH`
596db781477SPatrick Sanan           `ISLocalToGlobalMappingSetType()`, `ISLocalToGlobalMappingType`
5972362add9SBarry Smith @*/
598d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingCreate(MPI_Comm comm, PetscInt bs, PetscInt n, const PetscInt indices[], PetscCopyMode mode, ISLocalToGlobalMapping *mapping)
599d71ae5a4SJacob Faibussowitsch {
60032dcc486SBarry Smith   PetscInt *in;
601b46b645bSBarry Smith 
602b46b645bSBarry Smith   PetscFunctionBegin;
603064a246eSJacob Faibussowitsch   if (n) PetscValidIntPointer(indices, 4);
604064a246eSJacob Faibussowitsch   PetscValidPointer(mapping, 6);
605b46b645bSBarry Smith 
6060298fd71SBarry Smith   *mapping = NULL;
6079566063dSJacob Faibussowitsch   PetscCall(ISInitializePackage());
6082362add9SBarry Smith 
6099566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(*mapping, IS_LTOGM_CLASSID, "ISLocalToGlobalMapping", "Local to global mapping", "IS", comm, ISLocalToGlobalMappingDestroy, ISLocalToGlobalMappingView));
610d4bb536fSBarry Smith   (*mapping)->n  = n;
611f0413b6fSBarry Smith   (*mapping)->bs = bs;
612d5ad8652SBarry Smith   if (mode == PETSC_COPY_VALUES) {
6139566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &in));
6149566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(in, indices, n));
615d5ad8652SBarry Smith     (*mapping)->indices         = in;
61671910c26SVaclav Hapla     (*mapping)->dealloc_indices = PETSC_TRUE;
6176389a1a1SBarry Smith   } else if (mode == PETSC_OWN_POINTER) {
6186389a1a1SBarry Smith     (*mapping)->indices         = (PetscInt *)indices;
61971910c26SVaclav Hapla     (*mapping)->dealloc_indices = PETSC_TRUE;
62071910c26SVaclav Hapla   } else if (mode == PETSC_USE_POINTER) {
62171910c26SVaclav Hapla     (*mapping)->indices = (PetscInt *)indices;
6229371c9d4SSatish Balay   } else SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid mode %d", mode);
6233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6242362add9SBarry Smith }
6252362add9SBarry Smith 
626413f72f0SBarry Smith PetscFunctionList ISLocalToGlobalMappingList = NULL;
627413f72f0SBarry Smith 
62890f02eecSBarry Smith /*@
6297e99dc12SLawrence Mitchell    ISLocalToGlobalMappingSetFromOptions - Set mapping options from the options database.
6307e99dc12SLawrence Mitchell 
63120662ed9SBarry Smith    Not Collective
6327e99dc12SLawrence Mitchell 
633*2fe279fdSBarry Smith    Input Parameter:
6347e99dc12SLawrence Mitchell .  mapping - mapping data structure
6357e99dc12SLawrence Mitchell 
63620662ed9SBarry Smith    Options Database Key:
63720662ed9SBarry Smith .  -islocaltoglobalmapping_type - <basic,hash> nonscalable and scalable versions
63820662ed9SBarry Smith 
6397e99dc12SLawrence Mitchell    Level: advanced
6407e99dc12SLawrence Mitchell 
64120662ed9SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingSetFromOptions()`,
64220662ed9SBarry Smith           `ISLOCALTOGLOBALMAPPINGBASIC`, `ISLOCALTOGLOBALMAPPINGHASH`
643cab54364SBarry Smith           `ISLocalToGlobalMappingSetType()`, `ISLocalToGlobalMappingType`
6447e99dc12SLawrence Mitchell @*/
645d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingSetFromOptions(ISLocalToGlobalMapping mapping)
646d71ae5a4SJacob Faibussowitsch {
647413f72f0SBarry Smith   char                       type[256];
648413f72f0SBarry Smith   ISLocalToGlobalMappingType defaulttype = "Not set";
6497e99dc12SLawrence Mitchell   PetscBool                  flg;
6507e99dc12SLawrence Mitchell 
6517e99dc12SLawrence Mitchell   PetscFunctionBegin;
6527e99dc12SLawrence Mitchell   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
6539566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRegisterAll());
654d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)mapping);
6559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-islocaltoglobalmapping_type", "ISLocalToGlobalMapping method", "ISLocalToGlobalMappingSetType", ISLocalToGlobalMappingList, (char *)(((PetscObject)mapping)->type_name) ? ((PetscObject)mapping)->type_name : defaulttype, type, 256, &flg));
6561baa6e33SBarry Smith   if (flg) PetscCall(ISLocalToGlobalMappingSetType(mapping, type));
657d0609cedSBarry Smith   PetscOptionsEnd();
6583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6597e99dc12SLawrence Mitchell }
6607e99dc12SLawrence Mitchell 
6617e99dc12SLawrence Mitchell /*@
66290f02eecSBarry Smith    ISLocalToGlobalMappingDestroy - Destroys a mapping between a local (0 to n)
66390f02eecSBarry Smith    ordering and a global parallel ordering.
66490f02eecSBarry Smith 
66520662ed9SBarry Smith    Not Collective
666b9cd556bSLois Curfman McInnes 
667*2fe279fdSBarry Smith    Input Parameter:
66890f02eecSBarry Smith .  mapping - mapping data structure
66990f02eecSBarry Smith 
670a997ad1aSLois Curfman McInnes    Level: advanced
671a997ad1aSLois Curfman McInnes 
672cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingCreate()`
67390f02eecSBarry Smith @*/
674d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping *mapping)
675d71ae5a4SJacob Faibussowitsch {
6763a40ed3dSBarry Smith   PetscFunctionBegin;
6773ba16761SJacob Faibussowitsch   if (!*mapping) PetscFunctionReturn(PETSC_SUCCESS);
6786bf464f9SBarry Smith   PetscValidHeaderSpecific((*mapping), IS_LTOGM_CLASSID, 1);
6799371c9d4SSatish Balay   if (--((PetscObject)(*mapping))->refct > 0) {
6809371c9d4SSatish Balay     *mapping = NULL;
6813ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
68271910c26SVaclav Hapla   }
68348a46eb9SPierre Jolivet   if ((*mapping)->dealloc_indices) PetscCall(PetscFree((*mapping)->indices));
6849566063dSJacob Faibussowitsch   PetscCall(PetscFree((*mapping)->info_procs));
6859566063dSJacob Faibussowitsch   PetscCall(PetscFree((*mapping)->info_numprocs));
686268a049cSStefano Zampini   if ((*mapping)->info_indices) {
687268a049cSStefano Zampini     PetscInt i;
688268a049cSStefano Zampini 
6899566063dSJacob Faibussowitsch     PetscCall(PetscFree(((*mapping)->info_indices)[0]));
69048a46eb9SPierre Jolivet     for (i = 1; i < (*mapping)->info_nproc; i++) PetscCall(PetscFree(((*mapping)->info_indices)[i]));
6919566063dSJacob Faibussowitsch     PetscCall(PetscFree((*mapping)->info_indices));
692268a049cSStefano Zampini   }
69348a46eb9SPierre Jolivet   if ((*mapping)->info_nodei) PetscCall(PetscFree(((*mapping)->info_nodei)[0]));
6949566063dSJacob Faibussowitsch   PetscCall(PetscFree2((*mapping)->info_nodec, (*mapping)->info_nodei));
69548a46eb9SPierre Jolivet   if ((*mapping)->ops->destroy) PetscCall((*(*mapping)->ops->destroy)(*mapping));
6969566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(mapping));
6974c8fdceaSLisandro Dalcin   *mapping = NULL;
6983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
69990f02eecSBarry Smith }
70090f02eecSBarry Smith 
70190f02eecSBarry Smith /*@
702cab54364SBarry Smith     ISLocalToGlobalMappingApplyIS - Creates from an `IS` in the local numbering
703cab54364SBarry Smith     a new index set using the global numbering defined in an `ISLocalToGlobalMapping`
7043acfe500SLois Curfman McInnes     context.
70590f02eecSBarry Smith 
706c3339decSBarry Smith     Collective
707b9cd556bSLois Curfman McInnes 
70890f02eecSBarry Smith     Input Parameters:
709b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
710b9cd556bSLois Curfman McInnes -   is - index set in local numbering
71190f02eecSBarry Smith 
712cab54364SBarry Smith     Output Parameter:
71390f02eecSBarry Smith .   newis - index set in global numbering
71490f02eecSBarry Smith 
715a997ad1aSLois Curfman McInnes     Level: advanced
716a997ad1aSLois Curfman McInnes 
717cab54364SBarry Smith     Note:
71820662ed9SBarry Smith     The output `IS` will have the same communicator as the input `IS`.
719cab54364SBarry Smith 
720cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingCreate()`,
721db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`, `ISGlobalToLocalMappingApply()`
72290f02eecSBarry Smith @*/
723d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping, IS is, IS *newis)
724d71ae5a4SJacob Faibussowitsch {
725e24637baSBarry Smith   PetscInt        n, *idxout;
7265d0c19d7SBarry Smith   const PetscInt *idxin;
7273a40ed3dSBarry Smith 
7283a40ed3dSBarry Smith   PetscFunctionBegin;
7290700a824SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
7300700a824SBarry Smith   PetscValidHeaderSpecific(is, IS_CLASSID, 2);
7314482741eSBarry Smith   PetscValidPointer(newis, 3);
73290f02eecSBarry Smith 
7339566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is, &n));
7349566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is, &idxin));
7359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &idxout));
7369566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(mapping, n, idxin, idxout));
7379566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is, &idxin));
7389566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)is), n, idxout, PETSC_OWN_POINTER, newis));
7393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
74090f02eecSBarry Smith }
74190f02eecSBarry Smith 
742b89cb25eSSatish Balay /*@
7433acfe500SLois Curfman McInnes    ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering
7443acfe500SLois Curfman McInnes    and converts them to the global numbering.
74590f02eecSBarry Smith 
74620662ed9SBarry Smith    Not Collective
747b9cd556bSLois Curfman McInnes 
748bb25748dSBarry Smith    Input Parameters:
749b9cd556bSLois Curfman McInnes +  mapping - the local to global mapping context
750bb25748dSBarry Smith .  N - number of integers
751b9cd556bSLois Curfman McInnes -  in - input indices in local numbering
752bb25748dSBarry Smith 
753bb25748dSBarry Smith    Output Parameter:
754bb25748dSBarry Smith .  out - indices in global numbering
755bb25748dSBarry Smith 
756a997ad1aSLois Curfman McInnes    Level: advanced
757a997ad1aSLois Curfman McInnes 
758cab54364SBarry Smith    Note:
75920662ed9SBarry Smith    The `in` and `out` array parameters may be identical.
760cab54364SBarry Smith 
761cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApplyBlock()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingDestroy()`,
762c2e3fba1SPatrick Sanan           `ISLocalToGlobalMappingApplyIS()`, `AOCreateBasic()`, `AOApplicationToPetsc()`,
763db781477SPatrick Sanan           `AOPetscToApplication()`, `ISGlobalToLocalMappingApply()`
764afcb2eb5SJed Brown @*/
765d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping, PetscInt N, const PetscInt in[], PetscInt out[])
766d71ae5a4SJacob Faibussowitsch {
767cbc1caf0SMatthew G. Knepley   PetscInt i, bs, Nmax;
76845b6f7e9SBarry Smith 
76945b6f7e9SBarry Smith   PetscFunctionBegin;
770cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
771cbc1caf0SMatthew G. Knepley   bs   = mapping->bs;
772cbc1caf0SMatthew G. Knepley   Nmax = bs * mapping->n;
77345b6f7e9SBarry Smith   if (bs == 1) {
774cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
77545b6f7e9SBarry Smith     for (i = 0; i < N; i++) {
77645b6f7e9SBarry Smith       if (in[i] < 0) {
77745b6f7e9SBarry Smith         out[i] = in[i];
77845b6f7e9SBarry Smith         continue;
77945b6f7e9SBarry Smith       }
78008401ef6SPierre 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);
78145b6f7e9SBarry Smith       out[i] = idx[in[i]];
78245b6f7e9SBarry Smith     }
78345b6f7e9SBarry Smith   } else {
784cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
78545b6f7e9SBarry Smith     for (i = 0; i < N; i++) {
78645b6f7e9SBarry Smith       if (in[i] < 0) {
78745b6f7e9SBarry Smith         out[i] = in[i];
78845b6f7e9SBarry Smith         continue;
78945b6f7e9SBarry Smith       }
79008401ef6SPierre 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);
79145b6f7e9SBarry Smith       out[i] = idx[in[i] / bs] * bs + (in[i] % bs);
79245b6f7e9SBarry Smith     }
79345b6f7e9SBarry Smith   }
7943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
79545b6f7e9SBarry Smith }
79645b6f7e9SBarry Smith 
79745b6f7e9SBarry Smith /*@
7986006e8d2SBarry Smith    ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local block numbering and converts them to the global block numbering
79945b6f7e9SBarry Smith 
80020662ed9SBarry Smith    Not Collective
80145b6f7e9SBarry Smith 
80245b6f7e9SBarry Smith    Input Parameters:
80345b6f7e9SBarry Smith +  mapping - the local to global mapping context
80445b6f7e9SBarry Smith .  N - number of integers
8056006e8d2SBarry Smith -  in - input indices in local block numbering
80645b6f7e9SBarry Smith 
80745b6f7e9SBarry Smith    Output Parameter:
8086006e8d2SBarry Smith .  out - indices in global block numbering
80945b6f7e9SBarry Smith 
8106006e8d2SBarry Smith    Example:
811cab54364SBarry 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
8126006e8d2SBarry Smith      (the first block) would produce 0 and the mapping applied to 1 (the second block) would produce 3.
8136006e8d2SBarry Smith 
81420662ed9SBarry Smith    Level: advanced
81520662ed9SBarry Smith 
81620662ed9SBarry Smith    Note:
81720662ed9SBarry Smith    The `in` and `out` array parameters may be identical.
81820662ed9SBarry Smith 
819cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingDestroy()`,
820c2e3fba1SPatrick Sanan           `ISLocalToGlobalMappingApplyIS()`, `AOCreateBasic()`, `AOApplicationToPetsc()`,
821db781477SPatrick Sanan           `AOPetscToApplication()`, `ISGlobalToLocalMappingApply()`
82245b6f7e9SBarry Smith @*/
823d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping, PetscInt N, const PetscInt in[], PetscInt out[])
824d71ae5a4SJacob Faibussowitsch {
8258a1f772fSStefano Zampini   PetscInt        i, Nmax;
8268a1f772fSStefano Zampini   const PetscInt *idx;
827d4bb536fSBarry Smith 
828a0d79125SStefano Zampini   PetscFunctionBegin;
829a0d79125SStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
8308a1f772fSStefano Zampini   Nmax = mapping->n;
8318a1f772fSStefano Zampini   idx  = mapping->indices;
832afcb2eb5SJed Brown   for (i = 0; i < N; i++) {
833afcb2eb5SJed Brown     if (in[i] < 0) {
834afcb2eb5SJed Brown       out[i] = in[i];
835afcb2eb5SJed Brown       continue;
836afcb2eb5SJed Brown     }
83708401ef6SPierre 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);
838afcb2eb5SJed Brown     out[i] = idx[in[i]];
839afcb2eb5SJed Brown   }
8403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
841afcb2eb5SJed Brown }
842d4bb536fSBarry Smith 
8437e99dc12SLawrence Mitchell /*@
844a997ad1aSLois Curfman McInnes     ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers
845a997ad1aSLois Curfman McInnes     specified with a global numbering.
846d4bb536fSBarry Smith 
84720662ed9SBarry Smith     Not Collective
848b9cd556bSLois Curfman McInnes 
849d4bb536fSBarry Smith     Input Parameters:
850b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
851cab54364SBarry Smith .   type - `IS_GTOLM_MASK` - maps global indices with no local value to -1 in the output list (i.e., mask them)
852cab54364SBarry Smith            `IS_GTOLM_DROP` - drops the indices with no local value from the output list
853d4bb536fSBarry Smith .   n - number of global indices to map
854b9cd556bSLois Curfman McInnes -   idx - global indices to map
855d4bb536fSBarry Smith 
856d4bb536fSBarry Smith     Output Parameters:
857cab54364SBarry Smith +   nout - number of indices in output array (if type == `IS_GTOLM_MASK` then nout = n)
858b9cd556bSLois Curfman McInnes -   idxout - local index of each global index, one must pass in an array long enough
859cab54364SBarry Smith              to hold all the indices. You can call `ISGlobalToLocalMappingApply()` with
8600298fd71SBarry Smith              idxout == NULL to determine the required length (returned in nout)
861cab54364SBarry Smith              and then allocate the required space and call `ISGlobalToLocalMappingApply()`
862e182c471SBarry Smith              a second time to set the values.
863d4bb536fSBarry Smith 
864cab54364SBarry Smith     Level: advanced
865cab54364SBarry Smith 
866b9cd556bSLois Curfman McInnes     Notes:
86720662ed9SBarry Smith     Either `nout` or `idxout` may be `NULL`. `idx` and `idxout` may be identical.
868d4bb536fSBarry Smith 
86920662ed9SBarry Smith     For "small" problems when using `ISGlobalToLocalMappingApply()` and `ISGlobalToLocalMappingApplyBlock()`, the `ISLocalToGlobalMappingType` of
87020662ed9SBarry Smith     `ISLOCALTOGLOBALMAPPINGBASIC` will be used;
871cab54364SBarry 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.
872cab54364SBarry Smith     Use `ISLocalToGlobalMappingSetType()` or call `ISLocalToGlobalMappingSetFromOptions()` with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
8730f5bd95cSBarry Smith 
874cab54364SBarry Smith     Developer Note:
87520662ed9SBarry Smith     The manual page states that `idx` and `idxout` may be identical but the calling
87620662ed9SBarry Smith     sequence declares `idx` as const so it cannot be the same as `idxout`.
87732fd6b96SBarry Smith 
878cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApply()`, `ISGlobalToLocalMappingApplyBlock()`, `ISLocalToGlobalMappingCreate()`,
879db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`
880d4bb536fSBarry Smith @*/
881d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, PetscInt n, const PetscInt idx[], PetscInt *nout, PetscInt idxout[])
882d71ae5a4SJacob Faibussowitsch {
8839d90f715SBarry Smith   PetscFunctionBegin;
8849d90f715SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
88548a46eb9SPierre Jolivet   if (!mapping->data) PetscCall(ISGlobalToLocalMappingSetUp(mapping));
886dbbe0bcdSBarry Smith   PetscUseTypeMethod(mapping, globaltolocalmappingapply, type, n, idx, nout, idxout);
8873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8889d90f715SBarry Smith }
8899d90f715SBarry Smith 
890d4fe737eSStefano Zampini /*@
891cab54364SBarry Smith     ISGlobalToLocalMappingApplyIS - Creates from an `IS` in the global numbering
892cab54364SBarry Smith     a new index set using the local numbering defined in an `ISLocalToGlobalMapping`
893d4fe737eSStefano Zampini     context.
894d4fe737eSStefano Zampini 
89520662ed9SBarry Smith     Not Collective
896d4fe737eSStefano Zampini 
897d4fe737eSStefano Zampini     Input Parameters:
898d4fe737eSStefano Zampini +   mapping - mapping between local and global numbering
899cab54364SBarry Smith .   type - `IS_GTOLM_MASK` - maps global indices with no local value to -1 in the output list (i.e., mask them)
900cab54364SBarry Smith            `IS_GTOLM_DROP` - drops the indices with no local value from the output list
901d4fe737eSStefano Zampini -   is - index set in global numbering
902d4fe737eSStefano Zampini 
903*2fe279fdSBarry Smith     Output Parameter:
904d4fe737eSStefano Zampini .   newis - index set in local numbering
905d4fe737eSStefano Zampini 
906d4fe737eSStefano Zampini     Level: advanced
907d4fe737eSStefano Zampini 
908cab54364SBarry Smith     Note:
909cab54364SBarry Smith     The output `IS` will be sequential, as it encodes a purely local operation
910cab54364SBarry Smith 
911cab54364SBarry Smith .seealso: [](sec_scatter), `ISGlobalToLocalMapping`, `ISGlobalToLocalMappingApply()`, `ISLocalToGlobalMappingCreate()`,
912db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`
913d4fe737eSStefano Zampini @*/
914d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, IS is, IS *newis)
915d71ae5a4SJacob Faibussowitsch {
916d4fe737eSStefano Zampini   PetscInt        n, nout, *idxout;
917d4fe737eSStefano Zampini   const PetscInt *idxin;
918d4fe737eSStefano Zampini 
919d4fe737eSStefano Zampini   PetscFunctionBegin;
920d4fe737eSStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
921d4fe737eSStefano Zampini   PetscValidHeaderSpecific(is, IS_CLASSID, 3);
922d4fe737eSStefano Zampini   PetscValidPointer(newis, 4);
923d4fe737eSStefano Zampini 
9249566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is, &n));
9259566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is, &idxin));
926d4fe737eSStefano Zampini   if (type == IS_GTOLM_MASK) {
9279566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &idxout));
928d4fe737eSStefano Zampini   } else {
9299566063dSJacob Faibussowitsch     PetscCall(ISGlobalToLocalMappingApply(mapping, type, n, idxin, &nout, NULL));
9309566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nout, &idxout));
931d4fe737eSStefano Zampini   }
9329566063dSJacob Faibussowitsch   PetscCall(ISGlobalToLocalMappingApply(mapping, type, n, idxin, &nout, idxout));
9339566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is, &idxin));
9349566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nout, idxout, PETSC_OWN_POINTER, newis));
9353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
936d4fe737eSStefano Zampini }
937d4fe737eSStefano Zampini 
9389d90f715SBarry Smith /*@
9399d90f715SBarry Smith     ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers
9409d90f715SBarry Smith     specified with a block global numbering.
9419d90f715SBarry Smith 
94220662ed9SBarry Smith     Not Collective
9439d90f715SBarry Smith 
9449d90f715SBarry Smith     Input Parameters:
9459d90f715SBarry Smith +   mapping - mapping between local and global numbering
946cab54364SBarry Smith .   type - `IS_GTOLM_MASK` - maps global indices with no local value to -1 in the output list (i.e., mask them)
947cab54364SBarry Smith            `IS_GTOLM_DROP` - drops the indices with no local value from the output list
9489d90f715SBarry Smith .   n - number of global indices to map
9499d90f715SBarry Smith -   idx - global indices to map
9509d90f715SBarry Smith 
9519d90f715SBarry Smith     Output Parameters:
952cab54364SBarry Smith +   nout - number of indices in output array (if type == `IS_GTOLM_MASK` then nout = n)
9539d90f715SBarry Smith -   idxout - local index of each global index, one must pass in an array long enough
954cab54364SBarry Smith              to hold all the indices. You can call `ISGlobalToLocalMappingApplyBlock()` with
9559d90f715SBarry Smith              idxout == NULL to determine the required length (returned in nout)
956cab54364SBarry Smith              and then allocate the required space and call `ISGlobalToLocalMappingApplyBlock()`
9579d90f715SBarry Smith              a second time to set the values.
9589d90f715SBarry Smith 
959cab54364SBarry Smith     Level: advanced
960cab54364SBarry Smith 
9619d90f715SBarry Smith     Notes:
96220662ed9SBarry Smith     Either `nout` or `idxout` may be `NULL`. `idx` and `idxout` may be identical.
9639d90f715SBarry Smith 
96420662ed9SBarry Smith     For "small" problems when using `ISGlobalToLocalMappingApply()` and `ISGlobalToLocalMappingApplyBlock()`, the `ISLocalToGlobalMappingType` of
96520662ed9SBarry Smith     `ISLOCALTOGLOBALMAPPINGBASIC` will be used;
966cab54364SBarry 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.
967cab54364SBarry Smith     Use `ISLocalToGlobalMappingSetType()` or call `ISLocalToGlobalMappingSetFromOptions()` with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
9689d90f715SBarry Smith 
969cab54364SBarry Smith     Developer Note:
97020662ed9SBarry Smith     The manual page states that `idx` and `idxout` may be identical but the calling
97120662ed9SBarry Smith     sequence declares `idx` as const so it cannot be the same as `idxout`.
9729d90f715SBarry Smith 
973cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApply()`, `ISGlobalToLocalMappingApply()`, `ISLocalToGlobalMappingCreate()`,
974db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`
9759d90f715SBarry Smith @*/
976d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, PetscInt n, const PetscInt idx[], PetscInt *nout, PetscInt idxout[])
977d71ae5a4SJacob Faibussowitsch {
9783a40ed3dSBarry Smith   PetscFunctionBegin;
9790700a824SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
98048a46eb9SPierre Jolivet   if (!mapping->data) PetscCall(ISGlobalToLocalMappingSetUp(mapping));
981dbbe0bcdSBarry Smith   PetscUseTypeMethod(mapping, globaltolocalmappingapplyblock, type, n, idx, nout, idxout);
9823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
983d4bb536fSBarry Smith }
98490f02eecSBarry Smith 
98589d82c54SBarry Smith /*@C
9866a818285SBarry Smith     ISLocalToGlobalMappingGetBlockInfo - Gets the neighbor information for each processor and
98789d82c54SBarry Smith      each index shared by more than one processor
98889d82c54SBarry Smith 
989c3339decSBarry Smith     Collective
99089d82c54SBarry Smith 
991f899ff85SJose E. Roman     Input Parameter:
99289d82c54SBarry Smith .   mapping - the mapping from local to global indexing
99389d82c54SBarry Smith 
994d8d19677SJose E. Roman     Output Parameters:
99589d82c54SBarry Smith +   nproc - number of processors that are connected to this one
99689d82c54SBarry Smith .   proc - neighboring processors
99707b52d57SBarry Smith .   numproc - number of indices for each subdomain (processor)
9983463a7baSJed Brown -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
99989d82c54SBarry Smith 
100089d82c54SBarry Smith     Level: advanced
100189d82c54SBarry Smith 
10022cfcea29SBarry Smith     Fortran Usage:
1003cab54364SBarry Smith .vb
10042cfcea29SBarry Smith         PetscInt indices[nproc][numprocmax],ierr)
1005cab54364SBarry Smith         ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
1006cab54364SBarry Smith         ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
1007cab54364SBarry Smith .ve
1008cab54364SBarry Smith 
1009cab54364SBarry Smith    Fortran Note:
101020662ed9SBarry Smith         There is no `ISLocalToGlobalMappingRestoreInfo()` in Fortran. You must make sure that `procs`[], `numprocs`[] and
101120662ed9SBarry Smith         `indices`[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
10122cfcea29SBarry Smith 
1013cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1014db781477SPatrick Sanan           `ISLocalToGlobalMappingRestoreInfo()`
101589d82c54SBarry Smith @*/
1016d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[])
1017d71ae5a4SJacob Faibussowitsch {
1018268a049cSStefano Zampini   PetscFunctionBegin;
1019268a049cSStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
1020268a049cSStefano Zampini   if (mapping->info_cached) {
1021268a049cSStefano Zampini     *nproc    = mapping->info_nproc;
1022268a049cSStefano Zampini     *procs    = mapping->info_procs;
1023268a049cSStefano Zampini     *numprocs = mapping->info_numprocs;
1024268a049cSStefano Zampini     *indices  = mapping->info_indices;
1025268a049cSStefano Zampini   } else {
10269566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetBlockInfo_Private(mapping, nproc, procs, numprocs, indices));
1027268a049cSStefano Zampini   }
10283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1029268a049cSStefano Zampini }
1030268a049cSStefano Zampini 
1031d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[])
1032d71ae5a4SJacob Faibussowitsch {
103397f1f81fSBarry Smith   PetscMPIInt  size, rank, tag1, tag2, tag3, *len, *source, imdex;
103432dcc486SBarry Smith   PetscInt     i, n = mapping->n, Ng, ng, max = 0, *lindices = mapping->indices;
103532dcc486SBarry Smith   PetscInt    *nprocs, *owner, nsends, *sends, j, *starts, nmax, nrecvs, *recvs, proc;
1036c599c493SJunchao Zhang   PetscInt     cnt, scale, *ownedsenders, *nownedsenders, rstart;
103732dcc486SBarry Smith   PetscInt     node, nownedm, nt, *sends2, nsends2, *starts2, *lens2, *dest, nrecvs2, *starts3, *recvs2, k, *bprocs, *tmp;
103832dcc486SBarry Smith   PetscInt     first_procs, first_numprocs, *first_indices;
103989d82c54SBarry Smith   MPI_Request *recv_waits, *send_waits;
104030dcb7c9SBarry Smith   MPI_Status   recv_status, *send_status, *recv_statuses;
1041ce94432eSBarry Smith   MPI_Comm     comm;
1042ace3abfcSBarry Smith   PetscBool    debug = PETSC_FALSE;
104389d82c54SBarry Smith 
104489d82c54SBarry Smith   PetscFunctionBegin;
10459566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mapping, &comm));
10469566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
10479566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
104824cf384cSBarry Smith   if (size == 1) {
104924cf384cSBarry Smith     *nproc = 0;
10500298fd71SBarry Smith     *procs = NULL;
10519566063dSJacob Faibussowitsch     PetscCall(PetscNew(numprocs));
10521e2105dcSBarry Smith     (*numprocs)[0] = 0;
10539566063dSJacob Faibussowitsch     PetscCall(PetscNew(indices));
10540298fd71SBarry Smith     (*indices)[0] = NULL;
1055268a049cSStefano Zampini     /* save info for reuse */
1056268a049cSStefano Zampini     mapping->info_nproc    = *nproc;
1057268a049cSStefano Zampini     mapping->info_procs    = *procs;
1058268a049cSStefano Zampini     mapping->info_numprocs = *numprocs;
1059268a049cSStefano Zampini     mapping->info_indices  = *indices;
1060268a049cSStefano Zampini     mapping->info_cached   = PETSC_TRUE;
10613ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
106224cf384cSBarry Smith   }
106324cf384cSBarry Smith 
10649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)mapping)->options, NULL, "-islocaltoglobalmappinggetinfo_debug", &debug, NULL));
106507b52d57SBarry Smith 
10663677ff5aSBarry Smith   /*
10676a818285SBarry Smith     Notes on ISLocalToGlobalMappingGetBlockInfo
10683677ff5aSBarry Smith 
10693677ff5aSBarry Smith     globally owned node - the nodes that have been assigned to this processor in global
10703677ff5aSBarry Smith            numbering, just for this routine.
10713677ff5aSBarry Smith 
10723677ff5aSBarry Smith     nontrivial globally owned node - node assigned to this processor that is on a subdomain
10733677ff5aSBarry Smith            boundary (i.e. is has more than one local owner)
10743677ff5aSBarry Smith 
10753677ff5aSBarry Smith     locally owned node - node that exists on this processors subdomain
10763677ff5aSBarry Smith 
10773677ff5aSBarry Smith     nontrivial locally owned node - node that is not in the interior (i.e. has more than one
10783677ff5aSBarry Smith            local subdomain
10793677ff5aSBarry Smith   */
10809566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetNewTag((PetscObject)mapping, &tag1));
10819566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetNewTag((PetscObject)mapping, &tag2));
10829566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetNewTag((PetscObject)mapping, &tag3));
108389d82c54SBarry Smith 
108489d82c54SBarry Smith   for (i = 0; i < n; i++) {
108589d82c54SBarry Smith     if (lindices[i] > max) max = lindices[i];
108689d82c54SBarry Smith   }
10871c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&max, &Ng, 1, MPIU_INT, MPI_MAX, comm));
108878058e43SBarry Smith   Ng++;
10899566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
10909566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1091bc8ff85bSBarry Smith   scale = Ng / size + 1;
10929371c9d4SSatish Balay   ng    = scale;
10939371c9d4SSatish Balay   if (rank == size - 1) ng = Ng - scale * (size - 1);
10949371c9d4SSatish Balay   ng     = PetscMax(1, ng);
1095caba0dd0SBarry Smith   rstart = scale * rank;
109689d82c54SBarry Smith 
109789d82c54SBarry Smith   /* determine ownership ranges of global indices */
10989566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * size, &nprocs));
10999566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(nprocs, 2 * size));
110089d82c54SBarry Smith 
110189d82c54SBarry Smith   /* determine owners of each local node  */
11029566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &owner));
110389d82c54SBarry Smith   for (i = 0; i < n; i++) {
11043677ff5aSBarry Smith     proc                 = lindices[i] / scale; /* processor that globally owns this index */
110527c402fcSBarry Smith     nprocs[2 * proc + 1] = 1;                   /* processor globally owns at least one of ours */
11063677ff5aSBarry Smith     owner[i]             = proc;
110727c402fcSBarry Smith     nprocs[2 * proc]++; /* count of how many that processor globally owns of ours */
110889d82c54SBarry Smith   }
11099371c9d4SSatish Balay   nsends = 0;
11109371c9d4SSatish Balay   for (i = 0; i < size; i++) nsends += nprocs[2 * i + 1];
11119566063dSJacob Faibussowitsch   PetscCall(PetscInfo(mapping, "Number of global owners for my local data %" PetscInt_FMT "\n", nsends));
111289d82c54SBarry Smith 
111389d82c54SBarry Smith   /* inform other processors of number of messages and max length*/
11149566063dSJacob Faibussowitsch   PetscCall(PetscMaxSum(comm, nprocs, &nmax, &nrecvs));
11159566063dSJacob Faibussowitsch   PetscCall(PetscInfo(mapping, "Number of local owners for my global data %" PetscInt_FMT "\n", nrecvs));
111689d82c54SBarry Smith 
111789d82c54SBarry Smith   /* post receives for owned rows */
11189566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((2 * nrecvs + 1) * (nmax + 1), &recvs));
11199566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs + 1, &recv_waits));
112048a46eb9SPierre 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));
112189d82c54SBarry Smith 
112289d82c54SBarry Smith   /* pack messages containing lists of local nodes to owners */
11239566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * n + 1, &sends));
11249566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size + 1, &starts));
112589d82c54SBarry Smith   starts[0] = 0;
1126f6e5521dSKarl Rupp   for (i = 1; i < size; i++) starts[i] = starts[i - 1] + 2 * nprocs[2 * i - 2];
112789d82c54SBarry Smith   for (i = 0; i < n; i++) {
112889d82c54SBarry Smith     sends[starts[owner[i]]++] = lindices[i];
112930dcb7c9SBarry Smith     sends[starts[owner[i]]++] = i;
113089d82c54SBarry Smith   }
11319566063dSJacob Faibussowitsch   PetscCall(PetscFree(owner));
113289d82c54SBarry Smith   starts[0] = 0;
1133f6e5521dSKarl Rupp   for (i = 1; i < size; i++) starts[i] = starts[i - 1] + 2 * nprocs[2 * i - 2];
113489d82c54SBarry Smith 
113589d82c54SBarry Smith   /* send the messages */
11369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends + 1, &send_waits));
11379566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends + 1, &dest));
113889d82c54SBarry Smith   cnt = 0;
113989d82c54SBarry Smith   for (i = 0; i < size; i++) {
114027c402fcSBarry Smith     if (nprocs[2 * i]) {
11419566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Isend(sends + starts[i], 2 * nprocs[2 * i], MPIU_INT, i, tag1, comm, send_waits + cnt));
114230dcb7c9SBarry Smith       dest[cnt] = i;
114389d82c54SBarry Smith       cnt++;
114489d82c54SBarry Smith     }
114589d82c54SBarry Smith   }
11469566063dSJacob Faibussowitsch   PetscCall(PetscFree(starts));
114789d82c54SBarry Smith 
114889d82c54SBarry Smith   /* wait on receives */
11499566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs + 1, &source));
11509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs + 1, &len));
115189d82c54SBarry Smith   cnt = nrecvs;
11529566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(ng + 1, &nownedsenders));
115389d82c54SBarry Smith   while (cnt) {
11549566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitany(nrecvs, recv_waits, &imdex, &recv_status));
115589d82c54SBarry Smith     /* unpack receives into our local space */
11569566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Get_count(&recv_status, MPIU_INT, &len[imdex]));
115789d82c54SBarry Smith     source[imdex] = recv_status.MPI_SOURCE;
115830dcb7c9SBarry Smith     len[imdex]    = len[imdex] / 2;
1159caba0dd0SBarry Smith     /* count how many local owners for each of my global owned indices */
116030dcb7c9SBarry Smith     for (i = 0; i < len[imdex]; i++) nownedsenders[recvs[2 * imdex * nmax + 2 * i] - rstart]++;
116189d82c54SBarry Smith     cnt--;
116289d82c54SBarry Smith   }
11639566063dSJacob Faibussowitsch   PetscCall(PetscFree(recv_waits));
116489d82c54SBarry Smith 
116530dcb7c9SBarry Smith   /* count how many globally owned indices are on an edge multiplied by how many processors own them. */
1166bc8ff85bSBarry Smith   nownedm = 0;
1167bc8ff85bSBarry Smith   for (i = 0; i < ng; i++) {
1168c599c493SJunchao Zhang     if (nownedsenders[i] > 1) nownedm += nownedsenders[i];
1169bc8ff85bSBarry Smith   }
1170bc8ff85bSBarry Smith 
1171bc8ff85bSBarry Smith   /* create single array to contain rank of all local owners of each globally owned index */
11729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nownedm + 1, &ownedsenders));
11739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ng + 1, &starts));
1174bc8ff85bSBarry Smith   starts[0] = 0;
1175bc8ff85bSBarry Smith   for (i = 1; i < ng; i++) {
1176bc8ff85bSBarry Smith     if (nownedsenders[i - 1] > 1) starts[i] = starts[i - 1] + nownedsenders[i - 1];
1177bc8ff85bSBarry Smith     else starts[i] = starts[i - 1];
1178bc8ff85bSBarry Smith   }
1179bc8ff85bSBarry Smith 
11806aad120cSJose E. Roman   /* for each nontrivial globally owned node list all arriving processors */
1181bc8ff85bSBarry Smith   for (i = 0; i < nrecvs; i++) {
1182bc8ff85bSBarry Smith     for (j = 0; j < len[i]; j++) {
118330dcb7c9SBarry Smith       node = recvs[2 * i * nmax + 2 * j] - rstart;
1184f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) ownedsenders[starts[node]++] = source[i];
1185bc8ff85bSBarry Smith     }
1186bc8ff85bSBarry Smith   }
1187bc8ff85bSBarry Smith 
118807b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
118930dcb7c9SBarry Smith     starts[0] = 0;
119030dcb7c9SBarry Smith     for (i = 1; i < ng; i++) {
119130dcb7c9SBarry Smith       if (nownedsenders[i - 1] > 1) starts[i] = starts[i - 1] + nownedsenders[i - 1];
119230dcb7c9SBarry Smith       else starts[i] = starts[i - 1];
119330dcb7c9SBarry Smith     }
119430dcb7c9SBarry Smith     for (i = 0; i < ng; i++) {
119530dcb7c9SBarry Smith       if (nownedsenders[i] > 1) {
11969566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(comm, "[%d] global node %" PetscInt_FMT " local owner processors: ", rank, i + rstart));
119748a46eb9SPierre Jolivet         for (j = 0; j < nownedsenders[i]; j++) PetscCall(PetscSynchronizedPrintf(comm, "%" PetscInt_FMT " ", ownedsenders[starts[i] + j]));
11989566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(comm, "\n"));
119930dcb7c9SBarry Smith       }
120030dcb7c9SBarry Smith     }
12019566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
120207b52d57SBarry Smith   } /* -----------------------------------  */
120330dcb7c9SBarry Smith 
12043677ff5aSBarry Smith   /* wait on original sends */
12053a96401aSBarry Smith   if (nsends) {
12069566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nsends, &send_status));
12079566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nsends, send_waits, send_status));
12089566063dSJacob Faibussowitsch     PetscCall(PetscFree(send_status));
12093a96401aSBarry Smith   }
12109566063dSJacob Faibussowitsch   PetscCall(PetscFree(send_waits));
12119566063dSJacob Faibussowitsch   PetscCall(PetscFree(sends));
12129566063dSJacob Faibussowitsch   PetscCall(PetscFree(nprocs));
12133677ff5aSBarry Smith 
12143677ff5aSBarry Smith   /* pack messages to send back to local owners */
121530dcb7c9SBarry Smith   starts[0] = 0;
121630dcb7c9SBarry Smith   for (i = 1; i < ng; i++) {
121730dcb7c9SBarry Smith     if (nownedsenders[i - 1] > 1) starts[i] = starts[i - 1] + nownedsenders[i - 1];
121830dcb7c9SBarry Smith     else starts[i] = starts[i - 1];
121930dcb7c9SBarry Smith   }
122030dcb7c9SBarry Smith   nsends2 = nrecvs;
12219566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends2 + 1, &nprocs)); /* length of each message */
122230dcb7c9SBarry Smith   for (i = 0; i < nrecvs; i++) {
122330dcb7c9SBarry Smith     nprocs[i] = 1;
122430dcb7c9SBarry Smith     for (j = 0; j < len[i]; j++) {
122530dcb7c9SBarry Smith       node = recvs[2 * i * nmax + 2 * j] - rstart;
1226f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) nprocs[i] += 2 + nownedsenders[node];
122730dcb7c9SBarry Smith     }
122830dcb7c9SBarry Smith   }
1229f6e5521dSKarl Rupp   nt = 0;
1230f6e5521dSKarl Rupp   for (i = 0; i < nsends2; i++) nt += nprocs[i];
1231f6e5521dSKarl Rupp 
12329566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt + 1, &sends2));
12339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends2 + 1, &starts2));
1234f6e5521dSKarl Rupp 
1235f6e5521dSKarl Rupp   starts2[0] = 0;
1236f6e5521dSKarl Rupp   for (i = 1; i < nsends2; i++) starts2[i] = starts2[i - 1] + nprocs[i - 1];
123730dcb7c9SBarry Smith   /*
123830dcb7c9SBarry Smith      Each message is 1 + nprocs[i] long, and consists of
123930dcb7c9SBarry Smith        (0) the number of nodes being sent back
124030dcb7c9SBarry Smith        (1) the local node number,
124130dcb7c9SBarry Smith        (2) the number of processors sharing it,
124230dcb7c9SBarry Smith        (3) the processors sharing it
124330dcb7c9SBarry Smith   */
124430dcb7c9SBarry Smith   for (i = 0; i < nsends2; i++) {
124530dcb7c9SBarry Smith     cnt                = 1;
124630dcb7c9SBarry Smith     sends2[starts2[i]] = 0;
124730dcb7c9SBarry Smith     for (j = 0; j < len[i]; j++) {
124830dcb7c9SBarry Smith       node = recvs[2 * i * nmax + 2 * j] - rstart;
124930dcb7c9SBarry Smith       if (nownedsenders[node] > 1) {
125030dcb7c9SBarry Smith         sends2[starts2[i]]++;
125130dcb7c9SBarry Smith         sends2[starts2[i] + cnt++] = recvs[2 * i * nmax + 2 * j + 1];
125230dcb7c9SBarry Smith         sends2[starts2[i] + cnt++] = nownedsenders[node];
12539566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(&sends2[starts2[i] + cnt], &ownedsenders[starts[node]], nownedsenders[node]));
125430dcb7c9SBarry Smith         cnt += nownedsenders[node];
125530dcb7c9SBarry Smith       }
125630dcb7c9SBarry Smith     }
125730dcb7c9SBarry Smith   }
125830dcb7c9SBarry Smith 
125930dcb7c9SBarry Smith   /* receive the message lengths */
126030dcb7c9SBarry Smith   nrecvs2 = nsends;
12619566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs2 + 1, &lens2));
12629566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs2 + 1, &starts3));
12639566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs2 + 1, &recv_waits));
126448a46eb9SPierre Jolivet   for (i = 0; i < nrecvs2; i++) PetscCallMPI(MPI_Irecv(&lens2[i], 1, MPIU_INT, dest[i], tag2, comm, recv_waits + i));
1265d44834fbSBarry Smith 
12668a8e0b3aSBarry Smith   /* send the message lengths */
126748a46eb9SPierre Jolivet   for (i = 0; i < nsends2; i++) PetscCallMPI(MPI_Send(&nprocs[i], 1, MPIU_INT, source[i], tag2, comm));
12688a8e0b3aSBarry Smith 
1269d44834fbSBarry Smith   /* wait on receives of lens */
12700c468ba9SBarry Smith   if (nrecvs2) {
12719566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrecvs2, &recv_statuses));
12729566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrecvs2, recv_waits, recv_statuses));
12739566063dSJacob Faibussowitsch     PetscCall(PetscFree(recv_statuses));
12740c468ba9SBarry Smith   }
12759566063dSJacob Faibussowitsch   PetscCall(PetscFree(recv_waits));
1276d44834fbSBarry Smith 
127730dcb7c9SBarry Smith   starts3[0] = 0;
1278d44834fbSBarry Smith   nt         = 0;
127930dcb7c9SBarry Smith   for (i = 0; i < nrecvs2 - 1; i++) {
128030dcb7c9SBarry Smith     starts3[i + 1] = starts3[i] + lens2[i];
1281d44834fbSBarry Smith     nt += lens2[i];
128230dcb7c9SBarry Smith   }
128376466f69SStefano Zampini   if (nrecvs2) nt += lens2[nrecvs2 - 1];
1284d44834fbSBarry Smith 
12859566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt + 1, &recvs2));
12869566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs2 + 1, &recv_waits));
128748a46eb9SPierre Jolivet   for (i = 0; i < nrecvs2; i++) PetscCallMPI(MPI_Irecv(recvs2 + starts3[i], lens2[i], MPIU_INT, dest[i], tag3, comm, recv_waits + i));
128830dcb7c9SBarry Smith 
128930dcb7c9SBarry Smith   /* send the messages */
12909566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends2 + 1, &send_waits));
129148a46eb9SPierre Jolivet   for (i = 0; i < nsends2; i++) PetscCallMPI(MPI_Isend(sends2 + starts2[i], nprocs[i], MPIU_INT, source[i], tag3, comm, send_waits + i));
129230dcb7c9SBarry Smith 
129330dcb7c9SBarry Smith   /* wait on receives */
12940c468ba9SBarry Smith   if (nrecvs2) {
12959566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrecvs2, &recv_statuses));
12969566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrecvs2, recv_waits, recv_statuses));
12979566063dSJacob Faibussowitsch     PetscCall(PetscFree(recv_statuses));
12980c468ba9SBarry Smith   }
12999566063dSJacob Faibussowitsch   PetscCall(PetscFree(recv_waits));
13009566063dSJacob Faibussowitsch   PetscCall(PetscFree(nprocs));
130130dcb7c9SBarry Smith 
130207b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
130330dcb7c9SBarry Smith     cnt = 0;
130430dcb7c9SBarry Smith     for (i = 0; i < nrecvs2; i++) {
130530dcb7c9SBarry Smith       nt = recvs2[cnt++];
130630dcb7c9SBarry Smith       for (j = 0; j < nt; j++) {
13079566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(comm, "[%d] local node %" PetscInt_FMT " number of subdomains %" PetscInt_FMT ": ", rank, recvs2[cnt], recvs2[cnt + 1]));
130848a46eb9SPierre Jolivet         for (k = 0; k < recvs2[cnt + 1]; k++) PetscCall(PetscSynchronizedPrintf(comm, "%" PetscInt_FMT " ", recvs2[cnt + 2 + k]));
130930dcb7c9SBarry Smith         cnt += 2 + recvs2[cnt + 1];
13109566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(comm, "\n"));
131130dcb7c9SBarry Smith       }
131230dcb7c9SBarry Smith     }
13139566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
131407b52d57SBarry Smith   } /* -----------------------------------  */
131530dcb7c9SBarry Smith 
131630dcb7c9SBarry Smith   /* count number subdomains for each local node */
13179566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(size, &nprocs));
131830dcb7c9SBarry Smith   cnt = 0;
131930dcb7c9SBarry Smith   for (i = 0; i < nrecvs2; i++) {
132030dcb7c9SBarry Smith     nt = recvs2[cnt++];
132130dcb7c9SBarry Smith     for (j = 0; j < nt; j++) {
1322f6e5521dSKarl Rupp       for (k = 0; k < recvs2[cnt + 1]; k++) nprocs[recvs2[cnt + 2 + k]]++;
132330dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt + 1];
132430dcb7c9SBarry Smith     }
132530dcb7c9SBarry Smith   }
13269371c9d4SSatish Balay   nt = 0;
13279371c9d4SSatish Balay   for (i = 0; i < size; i++) nt += (nprocs[i] > 0);
132830dcb7c9SBarry Smith   *nproc = nt;
13299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt + 1, procs));
13309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt + 1, numprocs));
13319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt + 1, indices));
13320298fd71SBarry Smith   for (i = 0; i < nt + 1; i++) (*indices)[i] = NULL;
13339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size, &bprocs));
133430dcb7c9SBarry Smith   cnt = 0;
133530dcb7c9SBarry Smith   for (i = 0; i < size; i++) {
133630dcb7c9SBarry Smith     if (nprocs[i] > 0) {
133730dcb7c9SBarry Smith       bprocs[i]        = cnt;
133830dcb7c9SBarry Smith       (*procs)[cnt]    = i;
133930dcb7c9SBarry Smith       (*numprocs)[cnt] = nprocs[i];
13409566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nprocs[i], &(*indices)[cnt]));
134130dcb7c9SBarry Smith       cnt++;
134230dcb7c9SBarry Smith     }
134330dcb7c9SBarry Smith   }
134430dcb7c9SBarry Smith 
134530dcb7c9SBarry Smith   /* make the list of subdomains for each nontrivial local node */
13469566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(*numprocs, nt));
134730dcb7c9SBarry Smith   cnt = 0;
134830dcb7c9SBarry Smith   for (i = 0; i < nrecvs2; i++) {
134930dcb7c9SBarry Smith     nt = recvs2[cnt++];
135030dcb7c9SBarry Smith     for (j = 0; j < nt; j++) {
1351f6e5521dSKarl Rupp       for (k = 0; k < recvs2[cnt + 1]; k++) (*indices)[bprocs[recvs2[cnt + 2 + k]]][(*numprocs)[bprocs[recvs2[cnt + 2 + k]]]++] = recvs2[cnt];
135230dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt + 1];
135330dcb7c9SBarry Smith     }
135430dcb7c9SBarry Smith   }
13559566063dSJacob Faibussowitsch   PetscCall(PetscFree(bprocs));
13569566063dSJacob Faibussowitsch   PetscCall(PetscFree(recvs2));
135730dcb7c9SBarry Smith 
135807b52d57SBarry Smith   /* sort the node indexing by their global numbers */
135907b52d57SBarry Smith   nt = *nproc;
136007b52d57SBarry Smith   for (i = 0; i < nt; i++) {
13619566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1((*numprocs)[i], &tmp));
1362f6e5521dSKarl Rupp     for (j = 0; j < (*numprocs)[i]; j++) tmp[j] = lindices[(*indices)[i][j]];
13639566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithArray((*numprocs)[i], tmp, (*indices)[i]));
13649566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmp));
136507b52d57SBarry Smith   }
136607b52d57SBarry Smith 
136707b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
136830dcb7c9SBarry Smith     nt = *nproc;
136930dcb7c9SBarry Smith     for (i = 0; i < nt; i++) {
13709566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedPrintf(comm, "[%d] subdomain %" PetscInt_FMT " number of indices %" PetscInt_FMT ": ", rank, (*procs)[i], (*numprocs)[i]));
137148a46eb9SPierre Jolivet       for (j = 0; j < (*numprocs)[i]; j++) PetscCall(PetscSynchronizedPrintf(comm, "%" PetscInt_FMT " ", (*indices)[i][j]));
13729566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedPrintf(comm, "\n"));
137330dcb7c9SBarry Smith     }
13749566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
137507b52d57SBarry Smith   } /* -----------------------------------  */
137630dcb7c9SBarry Smith 
137730dcb7c9SBarry Smith   /* wait on sends */
137830dcb7c9SBarry Smith   if (nsends2) {
13799566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nsends2, &send_status));
13809566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nsends2, send_waits, send_status));
13819566063dSJacob Faibussowitsch     PetscCall(PetscFree(send_status));
138230dcb7c9SBarry Smith   }
138330dcb7c9SBarry Smith 
13849566063dSJacob Faibussowitsch   PetscCall(PetscFree(starts3));
13859566063dSJacob Faibussowitsch   PetscCall(PetscFree(dest));
13869566063dSJacob Faibussowitsch   PetscCall(PetscFree(send_waits));
13873677ff5aSBarry Smith 
13889566063dSJacob Faibussowitsch   PetscCall(PetscFree(nownedsenders));
13899566063dSJacob Faibussowitsch   PetscCall(PetscFree(ownedsenders));
13909566063dSJacob Faibussowitsch   PetscCall(PetscFree(starts));
13919566063dSJacob Faibussowitsch   PetscCall(PetscFree(starts2));
13929566063dSJacob Faibussowitsch   PetscCall(PetscFree(lens2));
139389d82c54SBarry Smith 
13949566063dSJacob Faibussowitsch   PetscCall(PetscFree(source));
13959566063dSJacob Faibussowitsch   PetscCall(PetscFree(len));
13969566063dSJacob Faibussowitsch   PetscCall(PetscFree(recvs));
13979566063dSJacob Faibussowitsch   PetscCall(PetscFree(nprocs));
13989566063dSJacob Faibussowitsch   PetscCall(PetscFree(sends2));
139924cf384cSBarry Smith 
140024cf384cSBarry Smith   /* put the information about myself as the first entry in the list */
140124cf384cSBarry Smith   first_procs    = (*procs)[0];
140224cf384cSBarry Smith   first_numprocs = (*numprocs)[0];
140324cf384cSBarry Smith   first_indices  = (*indices)[0];
140424cf384cSBarry Smith   for (i = 0; i < *nproc; i++) {
140524cf384cSBarry Smith     if ((*procs)[i] == rank) {
140624cf384cSBarry Smith       (*procs)[0]    = (*procs)[i];
140724cf384cSBarry Smith       (*numprocs)[0] = (*numprocs)[i];
140824cf384cSBarry Smith       (*indices)[0]  = (*indices)[i];
140924cf384cSBarry Smith       (*procs)[i]    = first_procs;
141024cf384cSBarry Smith       (*numprocs)[i] = first_numprocs;
141124cf384cSBarry Smith       (*indices)[i]  = first_indices;
141224cf384cSBarry Smith       break;
141324cf384cSBarry Smith     }
141424cf384cSBarry Smith   }
1415268a049cSStefano Zampini 
1416268a049cSStefano Zampini   /* save info for reuse */
1417268a049cSStefano Zampini   mapping->info_nproc    = *nproc;
1418268a049cSStefano Zampini   mapping->info_procs    = *procs;
1419268a049cSStefano Zampini   mapping->info_numprocs = *numprocs;
1420268a049cSStefano Zampini   mapping->info_indices  = *indices;
1421268a049cSStefano Zampini   mapping->info_cached   = PETSC_TRUE;
14223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
142389d82c54SBarry Smith }
142489d82c54SBarry Smith 
14256a818285SBarry Smith /*@C
1426cab54364SBarry Smith     ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by `ISLocalToGlobalMappingGetBlockInfo()`
14276a818285SBarry Smith 
1428c3339decSBarry Smith     Collective
14296a818285SBarry Smith 
1430f899ff85SJose E. Roman     Input Parameter:
14316a818285SBarry Smith .   mapping - the mapping from local to global indexing
14326a818285SBarry Smith 
1433d8d19677SJose E. Roman     Output Parameters:
14346a818285SBarry Smith +   nproc - number of processors that are connected to this one
14356a818285SBarry Smith .   proc - neighboring processors
14366a818285SBarry Smith .   numproc - number of indices for each processor
14376a818285SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
14386a818285SBarry Smith 
14396a818285SBarry Smith     Level: advanced
14406a818285SBarry Smith 
1441cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1442db781477SPatrick Sanan           `ISLocalToGlobalMappingGetInfo()`
14436a818285SBarry Smith @*/
1444d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[])
1445d71ae5a4SJacob Faibussowitsch {
14466a818285SBarry Smith   PetscFunctionBegin;
1447cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
1448268a049cSStefano Zampini   if (mapping->info_free) {
14499566063dSJacob Faibussowitsch     PetscCall(PetscFree(*numprocs));
14506a818285SBarry Smith     if (*indices) {
1451268a049cSStefano Zampini       PetscInt i;
1452268a049cSStefano Zampini 
14539566063dSJacob Faibussowitsch       PetscCall(PetscFree((*indices)[0]));
145448a46eb9SPierre Jolivet       for (i = 1; i < *nproc; i++) PetscCall(PetscFree((*indices)[i]));
14559566063dSJacob Faibussowitsch       PetscCall(PetscFree(*indices));
14566a818285SBarry Smith     }
1457268a049cSStefano Zampini   }
1458268a049cSStefano Zampini   *nproc    = 0;
1459268a049cSStefano Zampini   *procs    = NULL;
1460268a049cSStefano Zampini   *numprocs = NULL;
1461268a049cSStefano Zampini   *indices  = NULL;
14623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14636a818285SBarry Smith }
14646a818285SBarry Smith 
14656a818285SBarry Smith /*@C
14666a818285SBarry Smith     ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each processor and
14676a818285SBarry Smith      each index shared by more than one processor
14686a818285SBarry Smith 
1469c3339decSBarry Smith     Collective
14706a818285SBarry Smith 
1471f899ff85SJose E. Roman     Input Parameter:
14726a818285SBarry Smith .   mapping - the mapping from local to global indexing
14736a818285SBarry Smith 
1474d8d19677SJose E. Roman     Output Parameters:
14756a818285SBarry Smith +   nproc - number of processors that are connected to this one
14766a818285SBarry Smith .   proc - neighboring processors
14776a818285SBarry Smith .   numproc - number of indices for each subdomain (processor)
14786a818285SBarry Smith -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
14796a818285SBarry Smith 
14806a818285SBarry Smith     Level: advanced
14816a818285SBarry Smith 
1482cab54364SBarry Smith     Note:
1483cab54364SBarry Smith     The user needs to call `ISLocalToGlobalMappingRestoreInfo()` when the data is no longer needed.
14841bd0b88eSStefano Zampini 
14856a818285SBarry Smith     Fortran Usage:
1486cab54364SBarry Smith .vb
14876a818285SBarry Smith         PetscInt indices[nproc][numprocmax],ierr)
1488cab54364SBarry Smith         ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
1489cab54364SBarry Smith         ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
1490cab54364SBarry Smith .ve
1491cab54364SBarry Smith 
1492cab54364SBarry Smith     Fortran Note:
149320662ed9SBarry Smith         There is no `ISLocalToGlobalMappingRestoreInfo()` in Fortran. You must make sure that `procs`[], `numprocs`[] and
149420662ed9SBarry Smith         `indices`[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
14956a818285SBarry Smith 
1496cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1497db781477SPatrick Sanan           `ISLocalToGlobalMappingRestoreInfo()`
14986a818285SBarry Smith @*/
1499d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[])
1500d71ae5a4SJacob Faibussowitsch {
15018a1f772fSStefano Zampini   PetscInt **bindices = NULL, *bnumprocs = NULL, bs, i, j, k;
15026a818285SBarry Smith 
15036a818285SBarry Smith   PetscFunctionBegin;
15046a818285SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
15058a1f772fSStefano Zampini   bs = mapping->bs;
15069566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockInfo(mapping, nproc, procs, &bnumprocs, &bindices));
1507268a049cSStefano Zampini   if (bs > 1) { /* we need to expand the cached info */
15089566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(*nproc, &*indices));
15099566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(*nproc, &*numprocs));
15106a818285SBarry Smith     for (i = 0; i < *nproc; i++) {
15119566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(bs * bnumprocs[i], &(*indices)[i]));
1512268a049cSStefano Zampini       for (j = 0; j < bnumprocs[i]; j++) {
1513ad540459SPierre Jolivet         for (k = 0; k < bs; k++) (*indices)[i][j * bs + k] = bs * bindices[i][j] + k;
15146a818285SBarry Smith       }
1515268a049cSStefano Zampini       (*numprocs)[i] = bnumprocs[i] * bs;
15166a818285SBarry Smith     }
1517268a049cSStefano Zampini     mapping->info_free = PETSC_TRUE;
1518268a049cSStefano Zampini   } else {
1519268a049cSStefano Zampini     *numprocs = bnumprocs;
1520268a049cSStefano Zampini     *indices  = bindices;
15216a818285SBarry Smith   }
15223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15236a818285SBarry Smith }
15246a818285SBarry Smith 
152507b52d57SBarry Smith /*@C
1526cab54364SBarry Smith     ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by `ISLocalToGlobalMappingGetInfo()`
152789d82c54SBarry Smith 
1528c3339decSBarry Smith     Collective
152907b52d57SBarry Smith 
1530f899ff85SJose E. Roman     Input Parameter:
153107b52d57SBarry Smith .   mapping - the mapping from local to global indexing
153207b52d57SBarry Smith 
1533d8d19677SJose E. Roman     Output Parameters:
153407b52d57SBarry Smith +   nproc - number of processors that are connected to this one
153507b52d57SBarry Smith .   proc - neighboring processors
153607b52d57SBarry Smith .   numproc - number of indices for each processor
153707b52d57SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
153807b52d57SBarry Smith 
153907b52d57SBarry Smith     Level: advanced
154007b52d57SBarry Smith 
1541cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1542db781477SPatrick Sanan           `ISLocalToGlobalMappingGetInfo()`
154307b52d57SBarry Smith @*/
1544d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[])
1545d71ae5a4SJacob Faibussowitsch {
154607b52d57SBarry Smith   PetscFunctionBegin;
15479566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockInfo(mapping, nproc, procs, numprocs, indices));
15483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
154907b52d57SBarry Smith }
155086994e45SJed Brown 
155186994e45SJed Brown /*@C
1552cab54364SBarry Smith     ISLocalToGlobalMappingGetNodeInfo - Gets the neighbor information for each MPI rank
15531bd0b88eSStefano Zampini 
1554c3339decSBarry Smith     Collective
15551bd0b88eSStefano Zampini 
1556f899ff85SJose E. Roman     Input Parameter:
15571bd0b88eSStefano Zampini .   mapping - the mapping from local to global indexing
15581bd0b88eSStefano Zampini 
1559d8d19677SJose E. Roman     Output Parameters:
1560cab54364SBarry Smith +   nnodes - number of local nodes (same `ISLocalToGlobalMappingGetSize()`)
15611bd0b88eSStefano Zampini .   count - number of neighboring processors per node
15621bd0b88eSStefano Zampini -   indices - indices of processes sharing the node (sorted)
15631bd0b88eSStefano Zampini 
15641bd0b88eSStefano Zampini     Level: advanced
15651bd0b88eSStefano Zampini 
1566cab54364SBarry Smith     Note:
1567cab54364SBarry Smith     The user needs to call `ISLocalToGlobalMappingRestoreInfo()` when the data is no longer needed.
15681bd0b88eSStefano Zampini 
1569cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1570db781477SPatrick Sanan           `ISLocalToGlobalMappingGetInfo()`, `ISLocalToGlobalMappingRestoreNodeInfo()`
15711bd0b88eSStefano Zampini @*/
1572d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetNodeInfo(ISLocalToGlobalMapping mapping, PetscInt *nnodes, PetscInt *count[], PetscInt **indices[])
1573d71ae5a4SJacob Faibussowitsch {
15741bd0b88eSStefano Zampini   PetscInt n;
15751bd0b88eSStefano Zampini 
15761bd0b88eSStefano Zampini   PetscFunctionBegin;
15771bd0b88eSStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
15789566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetSize(mapping, &n));
15791bd0b88eSStefano Zampini   if (!mapping->info_nodec) {
15801bd0b88eSStefano Zampini     PetscInt i, m, n_neigh, *neigh, *n_shared, **shared;
15811bd0b88eSStefano Zampini 
15829566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(n + 1, &mapping->info_nodec, n, &mapping->info_nodei));
15839566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetInfo(mapping, &n_neigh, &neigh, &n_shared, &shared));
1584ad540459SPierre Jolivet     for (i = 0; i < n; i++) mapping->info_nodec[i] = 1;
1585071fcb05SBarry Smith     m                      = n;
1586071fcb05SBarry Smith     mapping->info_nodec[n] = 0;
15871bd0b88eSStefano Zampini     for (i = 1; i < n_neigh; i++) {
15881bd0b88eSStefano Zampini       PetscInt j;
15891bd0b88eSStefano Zampini 
15901bd0b88eSStefano Zampini       m += n_shared[i];
15911bd0b88eSStefano Zampini       for (j = 0; j < n_shared[i]; j++) mapping->info_nodec[shared[i][j]] += 1;
15921bd0b88eSStefano Zampini     }
15939566063dSJacob Faibussowitsch     if (n) PetscCall(PetscMalloc1(m, &mapping->info_nodei[0]));
15941bd0b88eSStefano Zampini     for (i = 1; i < n; i++) mapping->info_nodei[i] = mapping->info_nodei[i - 1] + mapping->info_nodec[i - 1];
15959566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(mapping->info_nodec, n));
15969371c9d4SSatish Balay     for (i = 0; i < n; i++) {
15979371c9d4SSatish Balay       mapping->info_nodec[i]    = 1;
15989371c9d4SSatish Balay       mapping->info_nodei[i][0] = neigh[0];
15999371c9d4SSatish Balay     }
16001bd0b88eSStefano Zampini     for (i = 1; i < n_neigh; i++) {
16011bd0b88eSStefano Zampini       PetscInt j;
16021bd0b88eSStefano Zampini 
16031bd0b88eSStefano Zampini       for (j = 0; j < n_shared[i]; j++) {
16041bd0b88eSStefano Zampini         PetscInt k = shared[i][j];
16051bd0b88eSStefano Zampini 
16061bd0b88eSStefano Zampini         mapping->info_nodei[k][mapping->info_nodec[k]] = neigh[i];
16071bd0b88eSStefano Zampini         mapping->info_nodec[k] += 1;
16081bd0b88eSStefano Zampini       }
16091bd0b88eSStefano Zampini     }
16109566063dSJacob Faibussowitsch     for (i = 0; i < n; i++) PetscCall(PetscSortRemoveDupsInt(&mapping->info_nodec[i], mapping->info_nodei[i]));
16119566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreInfo(mapping, &n_neigh, &neigh, &n_shared, &shared));
16121bd0b88eSStefano Zampini   }
16131bd0b88eSStefano Zampini   if (nnodes) *nnodes = n;
16141bd0b88eSStefano Zampini   if (count) *count = mapping->info_nodec;
16151bd0b88eSStefano Zampini   if (indices) *indices = mapping->info_nodei;
16163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16171bd0b88eSStefano Zampini }
16181bd0b88eSStefano Zampini 
16191bd0b88eSStefano Zampini /*@C
1620cab54364SBarry Smith     ISLocalToGlobalMappingRestoreNodeInfo - Frees the memory allocated by `ISLocalToGlobalMappingGetNodeInfo()`
16211bd0b88eSStefano Zampini 
1622c3339decSBarry Smith     Collective
16231bd0b88eSStefano Zampini 
1624f899ff85SJose E. Roman     Input Parameter:
16251bd0b88eSStefano Zampini .   mapping - the mapping from local to global indexing
16261bd0b88eSStefano Zampini 
1627d8d19677SJose E. Roman     Output Parameters:
16281bd0b88eSStefano Zampini +   nnodes - number of local nodes
16291bd0b88eSStefano Zampini .   count - number of neighboring processors per node
16301bd0b88eSStefano Zampini -   indices - indices of processes sharing the node (sorted)
16311bd0b88eSStefano Zampini 
16321bd0b88eSStefano Zampini     Level: advanced
16331bd0b88eSStefano Zampini 
1634cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1635db781477SPatrick Sanan           `ISLocalToGlobalMappingGetInfo()`
16361bd0b88eSStefano Zampini @*/
1637d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreNodeInfo(ISLocalToGlobalMapping mapping, PetscInt *nnodes, PetscInt *count[], PetscInt **indices[])
1638d71ae5a4SJacob Faibussowitsch {
16391bd0b88eSStefano Zampini   PetscFunctionBegin;
16401bd0b88eSStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
16411bd0b88eSStefano Zampini   if (nnodes) *nnodes = 0;
16421bd0b88eSStefano Zampini   if (count) *count = NULL;
16431bd0b88eSStefano Zampini   if (indices) *indices = NULL;
16443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16451bd0b88eSStefano Zampini }
16461bd0b88eSStefano Zampini 
16476ce40d10SBarry Smith /*MC
16486ce40d10SBarry Smith     ISLocalToGlobalMappingGetIndicesF90 - Get global indices for every local point that is mapped
16496ce40d10SBarry Smith 
16506ce40d10SBarry Smith     Synopsis:
16516ce40d10SBarry Smith     ISLocalToGlobalMappingGetIndicesF90(ISLocalToGlobalMapping ltog, PetscInt, pointer :: array(:)}, integer ierr)
16526ce40d10SBarry Smith 
16536ce40d10SBarry Smith     Not Collective
16546ce40d10SBarry Smith 
16556ce40d10SBarry Smith     Input Parameter:
16566ce40d10SBarry Smith .   A - the matrix
16576ce40d10SBarry Smith 
1658*2fe279fdSBarry Smith     Output Parameter:
16596ce40d10SBarry Smith .   array - array of indices, the length of this array may be obtained with `ISLocalToGlobalMappingGetSize()`
16606ce40d10SBarry Smith 
16616ce40d10SBarry Smith     Level: advanced
16626ce40d10SBarry Smith 
16636ce40d10SBarry Smith     Note:
16646ce40d10SBarry Smith     Use  `ISLocalToGlobalMappingGetIndicesF90()` when you no longer need access to the data
16656ce40d10SBarry Smith 
166620662ed9SBarry Smith .seealso: [](sec_fortranarrays), [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingGetIndices()`, `ISLocalToGlobalMappingRestoreIndices()`,
166720662ed9SBarry Smith           `ISLocalToGlobalMappingRestoreIndicesF90()`
16686ce40d10SBarry Smith M*/
16696ce40d10SBarry Smith 
16706ce40d10SBarry Smith /*MC
16716ce40d10SBarry Smith     ISLocalToGlobalMappingRestoreIndicesF90 - restores the global indices for every local point that is mapped obtained with `ISLocalToGlobalMappingGetIndicesF90()`
16726ce40d10SBarry Smith 
16736ce40d10SBarry Smith     Synopsis:
16746ce40d10SBarry Smith     ISLocalToGlobalMappingRestoreIndicesF90(ISLocalToGlobalMapping ltog, PetscInt, pointer :: array(:)}, integer ierr)
16756ce40d10SBarry Smith 
16766ce40d10SBarry Smith     Not Collective
16776ce40d10SBarry Smith 
16786ce40d10SBarry Smith     Input Parameters:
16796ce40d10SBarry Smith +   A - the matrix
16806ce40d10SBarry Smith -   array - array of indices, the length of this array may be obtained with `ISLocalToGlobalMappingGetSize()`
16816ce40d10SBarry Smith 
16826ce40d10SBarry Smith     Level: advanced
16836ce40d10SBarry Smith 
168420662ed9SBarry Smith .seealso: [](sec_fortranarrays), [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingGetIndices()`, `ISLocalToGlobalMappingRestoreIndices()`,
168520662ed9SBarry Smith           `ISLocalToGlobalMappingGetIndicesF90()`
16866ce40d10SBarry Smith M*/
16876ce40d10SBarry Smith 
16881bd0b88eSStefano Zampini /*@C
1689107e9a97SBarry Smith    ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped
169086994e45SJed Brown 
169186994e45SJed Brown    Not Collective
169286994e45SJed Brown 
16934165533cSJose E. Roman    Input Parameter:
169486994e45SJed Brown . ltog - local to global mapping
169586994e45SJed Brown 
16964165533cSJose E. Roman    Output Parameter:
1697cab54364SBarry Smith . array - array of indices, the length of this array may be obtained with `ISLocalToGlobalMappingGetSize()`
169886994e45SJed Brown 
169986994e45SJed Brown    Level: advanced
170086994e45SJed Brown 
1701cab54364SBarry Smith    Note:
1702cab54364SBarry Smith     `ISLocalToGlobalMappingGetSize()` returns the length the this array
1703107e9a97SBarry Smith 
170420662ed9SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingRestoreIndices()`,
170520662ed9SBarry Smith           `ISLocalToGlobalMappingGetBlockIndices()`, `ISLocalToGlobalMappingRestoreBlockIndices()`
170686994e45SJed Brown @*/
1707d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog, const PetscInt **array)
1708d71ae5a4SJacob Faibussowitsch {
170986994e45SJed Brown   PetscFunctionBegin;
171086994e45SJed Brown   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
171186994e45SJed Brown   PetscValidPointer(array, 2);
171245b6f7e9SBarry Smith   if (ltog->bs == 1) {
171386994e45SJed Brown     *array = ltog->indices;
171445b6f7e9SBarry Smith   } else {
171545b6f7e9SBarry Smith     PetscInt       *jj, k, i, j, n = ltog->n, bs = ltog->bs;
171645b6f7e9SBarry Smith     const PetscInt *ii;
171745b6f7e9SBarry Smith 
17189566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(bs * n, &jj));
171945b6f7e9SBarry Smith     *array = jj;
172045b6f7e9SBarry Smith     k      = 0;
172145b6f7e9SBarry Smith     ii     = ltog->indices;
172245b6f7e9SBarry Smith     for (i = 0; i < n; i++)
17239371c9d4SSatish Balay       for (j = 0; j < bs; j++) jj[k++] = bs * ii[i] + j;
172445b6f7e9SBarry Smith   }
17253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
172686994e45SJed Brown }
172786994e45SJed Brown 
172886994e45SJed Brown /*@C
1729cab54364SBarry Smith    ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with `ISLocalToGlobalMappingGetIndices()`
173086994e45SJed Brown 
173186994e45SJed Brown    Not Collective
173286994e45SJed Brown 
17334165533cSJose E. Roman    Input Parameters:
173486994e45SJed Brown + ltog - local to global mapping
173586994e45SJed Brown - array - array of indices
173686994e45SJed Brown 
173786994e45SJed Brown    Level: advanced
173886994e45SJed Brown 
1739cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingGetIndices()`
174086994e45SJed Brown @*/
1741d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog, const PetscInt **array)
1742d71ae5a4SJacob Faibussowitsch {
174386994e45SJed Brown   PetscFunctionBegin;
174486994e45SJed Brown   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
174586994e45SJed Brown   PetscValidPointer(array, 2);
1746c9cc58a2SBarry Smith   PetscCheck(ltog->bs != 1 || *array == ltog->indices, PETSC_COMM_SELF, PETSC_ERR_ARG_BADPTR, "Trying to return mismatched pointer");
174745b6f7e9SBarry Smith 
174848a46eb9SPierre Jolivet   if (ltog->bs > 1) PetscCall(PetscFree(*(void **)array));
17493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
175045b6f7e9SBarry Smith }
175145b6f7e9SBarry Smith 
175245b6f7e9SBarry Smith /*@C
175345b6f7e9SBarry Smith    ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block
175445b6f7e9SBarry Smith 
175545b6f7e9SBarry Smith    Not Collective
175645b6f7e9SBarry Smith 
17574165533cSJose E. Roman    Input Parameter:
175845b6f7e9SBarry Smith . ltog - local to global mapping
175945b6f7e9SBarry Smith 
17604165533cSJose E. Roman    Output Parameter:
176145b6f7e9SBarry Smith . array - array of indices
176245b6f7e9SBarry Smith 
176345b6f7e9SBarry Smith    Level: advanced
176445b6f7e9SBarry Smith 
1765cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingRestoreBlockIndices()`
176645b6f7e9SBarry Smith @*/
1767d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog, const PetscInt **array)
1768d71ae5a4SJacob Faibussowitsch {
176945b6f7e9SBarry Smith   PetscFunctionBegin;
177045b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
177145b6f7e9SBarry Smith   PetscValidPointer(array, 2);
177245b6f7e9SBarry Smith   *array = ltog->indices;
17733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
177445b6f7e9SBarry Smith }
177545b6f7e9SBarry Smith 
177645b6f7e9SBarry Smith /*@C
1777cab54364SBarry Smith    ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with `ISLocalToGlobalMappingGetBlockIndices()`
177845b6f7e9SBarry Smith 
177945b6f7e9SBarry Smith    Not Collective
178045b6f7e9SBarry Smith 
17814165533cSJose E. Roman    Input Parameters:
178245b6f7e9SBarry Smith + ltog - local to global mapping
178345b6f7e9SBarry Smith - array - array of indices
178445b6f7e9SBarry Smith 
178545b6f7e9SBarry Smith    Level: advanced
178645b6f7e9SBarry Smith 
1787cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingGetIndices()`
178845b6f7e9SBarry Smith @*/
1789d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog, const PetscInt **array)
1790d71ae5a4SJacob Faibussowitsch {
179145b6f7e9SBarry Smith   PetscFunctionBegin;
179245b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
179345b6f7e9SBarry Smith   PetscValidPointer(array, 2);
179408401ef6SPierre Jolivet   PetscCheck(*array == ltog->indices, PETSC_COMM_SELF, PETSC_ERR_ARG_BADPTR, "Trying to return mismatched pointer");
17950298fd71SBarry Smith   *array = NULL;
17963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
179786994e45SJed Brown }
1798f7efa3c7SJed Brown 
1799f7efa3c7SJed Brown /*@C
1800f7efa3c7SJed Brown    ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings
1801f7efa3c7SJed Brown 
1802f7efa3c7SJed Brown    Not Collective
1803f7efa3c7SJed Brown 
18044165533cSJose E. Roman    Input Parameters:
1805f7efa3c7SJed Brown + comm - communicator for the new mapping, must contain the communicator of every mapping to concatenate
1806f7efa3c7SJed Brown . n - number of mappings to concatenate
1807f7efa3c7SJed Brown - ltogs - local to global mappings
1808f7efa3c7SJed Brown 
18094165533cSJose E. Roman    Output Parameter:
1810f7efa3c7SJed Brown . ltogcat - new mapping
1811f7efa3c7SJed Brown 
1812f7efa3c7SJed Brown    Level: advanced
1813f7efa3c7SJed Brown 
1814cab54364SBarry Smith    Note:
1815cab54364SBarry Smith    This currently always returns a mapping with block size of 1
1816cab54364SBarry Smith 
1817cab54364SBarry Smith    Developer Note:
1818cab54364SBarry Smith    If all the input mapping have the same block size we could easily handle that as a special case
1819cab54364SBarry Smith 
1820cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingCreate()`
1821f7efa3c7SJed Brown @*/
1822d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm, PetscInt n, const ISLocalToGlobalMapping ltogs[], ISLocalToGlobalMapping *ltogcat)
1823d71ae5a4SJacob Faibussowitsch {
1824f7efa3c7SJed Brown   PetscInt i, cnt, m, *idx;
1825f7efa3c7SJed Brown 
1826f7efa3c7SJed Brown   PetscFunctionBegin;
182708401ef6SPierre Jolivet   PetscCheck(n >= 0, comm, PETSC_ERR_ARG_OUTOFRANGE, "Must have a non-negative number of mappings, given %" PetscInt_FMT, n);
1828f7efa3c7SJed Brown   if (n > 0) PetscValidPointer(ltogs, 3);
1829f7efa3c7SJed Brown   for (i = 0; i < n; i++) PetscValidHeaderSpecific(ltogs[i], IS_LTOGM_CLASSID, 3);
1830f7efa3c7SJed Brown   PetscValidPointer(ltogcat, 4);
1831f7efa3c7SJed Brown   for (cnt = 0, i = 0; i < n; i++) {
18329566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(ltogs[i], &m));
1833f7efa3c7SJed Brown     cnt += m;
1834f7efa3c7SJed Brown   }
18359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(cnt, &idx));
1836f7efa3c7SJed Brown   for (cnt = 0, i = 0; i < n; i++) {
1837f7efa3c7SJed Brown     const PetscInt *subidx;
18389566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(ltogs[i], &m));
18399566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetIndices(ltogs[i], &subidx));
18409566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(&idx[cnt], subidx, m));
18419566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreIndices(ltogs[i], &subidx));
1842f7efa3c7SJed Brown     cnt += m;
1843f7efa3c7SJed Brown   }
18449566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(comm, 1, cnt, idx, PETSC_OWN_POINTER, ltogcat));
18453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1846f7efa3c7SJed Brown }
184704a59952SBarry Smith 
1848413f72f0SBarry Smith /*MC
1849cab54364SBarry Smith       ISLOCALTOGLOBALMAPPINGBASIC - basic implementation of the `ISLocalToGlobalMapping` object. When `ISGlobalToLocalMappingApply()` is
1850413f72f0SBarry Smith                                     used this is good for only small and moderate size problems.
1851413f72f0SBarry Smith 
185220662ed9SBarry Smith    Options Database Key:
1853a2b725a8SWilliam Gropp .   -islocaltoglobalmapping_type basic - select this method
1854413f72f0SBarry Smith 
1855413f72f0SBarry Smith    Level: beginner
1856413f72f0SBarry Smith 
1857cab54364SBarry Smith    Developer Note:
1858cab54364SBarry Smith    This stores all the mapping information on each MPI rank.
1859cab54364SBarry Smith 
1860cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`, `ISLOCALTOGLOBALMAPPINGHASH`
1861413f72f0SBarry Smith M*/
1862d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Basic(ISLocalToGlobalMapping ltog)
1863d71ae5a4SJacob Faibussowitsch {
1864413f72f0SBarry Smith   PetscFunctionBegin;
1865413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Basic;
1866413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Basic;
1867413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Basic;
1868413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Basic;
18693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1870413f72f0SBarry Smith }
1871413f72f0SBarry Smith 
1872413f72f0SBarry Smith /*MC
1873cab54364SBarry Smith       ISLOCALTOGLOBALMAPPINGHASH - hash implementation of the `ISLocalToGlobalMapping` object. When `ISGlobalToLocalMappingApply()` is
1874ed56e8eaSBarry Smith                                     used this is good for large memory problems.
1875413f72f0SBarry Smith 
187620662ed9SBarry Smith    Options Database Key:
1877a2b725a8SWilliam Gropp .   -islocaltoglobalmapping_type hash - select this method
1878413f72f0SBarry Smith 
1879413f72f0SBarry Smith    Level: beginner
1880413f72f0SBarry Smith 
1881cab54364SBarry Smith    Note:
1882cab54364SBarry Smith     This is selected automatically for large problems if the user does not set the type.
1883cab54364SBarry Smith 
1884cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`, `ISLOCALTOGLOBALMAPPINGBASIC`
1885413f72f0SBarry Smith M*/
1886d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Hash(ISLocalToGlobalMapping ltog)
1887d71ae5a4SJacob Faibussowitsch {
1888413f72f0SBarry Smith   PetscFunctionBegin;
1889413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Hash;
1890413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Hash;
1891413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Hash;
1892413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Hash;
18933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1894413f72f0SBarry Smith }
1895413f72f0SBarry Smith 
1896413f72f0SBarry Smith /*@C
1897cab54364SBarry Smith     ISLocalToGlobalMappingRegister -  Registers a method for applying a global to local mapping with an `ISLocalToGlobalMapping`
1898413f72f0SBarry Smith 
1899413f72f0SBarry Smith    Not Collective
1900413f72f0SBarry Smith 
1901413f72f0SBarry Smith    Input Parameters:
1902413f72f0SBarry Smith +  sname - name of a new method
1903*2fe279fdSBarry Smith -  function - routine to create method context
1904413f72f0SBarry Smith 
1905413f72f0SBarry Smith    Sample usage:
1906413f72f0SBarry Smith .vb
1907413f72f0SBarry Smith    ISLocalToGlobalMappingRegister("my_mapper",MyCreate);
1908413f72f0SBarry Smith .ve
1909413f72f0SBarry Smith 
1910ed56e8eaSBarry Smith    Then, your mapping can be chosen with the procedural interface via
1911413f72f0SBarry Smith $     ISLocalToGlobalMappingSetType(ltog,"my_mapper")
1912413f72f0SBarry Smith    or at runtime via the option
1913ed56e8eaSBarry Smith $     -islocaltoglobalmapping_type my_mapper
1914413f72f0SBarry Smith 
1915413f72f0SBarry Smith    Level: advanced
1916413f72f0SBarry Smith 
1917cab54364SBarry Smith    Note:
1918cab54364SBarry Smith    `ISLocalToGlobalMappingRegister()` may be called multiple times to add several user-defined mappings.
1919413f72f0SBarry Smith 
1920cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingRegisterAll()`, `ISLocalToGlobalMappingRegisterDestroy()`, `ISLOCALTOGLOBALMAPPINGBASIC`,
1921cab54364SBarry Smith           `ISLOCALTOGLOBALMAPPINGHASH`, `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApply()`
1922413f72f0SBarry Smith @*/
1923d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRegister(const char sname[], PetscErrorCode (*function)(ISLocalToGlobalMapping))
1924d71ae5a4SJacob Faibussowitsch {
1925413f72f0SBarry Smith   PetscFunctionBegin;
19269566063dSJacob Faibussowitsch   PetscCall(ISInitializePackage());
19279566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&ISLocalToGlobalMappingList, sname, function));
19283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1929413f72f0SBarry Smith }
1930413f72f0SBarry Smith 
1931413f72f0SBarry Smith /*@C
1932cab54364SBarry Smith    ISLocalToGlobalMappingSetType - Sets the implementation type `ISLocalToGlobalMapping` will use
1933413f72f0SBarry Smith 
1934c3339decSBarry Smith    Logically Collective
1935413f72f0SBarry Smith 
1936413f72f0SBarry Smith    Input Parameters:
1937cab54364SBarry Smith +  ltog - the `ISLocalToGlobalMapping` object
1938413f72f0SBarry Smith -  type - a known method
1939413f72f0SBarry Smith 
1940413f72f0SBarry Smith    Options Database Key:
1941cab54364SBarry Smith .  -islocaltoglobalmapping_type  <method> - Sets the method; use -help for a list of available methods (for instance, basic or hash)
1942cab54364SBarry Smith 
1943cab54364SBarry Smith   Level: intermediate
1944413f72f0SBarry Smith 
1945413f72f0SBarry Smith    Notes:
194620662ed9SBarry Smith    See `ISLocalToGlobalMappingType` for available methods
1947413f72f0SBarry Smith 
1948cab54364SBarry Smith   Normally, it is best to use the `ISLocalToGlobalMappingSetFromOptions()` command and
1949cab54364SBarry Smith   then set the `ISLocalToGlobalMappingType` from the options database rather than by using
1950413f72f0SBarry Smith   this routine.
1951413f72f0SBarry Smith 
1952cab54364SBarry Smith   Developer Note:
1953cab54364SBarry Smith   `ISLocalToGlobalMappingRegister()` is used to add new types to `ISLocalToGlobalMappingList` from which they
1954cab54364SBarry Smith   are accessed by `ISLocalToGlobalMappingSetType()`.
1955413f72f0SBarry Smith 
1956cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingRegister()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingGetType()`
1957413f72f0SBarry Smith @*/
1958d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType type)
1959d71ae5a4SJacob Faibussowitsch {
1960413f72f0SBarry Smith   PetscBool match;
19615f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(ISLocalToGlobalMapping) = NULL;
1962413f72f0SBarry Smith 
1963413f72f0SBarry Smith   PetscFunctionBegin;
1964413f72f0SBarry Smith   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
1965a0d79125SStefano Zampini   if (type) PetscValidCharPointer(type, 2);
1966413f72f0SBarry Smith 
19679566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ltog, type, &match));
19683ba16761SJacob Faibussowitsch   if (match) PetscFunctionReturn(PETSC_SUCCESS);
1969413f72f0SBarry Smith 
1970a0d79125SStefano Zampini   /* L2G maps defer type setup at globaltolocal calls, allow passing NULL here */
1971a0d79125SStefano Zampini   if (type) {
19729566063dSJacob Faibussowitsch     PetscCall(PetscFunctionListFind(ISLocalToGlobalMappingList, type, &r));
1973a0d79125SStefano Zampini     PetscCheck(r, PetscObjectComm((PetscObject)ltog), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested ISLocalToGlobalMapping type %s", type);
1974a0d79125SStefano Zampini   }
1975413f72f0SBarry Smith   /* Destroy the previous private LTOG context */
1976dbbe0bcdSBarry Smith   PetscTryTypeMethod(ltog, destroy);
1977413f72f0SBarry Smith   ltog->ops->destroy = NULL;
1978dbbe0bcdSBarry Smith 
19799566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)ltog, type));
19809566063dSJacob Faibussowitsch   if (r) PetscCall((*r)(ltog));
19813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1982a0d79125SStefano Zampini }
1983a0d79125SStefano Zampini 
1984a0d79125SStefano Zampini /*@C
1985cab54364SBarry Smith    ISLocalToGlobalMappingGetType - Get the type of the `ISLocalToGlobalMapping`
1986a0d79125SStefano Zampini 
1987a0d79125SStefano Zampini    Not Collective
1988a0d79125SStefano Zampini 
1989a0d79125SStefano Zampini    Input Parameter:
1990cab54364SBarry Smith .  ltog - the `ISLocalToGlobalMapping` object
1991a0d79125SStefano Zampini 
1992a0d79125SStefano Zampini    Output Parameter:
1993a0d79125SStefano Zampini .  type - the type
1994a0d79125SStefano Zampini 
199549762cbcSSatish Balay    Level: intermediate
199649762cbcSSatish Balay 
1997cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingRegister()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`
1998a0d79125SStefano Zampini @*/
1999d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType *type)
2000d71ae5a4SJacob Faibussowitsch {
2001a0d79125SStefano Zampini   PetscFunctionBegin;
2002a0d79125SStefano Zampini   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
2003a0d79125SStefano Zampini   PetscValidPointer(type, 2);
2004a0d79125SStefano Zampini   *type = ((PetscObject)ltog)->type_name;
20053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2006413f72f0SBarry Smith }
2007413f72f0SBarry Smith 
2008413f72f0SBarry Smith PetscBool ISLocalToGlobalMappingRegisterAllCalled = PETSC_FALSE;
2009413f72f0SBarry Smith 
2010413f72f0SBarry Smith /*@C
2011cab54364SBarry Smith   ISLocalToGlobalMappingRegisterAll - Registers all of the local to global mapping components in the `IS` package.
2012413f72f0SBarry Smith 
2013413f72f0SBarry Smith   Not Collective
2014413f72f0SBarry Smith 
2015413f72f0SBarry Smith   Level: advanced
2016413f72f0SBarry Smith 
2017cab54364SBarry Smith .seealso: [](sec_scatter), `ISRegister()`, `ISLocalToGlobalRegister()`
2018413f72f0SBarry Smith @*/
2019d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRegisterAll(void)
2020d71ae5a4SJacob Faibussowitsch {
2021413f72f0SBarry Smith   PetscFunctionBegin;
20223ba16761SJacob Faibussowitsch   if (ISLocalToGlobalMappingRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
2023413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAllCalled = PETSC_TRUE;
20249566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGBASIC, ISLocalToGlobalMappingCreate_Basic));
20259566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGHASH, ISLocalToGlobalMappingCreate_Hash));
20263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2027413f72f0SBarry Smith }
2028