xref: /petsc/src/vec/is/utils/isltog.c (revision 7de13914194a88e5c240f442fe647525b445402a)
1af0996ceSBarry Smith #include <petsc/private/isimpl.h> /*I "petscis.h"  I*/
2e8f14785SLisandro Dalcin #include <petsc/private/hashmapi.h>
30c312b8eSJed Brown #include <petscsf.h>
4665c2dedSJed Brown #include <petscviewer.h>
5633354d9SStefano Zampini #include <petscbt.h>
62362add9SBarry Smith 
77087cfbeSBarry Smith PetscClassId          IS_LTOGM_CLASSID;
8633354d9SStefano Zampini static PetscErrorCode ISLocalToGlobalMappingSetUpBlockInfo_Private(ISLocalToGlobalMapping);
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 
221*7de13914SStefano Zampini static PetscErrorCode ISLocalToGlobalMappingResetBlockInfo_Private(ISLocalToGlobalMapping mapping)
222*7de13914SStefano Zampini {
223*7de13914SStefano Zampini   PetscFunctionBegin;
224*7de13914SStefano Zampini   PetscCall(PetscFree(mapping->info_procs));
225*7de13914SStefano Zampini   PetscCall(PetscFree(mapping->info_numprocs));
226*7de13914SStefano Zampini   if (mapping->info_indices) {
227*7de13914SStefano Zampini     for (PetscInt i = 0; i < mapping->info_nproc; i++) PetscCall(PetscFree(mapping->info_indices[i]));
228*7de13914SStefano Zampini     PetscCall(PetscFree(mapping->info_indices));
229*7de13914SStefano Zampini   }
230*7de13914SStefano Zampini   if (mapping->info_nodei) PetscCall(PetscFree(mapping->info_nodei[0]));
231*7de13914SStefano Zampini   PetscCall(PetscFree2(mapping->info_nodec, mapping->info_nodei));
232*7de13914SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
233*7de13914SStefano Zampini }
234*7de13914SStefano Zampini 
235413f72f0SBarry Smith #define GTOLTYPE _Basic
236413f72f0SBarry Smith #define GTOLNAME _Basic
237541bf97eSAdrian Croucher #define GTOLBS   mapping->bs
2389371c9d4SSatish Balay #define GTOL(g, local) \
2399371c9d4SSatish Balay   do { \
240413f72f0SBarry Smith     local = map->globals[g / bs - start]; \
2410040bde1SJunchao Zhang     if (local >= 0) local = bs * local + (g % bs); \
242413f72f0SBarry Smith   } while (0)
243541bf97eSAdrian Croucher 
244413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
245413f72f0SBarry Smith 
246413f72f0SBarry Smith #define GTOLTYPE _Basic
247413f72f0SBarry Smith #define GTOLNAME Block_Basic
248541bf97eSAdrian Croucher #define GTOLBS   1
2499371c9d4SSatish Balay #define GTOL(g, local) \
250d71ae5a4SJacob Faibussowitsch   do { \
251d71ae5a4SJacob Faibussowitsch     local = map->globals[g - start]; \
252d71ae5a4SJacob Faibussowitsch   } while (0)
253413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
254413f72f0SBarry Smith 
255413f72f0SBarry Smith #define GTOLTYPE _Hash
256413f72f0SBarry Smith #define GTOLNAME _Hash
257541bf97eSAdrian Croucher #define GTOLBS   mapping->bs
2589371c9d4SSatish Balay #define GTOL(g, local) \
2599371c9d4SSatish Balay   do { \
260e8f14785SLisandro Dalcin     (void)PetscHMapIGet(map->globalht, g / bs, &local); \
2610040bde1SJunchao Zhang     if (local >= 0) local = bs * local + (g % bs); \
262413f72f0SBarry Smith   } while (0)
263413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
264413f72f0SBarry Smith 
265413f72f0SBarry Smith #define GTOLTYPE _Hash
266413f72f0SBarry Smith #define GTOLNAME Block_Hash
267541bf97eSAdrian Croucher #define GTOLBS   1
2689371c9d4SSatish Balay #define GTOL(g, local) \
269d71ae5a4SJacob Faibussowitsch   do { \
270d71ae5a4SJacob Faibussowitsch     (void)PetscHMapIGet(map->globalht, g, &local); \
271d71ae5a4SJacob Faibussowitsch   } while (0)
272413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
273413f72f0SBarry Smith 
2746658fb44Sstefano_zampini /*@
2756658fb44Sstefano_zampini   ISLocalToGlobalMappingDuplicate - Duplicates the local to global mapping object
2766658fb44Sstefano_zampini 
2776658fb44Sstefano_zampini   Not Collective
2786658fb44Sstefano_zampini 
2796658fb44Sstefano_zampini   Input Parameter:
2806658fb44Sstefano_zampini . ltog - local to global mapping
2816658fb44Sstefano_zampini 
2826658fb44Sstefano_zampini   Output Parameter:
2836658fb44Sstefano_zampini . nltog - the duplicated local to global mapping
2846658fb44Sstefano_zampini 
2856658fb44Sstefano_zampini   Level: advanced
2866658fb44Sstefano_zampini 
287cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`
2886658fb44Sstefano_zampini @*/
289d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingDuplicate(ISLocalToGlobalMapping ltog, ISLocalToGlobalMapping *nltog)
290d71ae5a4SJacob Faibussowitsch {
291a0d79125SStefano Zampini   ISLocalToGlobalMappingType l2gtype;
2926658fb44Sstefano_zampini 
2936658fb44Sstefano_zampini   PetscFunctionBegin;
2946658fb44Sstefano_zampini   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
2959566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)ltog), ltog->bs, ltog->n, ltog->indices, PETSC_COPY_VALUES, nltog));
2969566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetType(ltog, &l2gtype));
2979566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingSetType(*nltog, l2gtype));
2983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2996658fb44Sstefano_zampini }
3006658fb44Sstefano_zampini 
301565245c5SBarry Smith /*@
302107e9a97SBarry Smith   ISLocalToGlobalMappingGetSize - Gets the local size of a local to global mapping
3033b9aefa3SBarry Smith 
3043b9aefa3SBarry Smith   Not Collective
3053b9aefa3SBarry Smith 
3063b9aefa3SBarry Smith   Input Parameter:
30738b5cf2dSJacob Faibussowitsch . mapping - local to global mapping
3083b9aefa3SBarry Smith 
3093b9aefa3SBarry Smith   Output Parameter:
310cab54364SBarry Smith . n - the number of entries in the local mapping, `ISLocalToGlobalMappingGetIndices()` returns an array of this length
3113b9aefa3SBarry Smith 
3123b9aefa3SBarry Smith   Level: advanced
3133b9aefa3SBarry Smith 
314cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`
3153b9aefa3SBarry Smith @*/
316d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping mapping, PetscInt *n)
317d71ae5a4SJacob Faibussowitsch {
3183b9aefa3SBarry Smith   PetscFunctionBegin;
3190700a824SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
3204f572ea9SToby Isaac   PetscAssertPointer(n, 2);
321107e9a97SBarry Smith   *n = mapping->bs * mapping->n;
3223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3233b9aefa3SBarry Smith }
3243b9aefa3SBarry Smith 
3255a5d4f66SBarry Smith /*@C
326cab54364SBarry Smith   ISLocalToGlobalMappingViewFromOptions - View an `ISLocalToGlobalMapping` based on values in the options database
327fe2efc57SMark 
328c3339decSBarry Smith   Collective
329fe2efc57SMark 
330fe2efc57SMark   Input Parameters:
331fe2efc57SMark + A    - the local to global mapping object
33220662ed9SBarry Smith . obj  - Optional object that provides the options prefix used for the options database query
333736c3998SJose E. Roman - name - command line option
334fe2efc57SMark 
335fe2efc57SMark   Level: intermediate
336cab54364SBarry Smith 
33720662ed9SBarry Smith   Note:
33820662ed9SBarry Smith   See `PetscObjectViewFromOptions()` for the available `PetscViewer` and `PetscViewerFormat`
33920662ed9SBarry Smith 
34020662ed9SBarry Smith .seealso: [](sec_scatter), `PetscViewer`, ``ISLocalToGlobalMapping`, `ISLocalToGlobalMappingView`, `PetscObjectViewFromOptions()`, `ISLocalToGlobalMappingCreate()`
341fe2efc57SMark @*/
342d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingViewFromOptions(ISLocalToGlobalMapping A, PetscObject obj, const char name[])
343d71ae5a4SJacob Faibussowitsch {
344fe2efc57SMark   PetscFunctionBegin;
345fe2efc57SMark   PetscValidHeaderSpecific(A, IS_LTOGM_CLASSID, 1);
3469566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
3473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
348fe2efc57SMark }
349fe2efc57SMark 
350fe2efc57SMark /*@C
3515a5d4f66SBarry Smith   ISLocalToGlobalMappingView - View a local to global mapping
3525a5d4f66SBarry Smith 
353*7de13914SStefano Zampini   Collective on viewer
354b9cd556bSLois Curfman McInnes 
3555a5d4f66SBarry Smith   Input Parameters:
35638b5cf2dSJacob Faibussowitsch + mapping - local to global mapping
3573b9aefa3SBarry Smith - viewer  - viewer
3585a5d4f66SBarry Smith 
359*7de13914SStefano Zampini   Level: intermediate
360a997ad1aSLois Curfman McInnes 
36120662ed9SBarry Smith .seealso: [](sec_scatter), `PetscViewer`, `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`
3625a5d4f66SBarry Smith @*/
363d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingView(ISLocalToGlobalMapping mapping, PetscViewer viewer)
364d71ae5a4SJacob Faibussowitsch {
365*7de13914SStefano Zampini   PetscBool         iascii, isbinary;
366*7de13914SStefano Zampini   PetscViewerFormat format;
3675a5d4f66SBarry Smith 
3685a5d4f66SBarry Smith   PetscFunctionBegin;
3690700a824SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
37048a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mapping), &viewer));
3710700a824SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
3725a5d4f66SBarry Smith 
3739566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
374*7de13914SStefano Zampini   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
375*7de13914SStefano Zampini   PetscCall(PetscViewerGetFormat(viewer, &format));
37632077d6dSBarry Smith   if (iascii) {
377*7de13914SStefano Zampini     if (format == PETSC_VIEWER_ASCII_MATLAB) {
378*7de13914SStefano Zampini       const PetscInt *idxs;
379*7de13914SStefano Zampini       IS              is;
380*7de13914SStefano Zampini       const char     *name = ((PetscObject)mapping)->name;
381*7de13914SStefano Zampini       char            iname[PETSC_MAX_PATH_LEN];
382*7de13914SStefano Zampini 
383*7de13914SStefano Zampini       PetscCall(PetscSNPrintf(iname, sizeof(iname), "%sl2g", name ? name : ""));
384*7de13914SStefano Zampini       PetscCall(ISLocalToGlobalMappingGetIndices(mapping, &idxs));
385*7de13914SStefano Zampini       PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)viewer), mapping->n * mapping->bs, idxs, PETSC_USE_POINTER, &is));
386*7de13914SStefano Zampini       PetscCall(PetscObjectSetName((PetscObject)is, iname));
387*7de13914SStefano Zampini       PetscCall(ISView(is, viewer));
388*7de13914SStefano Zampini       PetscCall(ISLocalToGlobalMappingRestoreIndices(mapping, &idxs));
389*7de13914SStefano Zampini       PetscCall(ISDestroy(&is));
390*7de13914SStefano Zampini     } else {
391*7de13914SStefano Zampini       PetscMPIInt rank;
392*7de13914SStefano Zampini 
393*7de13914SStefano Zampini       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mapping), &rank));
3949566063dSJacob Faibussowitsch       PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)mapping, viewer));
3959566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
396*7de13914SStefano Zampini       for (PetscInt i = 0; i < mapping->n; i++) {
397f2c6b1a2SJed Brown         PetscInt bs = mapping->bs, g = mapping->indices[i];
398f2c6b1a2SJed Brown         if (bs == 1) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %" PetscInt_FMT " %" PetscInt_FMT "\n", rank, i, g));
399f2c6b1a2SJed 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));
400f2c6b1a2SJed Brown       }
4019566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
4029566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
4031575c14dSBarry Smith     }
404*7de13914SStefano Zampini   } else if (isbinary) {
405*7de13914SStefano Zampini     PetscBool skipHeader;
406*7de13914SStefano Zampini 
407*7de13914SStefano Zampini     PetscCall(PetscViewerSetUp(viewer));
408*7de13914SStefano Zampini     PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &skipHeader));
409*7de13914SStefano Zampini     if (!skipHeader) {
410*7de13914SStefano Zampini       PetscMPIInt size;
411*7de13914SStefano Zampini       PetscInt    tr[3];
412*7de13914SStefano Zampini 
413*7de13914SStefano Zampini       PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)viewer), &size));
414*7de13914SStefano Zampini       tr[0] = IS_LTOGM_FILE_CLASSID;
415*7de13914SStefano Zampini       tr[1] = mapping->bs;
416*7de13914SStefano Zampini       tr[2] = size;
417*7de13914SStefano Zampini       PetscCall(PetscViewerBinaryWrite(viewer, tr, 3, PETSC_INT));
418*7de13914SStefano Zampini       PetscCall(PetscViewerBinaryWriteAll(viewer, &mapping->n, 1, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_INT));
419*7de13914SStefano Zampini     }
420*7de13914SStefano Zampini     /* write block indices */
421*7de13914SStefano Zampini     PetscCall(PetscViewerBinaryWriteAll(viewer, mapping->indices, mapping->n, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_INT));
422*7de13914SStefano Zampini   }
423*7de13914SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
424*7de13914SStefano Zampini }
425*7de13914SStefano Zampini 
426*7de13914SStefano Zampini /*@
427*7de13914SStefano Zampini   ISLocalToGlobalMappingLoad - Loads a local-to-global mapping that has been stored in binary format.
428*7de13914SStefano Zampini 
429*7de13914SStefano Zampini   Collective on viewer
430*7de13914SStefano Zampini 
431*7de13914SStefano Zampini   Input Parameters:
432*7de13914SStefano Zampini + mapping - the newly loaded map, this needs to have been created with `ISLocalToGlobalMappingCreate()` or some related function before a call to `ISLocalToGlobalMappingLoad()`
433*7de13914SStefano Zampini - viewer  - binary file viewer, obtained from `PetscViewerBinaryOpen()`
434*7de13914SStefano Zampini 
435*7de13914SStefano Zampini   Level: intermediate
436*7de13914SStefano Zampini 
437*7de13914SStefano Zampini .seealso: [](sec_scatter), `PetscViewer`, `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingView()`, `ISLocalToGlobalMappingCreate()`
438*7de13914SStefano Zampini @*/
439*7de13914SStefano Zampini PetscErrorCode ISLocalToGlobalMappingLoad(ISLocalToGlobalMapping mapping, PetscViewer viewer)
440*7de13914SStefano Zampini {
441*7de13914SStefano Zampini   PetscBool isbinary, skipHeader;
442*7de13914SStefano Zampini 
443*7de13914SStefano Zampini   PetscFunctionBegin;
444*7de13914SStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
445*7de13914SStefano Zampini   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
446*7de13914SStefano Zampini   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
447*7de13914SStefano Zampini   PetscCheck(isbinary, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Invalid viewer of type %s", ((PetscObject)viewer)->type_name);
448*7de13914SStefano Zampini 
449*7de13914SStefano Zampini   /* reset previous data */
450*7de13914SStefano Zampini   PetscCall(ISLocalToGlobalMappingResetBlockInfo_Private(mapping));
451*7de13914SStefano Zampini 
452*7de13914SStefano Zampini   PetscCall(PetscViewerSetUp(viewer));
453*7de13914SStefano Zampini   PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &skipHeader));
454*7de13914SStefano Zampini 
455*7de13914SStefano Zampini   /* When skipping header, it assumes bs and n have been already set */
456*7de13914SStefano Zampini   if (!skipHeader) {
457*7de13914SStefano Zampini     MPI_Comm comm = PetscObjectComm((PetscObject)viewer);
458*7de13914SStefano Zampini     PetscInt tr[3], nold = mapping->n, *sizes, nmaps = PETSC_DECIDE, st = 0;
459*7de13914SStefano Zampini 
460*7de13914SStefano Zampini     PetscCall(PetscViewerBinaryRead(viewer, tr, 3, NULL, PETSC_INT));
461*7de13914SStefano Zampini     PetscCheck(tr[0] == IS_LTOGM_FILE_CLASSID, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Not a local-to-global map next in file");
462*7de13914SStefano Zampini 
463*7de13914SStefano Zampini     mapping->bs = tr[1];
464*7de13914SStefano Zampini     PetscCall(PetscMalloc1(tr[2], &sizes));
465*7de13914SStefano Zampini     PetscCall(PetscViewerBinaryRead(viewer, sizes, tr[2], NULL, PETSC_INT));
466*7de13914SStefano Zampini 
467*7de13914SStefano Zampini     /* consume the input, read multiple maps per process if needed */
468*7de13914SStefano Zampini     PetscCall(PetscSplitOwnership(comm, &nmaps, &tr[2]));
469*7de13914SStefano Zampini     PetscCallMPI(MPI_Exscan(&nmaps, &st, 1, MPIU_INT, MPI_SUM, comm));
470*7de13914SStefano Zampini     mapping->n = 0;
471*7de13914SStefano Zampini     for (PetscInt i = st; i < st + nmaps; i++) mapping->n += sizes[i];
472*7de13914SStefano Zampini     PetscCall(PetscFree(sizes));
473*7de13914SStefano Zampini 
474*7de13914SStefano Zampini     if (nold != mapping->n) {
475*7de13914SStefano Zampini       if (mapping->dealloc_indices) PetscCall(PetscFree(mapping->indices));
476*7de13914SStefano Zampini       mapping->indices = NULL;
477*7de13914SStefano Zampini     }
478*7de13914SStefano Zampini   }
479*7de13914SStefano Zampini 
480*7de13914SStefano Zampini   /* read indices */
481*7de13914SStefano Zampini   if (mapping->n && !mapping->indices) {
482*7de13914SStefano Zampini     PetscCall(PetscMalloc1(mapping->n, &mapping->indices));
483*7de13914SStefano Zampini     mapping->dealloc_indices = PETSC_TRUE;
484*7de13914SStefano Zampini   }
485*7de13914SStefano Zampini   PetscCall(PetscViewerBinaryReadAll(viewer, mapping->indices, mapping->n, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_INT));
4863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4875a5d4f66SBarry Smith }
4885a5d4f66SBarry Smith 
4891f428162SBarry Smith /*@
4902bdab257SBarry Smith   ISLocalToGlobalMappingCreateIS - Creates a mapping between a local (0 to n)
4912bdab257SBarry Smith   ordering and a global parallel ordering.
4922bdab257SBarry Smith 
49320662ed9SBarry Smith   Not Collective
494b9cd556bSLois Curfman McInnes 
495a997ad1aSLois Curfman McInnes   Input Parameter:
4968c03b21aSDmitry Karpeev . is - index set containing the global numbers for each local number
4972bdab257SBarry Smith 
498a997ad1aSLois Curfman McInnes   Output Parameter:
4992bdab257SBarry Smith . mapping - new mapping data structure
5002bdab257SBarry Smith 
501a997ad1aSLois Curfman McInnes   Level: advanced
502a997ad1aSLois Curfman McInnes 
503cab54364SBarry Smith   Note:
504cab54364SBarry Smith   the block size of the `IS` determines the block size of the mapping
505cab54364SBarry Smith 
506cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetFromOptions()`
5072bdab257SBarry Smith @*/
508d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingCreateIS(IS is, ISLocalToGlobalMapping *mapping)
509d71ae5a4SJacob Faibussowitsch {
5103bbf0e92SBarry Smith   PetscInt        n, bs;
5115d0c19d7SBarry Smith   const PetscInt *indices;
5122bdab257SBarry Smith   MPI_Comm        comm;
5133bbf0e92SBarry Smith   PetscBool       isblock;
5143a40ed3dSBarry Smith 
5153a40ed3dSBarry Smith   PetscFunctionBegin;
5160700a824SBarry Smith   PetscValidHeaderSpecific(is, IS_CLASSID, 1);
5174f572ea9SToby Isaac   PetscAssertPointer(mapping, 2);
5182bdab257SBarry Smith 
5199566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)is, &comm));
5209566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is, &n));
5219566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)is, ISBLOCK, &isblock));
5226006e8d2SBarry Smith   if (!isblock) {
5239566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(is, &indices));
5249566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(comm, 1, n, indices, PETSC_COPY_VALUES, mapping));
5259566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(is, &indices));
5266006e8d2SBarry Smith   } else {
5279566063dSJacob Faibussowitsch     PetscCall(ISGetBlockSize(is, &bs));
5289566063dSJacob Faibussowitsch     PetscCall(ISBlockGetIndices(is, &indices));
5299566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(comm, bs, n / bs, indices, PETSC_COPY_VALUES, mapping));
5309566063dSJacob Faibussowitsch     PetscCall(ISBlockRestoreIndices(is, &indices));
5316006e8d2SBarry Smith   }
5323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5332bdab257SBarry Smith }
5345a5d4f66SBarry Smith 
535a4d96a55SJed Brown /*@C
536a4d96a55SJed Brown   ISLocalToGlobalMappingCreateSF - Creates a mapping between a local (0 to n)
537a4d96a55SJed Brown   ordering and a global parallel ordering.
538a4d96a55SJed Brown 
539a4d96a55SJed Brown   Collective
540a4d96a55SJed Brown 
541d8d19677SJose E. Roman   Input Parameters:
542a4d96a55SJed Brown + sf    - star forest mapping contiguous local indices to (rank, offset)
543cab54364SBarry Smith - start - first global index on this process, or `PETSC_DECIDE` to compute contiguous global numbering automatically
544a4d96a55SJed Brown 
545a4d96a55SJed Brown   Output Parameter:
546a4d96a55SJed Brown . mapping - new mapping data structure
547a4d96a55SJed Brown 
548a4d96a55SJed Brown   Level: advanced
549a4d96a55SJed Brown 
55020662ed9SBarry Smith   Note:
551633354d9SStefano Zampini   If a process calls this function with `start` = `PETSC_DECIDE` then all processes must, otherwise the program will hang.
5529a535bafSVaclav Hapla 
553cab54364SBarry Smith .seealso: [](sec_scatter), `PetscSF`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingSetFromOptions()`
554a4d96a55SJed Brown @*/
555d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF sf, PetscInt start, ISLocalToGlobalMapping *mapping)
556d71ae5a4SJacob Faibussowitsch {
557a4d96a55SJed Brown   PetscInt i, maxlocal, nroots, nleaves, *globals, *ltog;
558a4d96a55SJed Brown   MPI_Comm comm;
559a4d96a55SJed Brown 
560a4d96a55SJed Brown   PetscFunctionBegin;
561a4d96a55SJed Brown   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
5624f572ea9SToby Isaac   PetscAssertPointer(mapping, 3);
5639566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
56441f4c31fSVaclav Hapla   PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL));
5659a535bafSVaclav Hapla   if (start == PETSC_DECIDE) {
5669a535bafSVaclav Hapla     start = 0;
5679566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Exscan(&nroots, &start, 1, MPIU_INT, MPI_SUM, comm));
56841f4c31fSVaclav Hapla   } else PetscCheck(start >= 0, comm, PETSC_ERR_ARG_OUTOFRANGE, "start must be nonnegative or PETSC_DECIDE");
56941f4c31fSVaclav Hapla   PetscCall(PetscSFGetLeafRange(sf, NULL, &maxlocal));
57041f4c31fSVaclav Hapla   ++maxlocal;
5719566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nroots, &globals));
5729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxlocal, &ltog));
573a4d96a55SJed Brown   for (i = 0; i < nroots; i++) globals[i] = start + i;
574a4d96a55SJed Brown   for (i = 0; i < maxlocal; i++) ltog[i] = -1;
5759566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sf, MPIU_INT, globals, ltog, MPI_REPLACE));
5769566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sf, MPIU_INT, globals, ltog, MPI_REPLACE));
5779566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(comm, 1, maxlocal, ltog, PETSC_OWN_POINTER, mapping));
5789566063dSJacob Faibussowitsch   PetscCall(PetscFree(globals));
5793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
580a4d96a55SJed Brown }
581b46b645bSBarry Smith 
58263fa5c83Sstefano_zampini /*@
58363fa5c83Sstefano_zampini   ISLocalToGlobalMappingSetBlockSize - Sets the blocksize of the mapping
58463fa5c83Sstefano_zampini 
58520662ed9SBarry Smith   Not Collective
58663fa5c83Sstefano_zampini 
58763fa5c83Sstefano_zampini   Input Parameters:
588a2b725a8SWilliam Gropp + mapping - mapping data structure
589a2b725a8SWilliam Gropp - bs      - the blocksize
59063fa5c83Sstefano_zampini 
59163fa5c83Sstefano_zampini   Level: advanced
59263fa5c83Sstefano_zampini 
593cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`
59463fa5c83Sstefano_zampini @*/
595d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingSetBlockSize(ISLocalToGlobalMapping mapping, PetscInt bs)
596d71ae5a4SJacob Faibussowitsch {
597a59f3c4dSstefano_zampini   PetscInt       *nid;
598a59f3c4dSstefano_zampini   const PetscInt *oid;
599a59f3c4dSstefano_zampini   PetscInt        i, cn, on, obs, nn;
60063fa5c83Sstefano_zampini 
60163fa5c83Sstefano_zampini   PetscFunctionBegin;
60263fa5c83Sstefano_zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
60308401ef6SPierre Jolivet   PetscCheck(bs >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid block size %" PetscInt_FMT, bs);
6043ba16761SJacob Faibussowitsch   if (bs == mapping->bs) PetscFunctionReturn(PETSC_SUCCESS);
60563fa5c83Sstefano_zampini   on  = mapping->n;
60663fa5c83Sstefano_zampini   obs = mapping->bs;
60763fa5c83Sstefano_zampini   oid = mapping->indices;
60863fa5c83Sstefano_zampini   nn  = (on * obs) / bs;
60908401ef6SPierre 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);
610a59f3c4dSstefano_zampini 
6119566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nn, &nid));
6129566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetIndices(mapping, &oid));
613a59f3c4dSstefano_zampini   for (i = 0; i < nn; i++) {
614a59f3c4dSstefano_zampini     PetscInt j;
615a59f3c4dSstefano_zampini     for (j = 0, cn = 0; j < bs - 1; j++) {
6169371c9d4SSatish Balay       if (oid[i * bs + j] < 0) {
6179371c9d4SSatish Balay         cn++;
6189371c9d4SSatish Balay         continue;
6199371c9d4SSatish Balay       }
62008401ef6SPierre 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]);
621a59f3c4dSstefano_zampini     }
622a59f3c4dSstefano_zampini     if (oid[i * bs + j] < 0) cn++;
6238b7cb0e6Sstefano_zampini     if (cn) {
62408401ef6SPierre 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);
625a59f3c4dSstefano_zampini       nid[i] = -1;
6268b7cb0e6Sstefano_zampini     } else {
627a59f3c4dSstefano_zampini       nid[i] = oid[i * bs] / bs;
62863fa5c83Sstefano_zampini     }
62963fa5c83Sstefano_zampini   }
6309566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreIndices(mapping, &oid));
631a59f3c4dSstefano_zampini 
63263fa5c83Sstefano_zampini   mapping->n  = nn;
63363fa5c83Sstefano_zampini   mapping->bs = bs;
6349566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping->indices));
63563fa5c83Sstefano_zampini   mapping->indices     = nid;
636c9345713Sstefano_zampini   mapping->globalstart = 0;
637c9345713Sstefano_zampini   mapping->globalend   = 0;
6381bd0b88eSStefano Zampini 
6391bd0b88eSStefano Zampini   /* reset the cached information */
640633354d9SStefano Zampini   PetscCall(ISLocalToGlobalMappingResetBlockInfo_Private(mapping));
641dbbe0bcdSBarry Smith   PetscTryTypeMethod(mapping, destroy);
6423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
64363fa5c83Sstefano_zampini }
64463fa5c83Sstefano_zampini 
64545b6f7e9SBarry Smith /*@
64645b6f7e9SBarry Smith   ISLocalToGlobalMappingGetBlockSize - Gets the blocksize of the mapping
64745b6f7e9SBarry Smith   ordering and a global parallel ordering.
64845b6f7e9SBarry Smith 
64945b6f7e9SBarry Smith   Not Collective
65045b6f7e9SBarry Smith 
6512fe279fdSBarry Smith   Input Parameter:
65245b6f7e9SBarry Smith . mapping - mapping data structure
65345b6f7e9SBarry Smith 
65445b6f7e9SBarry Smith   Output Parameter:
65545b6f7e9SBarry Smith . bs - the blocksize
65645b6f7e9SBarry Smith 
65745b6f7e9SBarry Smith   Level: advanced
65845b6f7e9SBarry Smith 
659cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`
66045b6f7e9SBarry Smith @*/
661d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping mapping, PetscInt *bs)
662d71ae5a4SJacob Faibussowitsch {
66345b6f7e9SBarry Smith   PetscFunctionBegin;
664cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
66545b6f7e9SBarry Smith   *bs = mapping->bs;
6663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
66745b6f7e9SBarry Smith }
66845b6f7e9SBarry Smith 
669ba5bb76aSSatish Balay /*@
67090f02eecSBarry Smith   ISLocalToGlobalMappingCreate - Creates a mapping between a local (0 to n)
67190f02eecSBarry Smith   ordering and a global parallel ordering.
6722362add9SBarry Smith 
67389d82c54SBarry Smith   Not Collective, but communicator may have more than one process
674b9cd556bSLois Curfman McInnes 
6752362add9SBarry Smith   Input Parameters:
67689d82c54SBarry Smith + comm    - MPI communicator
677f0413b6fSBarry Smith . bs      - the block size
67828bc9809SBarry Smith . n       - the number of local elements divided by the block size, or equivalently the number of block indices
67928bc9809SBarry 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
680d5ad8652SBarry Smith - mode    - see PetscCopyMode
6812362add9SBarry Smith 
682a997ad1aSLois Curfman McInnes   Output Parameter:
68390f02eecSBarry Smith . mapping - new mapping data structure
6842362add9SBarry Smith 
685cab54364SBarry Smith   Level: advanced
686cab54364SBarry Smith 
68795452b02SPatrick Sanan   Notes:
68895452b02SPatrick Sanan   There is one integer value in indices per block and it represents the actual indices bs*idx + j, where j=0,..,bs-1
689413f72f0SBarry Smith 
690cab54364SBarry Smith   For "small" problems when using `ISGlobalToLocalMappingApply()` and `ISGlobalToLocalMappingApplyBlock()`, the `ISLocalToGlobalMappingType`
691cab54364SBarry Smith   of `ISLOCALTOGLOBALMAPPINGBASIC` will be used; this uses more memory but is faster; this approach is not scalable for extremely large mappings.
692413f72f0SBarry Smith 
693cab54364SBarry Smith   For large problems `ISLOCALTOGLOBALMAPPINGHASH` is used, this is scalable.
69420662ed9SBarry Smith   Use `ISLocalToGlobalMappingSetType()` or call `ISLocalToGlobalMappingSetFromOptions()` with the option
69520662ed9SBarry Smith   `-islocaltoglobalmapping_type` <`basic`,`hash`> to control which is used.
696a997ad1aSLois Curfman McInnes 
69720662ed9SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingSetFromOptions()`,
69820662ed9SBarry Smith           `ISLOCALTOGLOBALMAPPINGBASIC`, `ISLOCALTOGLOBALMAPPINGHASH`
699db781477SPatrick Sanan           `ISLocalToGlobalMappingSetType()`, `ISLocalToGlobalMappingType`
7002362add9SBarry Smith @*/
701d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingCreate(MPI_Comm comm, PetscInt bs, PetscInt n, const PetscInt indices[], PetscCopyMode mode, ISLocalToGlobalMapping *mapping)
702d71ae5a4SJacob Faibussowitsch {
70332dcc486SBarry Smith   PetscInt *in;
704b46b645bSBarry Smith 
705b46b645bSBarry Smith   PetscFunctionBegin;
7064f572ea9SToby Isaac   if (n) PetscAssertPointer(indices, 4);
7074f572ea9SToby Isaac   PetscAssertPointer(mapping, 6);
708b46b645bSBarry Smith 
7090298fd71SBarry Smith   *mapping = NULL;
7109566063dSJacob Faibussowitsch   PetscCall(ISInitializePackage());
7112362add9SBarry Smith 
7129566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(*mapping, IS_LTOGM_CLASSID, "ISLocalToGlobalMapping", "Local to global mapping", "IS", comm, ISLocalToGlobalMappingDestroy, ISLocalToGlobalMappingView));
713d4bb536fSBarry Smith   (*mapping)->n  = n;
714f0413b6fSBarry Smith   (*mapping)->bs = bs;
715d5ad8652SBarry Smith   if (mode == PETSC_COPY_VALUES) {
7169566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &in));
7179566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(in, indices, n));
718d5ad8652SBarry Smith     (*mapping)->indices         = in;
71971910c26SVaclav Hapla     (*mapping)->dealloc_indices = PETSC_TRUE;
7206389a1a1SBarry Smith   } else if (mode == PETSC_OWN_POINTER) {
7216389a1a1SBarry Smith     (*mapping)->indices         = (PetscInt *)indices;
72271910c26SVaclav Hapla     (*mapping)->dealloc_indices = PETSC_TRUE;
72371910c26SVaclav Hapla   } else if (mode == PETSC_USE_POINTER) {
72471910c26SVaclav Hapla     (*mapping)->indices = (PetscInt *)indices;
7259371c9d4SSatish Balay   } else SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid mode %d", mode);
7263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7272362add9SBarry Smith }
7282362add9SBarry Smith 
729413f72f0SBarry Smith PetscFunctionList ISLocalToGlobalMappingList = NULL;
730413f72f0SBarry Smith 
73190f02eecSBarry Smith /*@
7327e99dc12SLawrence Mitchell   ISLocalToGlobalMappingSetFromOptions - Set mapping options from the options database.
7337e99dc12SLawrence Mitchell 
73420662ed9SBarry Smith   Not Collective
7357e99dc12SLawrence Mitchell 
7362fe279fdSBarry Smith   Input Parameter:
7377e99dc12SLawrence Mitchell . mapping - mapping data structure
7387e99dc12SLawrence Mitchell 
73920662ed9SBarry Smith   Options Database Key:
74020662ed9SBarry Smith . -islocaltoglobalmapping_type - <basic,hash> nonscalable and scalable versions
74120662ed9SBarry Smith 
7427e99dc12SLawrence Mitchell   Level: advanced
7437e99dc12SLawrence Mitchell 
74442747ad1SJacob Faibussowitsch .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`,
74542747ad1SJacob Faibussowitsch `ISLocalToGlobalMappingCreateIS()`, `ISLOCALTOGLOBALMAPPINGBASIC`,
74642747ad1SJacob Faibussowitsch `ISLOCALTOGLOBALMAPPINGHASH`, `ISLocalToGlobalMappingSetType()`, `ISLocalToGlobalMappingType`
7477e99dc12SLawrence Mitchell @*/
748d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingSetFromOptions(ISLocalToGlobalMapping mapping)
749d71ae5a4SJacob Faibussowitsch {
750413f72f0SBarry Smith   char                       type[256];
751413f72f0SBarry Smith   ISLocalToGlobalMappingType defaulttype = "Not set";
7527e99dc12SLawrence Mitchell   PetscBool                  flg;
7537e99dc12SLawrence Mitchell 
7547e99dc12SLawrence Mitchell   PetscFunctionBegin;
7557e99dc12SLawrence Mitchell   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
7569566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRegisterAll());
757d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)mapping);
7589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-islocaltoglobalmapping_type", "ISLocalToGlobalMapping method", "ISLocalToGlobalMappingSetType", ISLocalToGlobalMappingList, (char *)(((PetscObject)mapping)->type_name) ? ((PetscObject)mapping)->type_name : defaulttype, type, 256, &flg));
7591baa6e33SBarry Smith   if (flg) PetscCall(ISLocalToGlobalMappingSetType(mapping, type));
760d0609cedSBarry Smith   PetscOptionsEnd();
7613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7627e99dc12SLawrence Mitchell }
7637e99dc12SLawrence Mitchell 
7647e99dc12SLawrence Mitchell /*@
76590f02eecSBarry Smith   ISLocalToGlobalMappingDestroy - Destroys a mapping between a local (0 to n)
76690f02eecSBarry Smith   ordering and a global parallel ordering.
76790f02eecSBarry Smith 
76820662ed9SBarry Smith   Not Collective
769b9cd556bSLois Curfman McInnes 
7702fe279fdSBarry Smith   Input Parameter:
77190f02eecSBarry Smith . mapping - mapping data structure
77290f02eecSBarry Smith 
773a997ad1aSLois Curfman McInnes   Level: advanced
774a997ad1aSLois Curfman McInnes 
775cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingCreate()`
77690f02eecSBarry Smith @*/
777d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping *mapping)
778d71ae5a4SJacob Faibussowitsch {
7793a40ed3dSBarry Smith   PetscFunctionBegin;
7803ba16761SJacob Faibussowitsch   if (!*mapping) PetscFunctionReturn(PETSC_SUCCESS);
781f4f49eeaSPierre Jolivet   PetscValidHeaderSpecific(*mapping, IS_LTOGM_CLASSID, 1);
782f4f49eeaSPierre Jolivet   if (--((PetscObject)*mapping)->refct > 0) {
7839371c9d4SSatish Balay     *mapping = NULL;
7843ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
78571910c26SVaclav Hapla   }
78648a46eb9SPierre Jolivet   if ((*mapping)->dealloc_indices) PetscCall(PetscFree((*mapping)->indices));
787633354d9SStefano Zampini   PetscCall(ISLocalToGlobalMappingResetBlockInfo_Private(*mapping));
788213acdd3SPierre Jolivet   PetscTryTypeMethod(*mapping, destroy);
7899566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(mapping));
7903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
79190f02eecSBarry Smith }
79290f02eecSBarry Smith 
79390f02eecSBarry Smith /*@
794cab54364SBarry Smith   ISLocalToGlobalMappingApplyIS - Creates from an `IS` in the local numbering
795cab54364SBarry Smith   a new index set using the global numbering defined in an `ISLocalToGlobalMapping`
7963acfe500SLois Curfman McInnes   context.
79790f02eecSBarry Smith 
798c3339decSBarry Smith   Collective
799b9cd556bSLois Curfman McInnes 
80090f02eecSBarry Smith   Input Parameters:
801b9cd556bSLois Curfman McInnes + mapping - mapping between local and global numbering
802b9cd556bSLois Curfman McInnes - is      - index set in local numbering
80390f02eecSBarry Smith 
804cab54364SBarry Smith   Output Parameter:
80590f02eecSBarry Smith . newis - index set in global numbering
80690f02eecSBarry Smith 
807a997ad1aSLois Curfman McInnes   Level: advanced
808a997ad1aSLois Curfman McInnes 
809cab54364SBarry Smith   Note:
81020662ed9SBarry Smith   The output `IS` will have the same communicator as the input `IS`.
811cab54364SBarry Smith 
812cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingCreate()`,
813db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`, `ISGlobalToLocalMappingApply()`
81490f02eecSBarry Smith @*/
815d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping, IS is, IS *newis)
816d71ae5a4SJacob Faibussowitsch {
817e24637baSBarry Smith   PetscInt        n, *idxout;
8185d0c19d7SBarry Smith   const PetscInt *idxin;
8193a40ed3dSBarry Smith 
8203a40ed3dSBarry Smith   PetscFunctionBegin;
8210700a824SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
8220700a824SBarry Smith   PetscValidHeaderSpecific(is, IS_CLASSID, 2);
8234f572ea9SToby Isaac   PetscAssertPointer(newis, 3);
82490f02eecSBarry Smith 
8259566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is, &n));
8269566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is, &idxin));
8279566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &idxout));
8289566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(mapping, n, idxin, idxout));
8299566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is, &idxin));
8309566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)is), n, idxout, PETSC_OWN_POINTER, newis));
8313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
83290f02eecSBarry Smith }
83390f02eecSBarry Smith 
834b89cb25eSSatish Balay /*@
8353acfe500SLois Curfman McInnes   ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering
8363acfe500SLois Curfman McInnes   and converts them to the global numbering.
83790f02eecSBarry Smith 
83820662ed9SBarry Smith   Not Collective
839b9cd556bSLois Curfman McInnes 
840bb25748dSBarry Smith   Input Parameters:
841b9cd556bSLois Curfman McInnes + mapping - the local to global mapping context
842bb25748dSBarry Smith . N       - number of integers
843b9cd556bSLois Curfman McInnes - in      - input indices in local numbering
844bb25748dSBarry Smith 
845bb25748dSBarry Smith   Output Parameter:
846bb25748dSBarry Smith . out - indices in global numbering
847bb25748dSBarry Smith 
848a997ad1aSLois Curfman McInnes   Level: advanced
849a997ad1aSLois Curfman McInnes 
850cab54364SBarry Smith   Note:
85120662ed9SBarry Smith   The `in` and `out` array parameters may be identical.
852cab54364SBarry Smith 
853cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApplyBlock()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingDestroy()`,
854c2e3fba1SPatrick Sanan           `ISLocalToGlobalMappingApplyIS()`, `AOCreateBasic()`, `AOApplicationToPetsc()`,
855db781477SPatrick Sanan           `AOPetscToApplication()`, `ISGlobalToLocalMappingApply()`
856afcb2eb5SJed Brown @*/
857d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping, PetscInt N, const PetscInt in[], PetscInt out[])
858d71ae5a4SJacob Faibussowitsch {
859cbc1caf0SMatthew G. Knepley   PetscInt i, bs, Nmax;
86045b6f7e9SBarry Smith 
86145b6f7e9SBarry Smith   PetscFunctionBegin;
862cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
863cbc1caf0SMatthew G. Knepley   bs   = mapping->bs;
864cbc1caf0SMatthew G. Knepley   Nmax = bs * mapping->n;
86545b6f7e9SBarry Smith   if (bs == 1) {
866cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
86745b6f7e9SBarry Smith     for (i = 0; i < N; i++) {
86845b6f7e9SBarry Smith       if (in[i] < 0) {
86945b6f7e9SBarry Smith         out[i] = in[i];
87045b6f7e9SBarry Smith         continue;
87145b6f7e9SBarry Smith       }
87208401ef6SPierre 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);
87345b6f7e9SBarry Smith       out[i] = idx[in[i]];
87445b6f7e9SBarry Smith     }
87545b6f7e9SBarry Smith   } else {
876cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
87745b6f7e9SBarry Smith     for (i = 0; i < N; i++) {
87845b6f7e9SBarry Smith       if (in[i] < 0) {
87945b6f7e9SBarry Smith         out[i] = in[i];
88045b6f7e9SBarry Smith         continue;
88145b6f7e9SBarry Smith       }
88208401ef6SPierre 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);
88345b6f7e9SBarry Smith       out[i] = idx[in[i] / bs] * bs + (in[i] % bs);
88445b6f7e9SBarry Smith     }
88545b6f7e9SBarry Smith   }
8863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
88745b6f7e9SBarry Smith }
88845b6f7e9SBarry Smith 
88945b6f7e9SBarry Smith /*@
8906006e8d2SBarry Smith   ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local block numbering and converts them to the global block numbering
89145b6f7e9SBarry Smith 
89220662ed9SBarry Smith   Not Collective
89345b6f7e9SBarry Smith 
89445b6f7e9SBarry Smith   Input Parameters:
89545b6f7e9SBarry Smith + mapping - the local to global mapping context
89645b6f7e9SBarry Smith . N       - number of integers
8976006e8d2SBarry Smith - in      - input indices in local block numbering
89845b6f7e9SBarry Smith 
89945b6f7e9SBarry Smith   Output Parameter:
9006006e8d2SBarry Smith . out - indices in global block numbering
90145b6f7e9SBarry Smith 
9026006e8d2SBarry Smith   Example:
903cab54364SBarry 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
9046006e8d2SBarry Smith   (the first block) would produce 0 and the mapping applied to 1 (the second block) would produce 3.
9056006e8d2SBarry Smith 
90620662ed9SBarry Smith   Level: advanced
90720662ed9SBarry Smith 
90820662ed9SBarry Smith   Note:
90920662ed9SBarry Smith   The `in` and `out` array parameters may be identical.
91020662ed9SBarry Smith 
911cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingDestroy()`,
912c2e3fba1SPatrick Sanan           `ISLocalToGlobalMappingApplyIS()`, `AOCreateBasic()`, `AOApplicationToPetsc()`,
913db781477SPatrick Sanan           `AOPetscToApplication()`, `ISGlobalToLocalMappingApply()`
91445b6f7e9SBarry Smith @*/
915d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping, PetscInt N, const PetscInt in[], PetscInt out[])
916d71ae5a4SJacob Faibussowitsch {
9178a1f772fSStefano Zampini   PetscInt        i, Nmax;
9188a1f772fSStefano Zampini   const PetscInt *idx;
919d4bb536fSBarry Smith 
920a0d79125SStefano Zampini   PetscFunctionBegin;
921a0d79125SStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
9228a1f772fSStefano Zampini   Nmax = mapping->n;
9238a1f772fSStefano Zampini   idx  = mapping->indices;
924afcb2eb5SJed Brown   for (i = 0; i < N; i++) {
925afcb2eb5SJed Brown     if (in[i] < 0) {
926afcb2eb5SJed Brown       out[i] = in[i];
927afcb2eb5SJed Brown       continue;
928afcb2eb5SJed Brown     }
92908401ef6SPierre 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);
930afcb2eb5SJed Brown     out[i] = idx[in[i]];
931afcb2eb5SJed Brown   }
9323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
933afcb2eb5SJed Brown }
934d4bb536fSBarry Smith 
9357e99dc12SLawrence Mitchell /*@
936a997ad1aSLois Curfman McInnes   ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers
937a997ad1aSLois Curfman McInnes   specified with a global numbering.
938d4bb536fSBarry Smith 
93920662ed9SBarry Smith   Not Collective
940b9cd556bSLois Curfman McInnes 
941d4bb536fSBarry Smith   Input Parameters:
942b9cd556bSLois Curfman McInnes + mapping - mapping between local and global numbering
943cab54364SBarry Smith . type    - `IS_GTOLM_MASK` - maps global indices with no local value to -1 in the output list (i.e., mask them)
944cab54364SBarry Smith            `IS_GTOLM_DROP` - drops the indices with no local value from the output list
945d4bb536fSBarry Smith . n       - number of global indices to map
946b9cd556bSLois Curfman McInnes - idx     - global indices to map
947d4bb536fSBarry Smith 
948d4bb536fSBarry Smith   Output Parameters:
949cab54364SBarry Smith + nout   - number of indices in output array (if type == `IS_GTOLM_MASK` then nout = n)
950b9cd556bSLois Curfman McInnes - idxout - local index of each global index, one must pass in an array long enough
951cab54364SBarry Smith              to hold all the indices. You can call `ISGlobalToLocalMappingApply()` with
9520298fd71SBarry Smith              idxout == NULL to determine the required length (returned in nout)
953cab54364SBarry Smith              and then allocate the required space and call `ISGlobalToLocalMappingApply()`
954e182c471SBarry Smith              a second time to set the values.
955d4bb536fSBarry Smith 
956cab54364SBarry Smith   Level: advanced
957cab54364SBarry Smith 
958b9cd556bSLois Curfman McInnes   Notes:
95920662ed9SBarry Smith   Either `nout` or `idxout` may be `NULL`. `idx` and `idxout` may be identical.
960d4bb536fSBarry Smith 
96120662ed9SBarry Smith   For "small" problems when using `ISGlobalToLocalMappingApply()` and `ISGlobalToLocalMappingApplyBlock()`, the `ISLocalToGlobalMappingType` of
96220662ed9SBarry Smith   `ISLOCALTOGLOBALMAPPINGBASIC` will be used;
963cab54364SBarry 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.
964cab54364SBarry Smith   Use `ISLocalToGlobalMappingSetType()` or call `ISLocalToGlobalMappingSetFromOptions()` with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
9650f5bd95cSBarry Smith 
96638b5cf2dSJacob Faibussowitsch   Developer Notes:
96720662ed9SBarry Smith   The manual page states that `idx` and `idxout` may be identical but the calling
96820662ed9SBarry Smith   sequence declares `idx` as const so it cannot be the same as `idxout`.
96932fd6b96SBarry Smith 
970cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApply()`, `ISGlobalToLocalMappingApplyBlock()`, `ISLocalToGlobalMappingCreate()`,
971db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`
972d4bb536fSBarry Smith @*/
973d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, PetscInt n, const PetscInt idx[], PetscInt *nout, PetscInt idxout[])
974d71ae5a4SJacob Faibussowitsch {
9759d90f715SBarry Smith   PetscFunctionBegin;
9769d90f715SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
97748a46eb9SPierre Jolivet   if (!mapping->data) PetscCall(ISGlobalToLocalMappingSetUp(mapping));
978dbbe0bcdSBarry Smith   PetscUseTypeMethod(mapping, globaltolocalmappingapply, type, n, idx, nout, idxout);
9793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9809d90f715SBarry Smith }
9819d90f715SBarry Smith 
982d4fe737eSStefano Zampini /*@
983cab54364SBarry Smith   ISGlobalToLocalMappingApplyIS - Creates from an `IS` in the global numbering
984cab54364SBarry Smith   a new index set using the local numbering defined in an `ISLocalToGlobalMapping`
985d4fe737eSStefano Zampini   context.
986d4fe737eSStefano Zampini 
98720662ed9SBarry Smith   Not Collective
988d4fe737eSStefano Zampini 
989d4fe737eSStefano Zampini   Input Parameters:
990d4fe737eSStefano Zampini + mapping - mapping between local and global numbering
991cab54364SBarry Smith . type    - `IS_GTOLM_MASK` - maps global indices with no local value to -1 in the output list (i.e., mask them)
992cab54364SBarry Smith            `IS_GTOLM_DROP` - drops the indices with no local value from the output list
993d4fe737eSStefano Zampini - is      - index set in global numbering
994d4fe737eSStefano Zampini 
9952fe279fdSBarry Smith   Output Parameter:
996d4fe737eSStefano Zampini . newis - index set in local numbering
997d4fe737eSStefano Zampini 
998d4fe737eSStefano Zampini   Level: advanced
999d4fe737eSStefano Zampini 
1000cab54364SBarry Smith   Note:
1001cab54364SBarry Smith   The output `IS` will be sequential, as it encodes a purely local operation
1002cab54364SBarry Smith 
1003cab54364SBarry Smith .seealso: [](sec_scatter), `ISGlobalToLocalMapping`, `ISGlobalToLocalMappingApply()`, `ISLocalToGlobalMappingCreate()`,
1004db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`
1005d4fe737eSStefano Zampini @*/
1006d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, IS is, IS *newis)
1007d71ae5a4SJacob Faibussowitsch {
1008d4fe737eSStefano Zampini   PetscInt        n, nout, *idxout;
1009d4fe737eSStefano Zampini   const PetscInt *idxin;
1010d4fe737eSStefano Zampini 
1011d4fe737eSStefano Zampini   PetscFunctionBegin;
1012d4fe737eSStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
1013d4fe737eSStefano Zampini   PetscValidHeaderSpecific(is, IS_CLASSID, 3);
10144f572ea9SToby Isaac   PetscAssertPointer(newis, 4);
1015d4fe737eSStefano Zampini 
10169566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is, &n));
10179566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is, &idxin));
1018d4fe737eSStefano Zampini   if (type == IS_GTOLM_MASK) {
10199566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &idxout));
1020d4fe737eSStefano Zampini   } else {
10219566063dSJacob Faibussowitsch     PetscCall(ISGlobalToLocalMappingApply(mapping, type, n, idxin, &nout, NULL));
10229566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nout, &idxout));
1023d4fe737eSStefano Zampini   }
10249566063dSJacob Faibussowitsch   PetscCall(ISGlobalToLocalMappingApply(mapping, type, n, idxin, &nout, idxout));
10259566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is, &idxin));
10269566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nout, idxout, PETSC_OWN_POINTER, newis));
10273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1028d4fe737eSStefano Zampini }
1029d4fe737eSStefano Zampini 
10309d90f715SBarry Smith /*@
10319d90f715SBarry Smith   ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers
10329d90f715SBarry Smith   specified with a block global numbering.
10339d90f715SBarry Smith 
103420662ed9SBarry Smith   Not Collective
10359d90f715SBarry Smith 
10369d90f715SBarry Smith   Input Parameters:
10379d90f715SBarry Smith + mapping - mapping between local and global numbering
1038cab54364SBarry Smith . type    - `IS_GTOLM_MASK` - maps global indices with no local value to -1 in the output list (i.e., mask them)
1039cab54364SBarry Smith            `IS_GTOLM_DROP` - drops the indices with no local value from the output list
10409d90f715SBarry Smith . n       - number of global indices to map
10419d90f715SBarry Smith - idx     - global indices to map
10429d90f715SBarry Smith 
10439d90f715SBarry Smith   Output Parameters:
1044cab54364SBarry Smith + nout   - number of indices in output array (if type == `IS_GTOLM_MASK` then nout = n)
10459d90f715SBarry Smith - idxout - local index of each global index, one must pass in an array long enough
1046cab54364SBarry Smith              to hold all the indices. You can call `ISGlobalToLocalMappingApplyBlock()` with
10479d90f715SBarry Smith              idxout == NULL to determine the required length (returned in nout)
1048cab54364SBarry Smith              and then allocate the required space and call `ISGlobalToLocalMappingApplyBlock()`
10499d90f715SBarry Smith              a second time to set the values.
10509d90f715SBarry Smith 
1051cab54364SBarry Smith   Level: advanced
1052cab54364SBarry Smith 
10539d90f715SBarry Smith   Notes:
105420662ed9SBarry Smith   Either `nout` or `idxout` may be `NULL`. `idx` and `idxout` may be identical.
10559d90f715SBarry Smith 
105620662ed9SBarry Smith   For "small" problems when using `ISGlobalToLocalMappingApply()` and `ISGlobalToLocalMappingApplyBlock()`, the `ISLocalToGlobalMappingType` of
105720662ed9SBarry Smith   `ISLOCALTOGLOBALMAPPINGBASIC` will be used;
1058cab54364SBarry 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.
1059cab54364SBarry Smith   Use `ISLocalToGlobalMappingSetType()` or call `ISLocalToGlobalMappingSetFromOptions()` with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
10609d90f715SBarry Smith 
106138b5cf2dSJacob Faibussowitsch   Developer Notes:
106220662ed9SBarry Smith   The manual page states that `idx` and `idxout` may be identical but the calling
106320662ed9SBarry Smith   sequence declares `idx` as const so it cannot be the same as `idxout`.
10649d90f715SBarry Smith 
1065cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApply()`, `ISGlobalToLocalMappingApply()`, `ISLocalToGlobalMappingCreate()`,
1066db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`
10679d90f715SBarry Smith @*/
1068d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, PetscInt n, const PetscInt idx[], PetscInt *nout, PetscInt idxout[])
1069d71ae5a4SJacob Faibussowitsch {
10703a40ed3dSBarry Smith   PetscFunctionBegin;
10710700a824SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
107248a46eb9SPierre Jolivet   if (!mapping->data) PetscCall(ISGlobalToLocalMappingSetUp(mapping));
1073dbbe0bcdSBarry Smith   PetscUseTypeMethod(mapping, globaltolocalmappingapplyblock, type, n, idx, nout, idxout);
10743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1075d4bb536fSBarry Smith }
107690f02eecSBarry Smith 
107789d82c54SBarry Smith /*@C
1078633354d9SStefano Zampini   ISLocalToGlobalMappingGetBlockInfo - Gets the neighbor information
107989d82c54SBarry Smith 
1080633354d9SStefano Zampini   Collective the first time it is called
108189d82c54SBarry Smith 
1082f899ff85SJose E. Roman   Input Parameter:
108389d82c54SBarry Smith . mapping - the mapping from local to global indexing
108489d82c54SBarry Smith 
1085d8d19677SJose E. Roman   Output Parameters:
1086633354d9SStefano Zampini + nproc    - number of processes that are connected to the calling process
1087633354d9SStefano Zampini . procs    - neighboring processes
1088633354d9SStefano Zampini . numprocs - number of block indices for each process
1089633354d9SStefano Zampini - indices  - block indices (in local numbering) shared with neighbors (sorted by global numbering)
109089d82c54SBarry Smith 
109189d82c54SBarry Smith   Level: advanced
109289d82c54SBarry Smith 
1093cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1094633354d9SStefano Zampini           `ISLocalToGlobalMappingRestoreBlockInfo()`
109589d82c54SBarry Smith @*/
1096d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[])
1097d71ae5a4SJacob Faibussowitsch {
1098268a049cSStefano Zampini   PetscFunctionBegin;
1099268a049cSStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
1100633354d9SStefano Zampini   PetscCall(ISLocalToGlobalMappingSetUpBlockInfo_Private(mapping));
1101633354d9SStefano Zampini   if (nproc) *nproc = mapping->info_nproc;
1102633354d9SStefano Zampini   if (procs) *procs = mapping->info_procs;
1103633354d9SStefano Zampini   if (numprocs) *numprocs = mapping->info_numprocs;
1104633354d9SStefano Zampini   if (indices) *indices = mapping->info_indices;
11053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1106268a049cSStefano Zampini }
1107268a049cSStefano Zampini 
1108633354d9SStefano Zampini /*@C
1109633354d9SStefano Zampini   ISLocalToGlobalMappingGetBlockNodeInfo - Gets the neighbor information for each local block index
1110633354d9SStefano Zampini 
1111633354d9SStefano Zampini   Collective the first time it is called
1112633354d9SStefano Zampini 
1113633354d9SStefano Zampini   Input Parameter:
1114633354d9SStefano Zampini . mapping - the mapping from local to global indexing
1115633354d9SStefano Zampini 
1116633354d9SStefano Zampini   Output Parameter:
1117633354d9SStefano Zampini + n       - number of local block nodes
1118633354d9SStefano Zampini . n_procs - an array storing the number of processes for each local block node (including self)
1119633354d9SStefano Zampini - procs   - the processes' rank for each local block node (sorted, self is first)
1120633354d9SStefano Zampini 
1121633354d9SStefano Zampini   Level: advanced
1122633354d9SStefano Zampini 
1123633354d9SStefano Zampini   Notes:
1124633354d9SStefano Zampini   The user needs to call `ISLocalToGlobalMappingRestoreBlockNodeInfo()` when the data is no longer needed.
1125633354d9SStefano Zampini   The information returned by this function complements that of `ISLocalToGlobalMappingGetBlockInfo()`.
1126633354d9SStefano Zampini   The latter only provides local information, and the neighboring information
1127633354d9SStefano Zampini   cannot be inferred in the general case, unless the mapping is locally one-to-one on each process.
1128633354d9SStefano Zampini 
1129633354d9SStefano Zampini .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1130633354d9SStefano Zampini           `ISLocalToGlobalMappingGetBlockInfo()`, `ISLocalToGlobalMappingRestoreBlockNodeInfo()`, `ISLocalToGlobalMappingGetNodeInfo()`
1131633354d9SStefano Zampini @*/
1132633354d9SStefano Zampini PetscErrorCode ISLocalToGlobalMappingGetBlockNodeInfo(ISLocalToGlobalMapping mapping, PetscInt *n, PetscInt *n_procs[], PetscInt **procs[])
1133d71ae5a4SJacob Faibussowitsch {
1134633354d9SStefano Zampini   PetscFunctionBegin;
1135633354d9SStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
1136633354d9SStefano Zampini   PetscCall(ISLocalToGlobalMappingSetUpBlockInfo_Private(mapping));
1137633354d9SStefano Zampini   if (n) *n = mapping->n;
1138633354d9SStefano Zampini   if (n_procs) *n_procs = mapping->info_nodec;
1139633354d9SStefano Zampini   if (procs) *procs = mapping->info_nodei;
1140633354d9SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1141633354d9SStefano Zampini }
1142633354d9SStefano Zampini 
1143633354d9SStefano Zampini /*@C
1144633354d9SStefano Zampini   ISLocalToGlobalMappingRestoreBlockNodeInfo - Frees the memory allocated by `ISLocalToGlobalMappingGetBlockNodeInfo()`
1145633354d9SStefano Zampini 
1146633354d9SStefano Zampini   Not Collective
1147633354d9SStefano Zampini 
1148633354d9SStefano Zampini   Input Parameters:
1149633354d9SStefano Zampini + mapping - the mapping from local to global indexing
1150633354d9SStefano Zampini . n       - number of local block nodes
1151633354d9SStefano Zampini . n_procs - an array storing the number of processes for each local block nodes (including self)
1152633354d9SStefano Zampini - procs   - the processes' rank for each local block node (sorted, self is first)
1153633354d9SStefano Zampini 
1154633354d9SStefano Zampini   Level: advanced
1155633354d9SStefano Zampini 
1156633354d9SStefano Zampini .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1157633354d9SStefano Zampini           `ISLocalToGlobalMappingGetBlockNodeInfo()`
1158633354d9SStefano Zampini @*/
1159633354d9SStefano Zampini PetscErrorCode ISLocalToGlobalMappingRestoreBlockNodeInfo(ISLocalToGlobalMapping mapping, PetscInt *n, PetscInt *n_procs[], PetscInt **procs[])
1160633354d9SStefano Zampini {
1161633354d9SStefano Zampini   PetscFunctionBegin;
1162633354d9SStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
1163633354d9SStefano Zampini   if (n) *n = 0;
1164633354d9SStefano Zampini   if (n_procs) *n_procs = NULL;
1165633354d9SStefano Zampini   if (procs) *procs = NULL;
1166633354d9SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1167633354d9SStefano Zampini }
1168633354d9SStefano Zampini 
1169633354d9SStefano Zampini static PetscErrorCode ISLocalToGlobalMappingSetUpBlockInfo_Private(ISLocalToGlobalMapping mapping)
1170633354d9SStefano Zampini {
1171633354d9SStefano Zampini   PetscSF            sf;
1172ce94432eSBarry Smith   MPI_Comm           comm;
1173633354d9SStefano Zampini   const PetscSFNode *sfnode;
1174633354d9SStefano Zampini   PetscSFNode       *newsfnode;
1175633354d9SStefano Zampini   PetscLayout        layout;
1176633354d9SStefano Zampini   PetscHMapI         neighs;
1177633354d9SStefano Zampini   PetscHashIter      iter;
1178633354d9SStefano Zampini   PetscBool          missing;
1179633354d9SStefano Zampini   const PetscInt    *gidxs, *rootdegree;
1180633354d9SStefano Zampini   PetscInt          *mask, *mrootdata, *leafdata, *newleafdata, *leafrd, *tmpg;
1181633354d9SStefano Zampini   PetscInt           nroots, nleaves, newnleaves, bs, i, j, m, mnroots, p;
1182633354d9SStefano Zampini   PetscMPIInt        rank, size;
118389d82c54SBarry Smith 
118489d82c54SBarry Smith   PetscFunctionBegin;
1185633354d9SStefano Zampini   if (mapping->info_numprocs) PetscFunctionReturn(PETSC_SUCCESS);
11869566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mapping, &comm));
11879566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
11889566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1189633354d9SStefano Zampini 
1190633354d9SStefano Zampini   /* Get mapping indices */
1191633354d9SStefano Zampini   PetscCall(ISLocalToGlobalMappingGetBlockSize(mapping, &bs));
1192633354d9SStefano Zampini   PetscCall(ISLocalToGlobalMappingGetBlockIndices(mapping, &gidxs));
1193633354d9SStefano Zampini   PetscCall(ISLocalToGlobalMappingGetSize(mapping, &nleaves));
1194633354d9SStefano Zampini   nleaves /= bs;
1195633354d9SStefano Zampini 
1196633354d9SStefano Zampini   /* Create layout for global indices */
1197633354d9SStefano Zampini   for (i = 0, m = 0; i < nleaves; i++) m = PetscMax(m, gidxs[i]);
1198633354d9SStefano Zampini   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &m, 1, MPIU_INT, MPI_MAX, comm));
1199633354d9SStefano Zampini   PetscCall(PetscLayoutCreate(comm, &layout));
1200633354d9SStefano Zampini   PetscCall(PetscLayoutSetSize(layout, m + 1));
1201633354d9SStefano Zampini   PetscCall(PetscLayoutSetUp(layout));
1202633354d9SStefano Zampini 
1203633354d9SStefano Zampini   /* Create SF to share global indices */
1204633354d9SStefano Zampini   PetscCall(PetscSFCreate(comm, &sf));
1205633354d9SStefano Zampini   PetscCall(PetscSFSetGraphLayout(sf, layout, nleaves, NULL, PETSC_OWN_POINTER, gidxs));
1206633354d9SStefano Zampini   PetscCall(PetscSFSetUp(sf));
1207633354d9SStefano Zampini   PetscCall(PetscLayoutDestroy(&layout));
1208633354d9SStefano Zampini 
1209633354d9SStefano Zampini   /* communicate root degree to leaves */
1210633354d9SStefano Zampini   PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, &sfnode));
1211633354d9SStefano Zampini   PetscCall(PetscSFComputeDegreeBegin(sf, &rootdegree));
1212633354d9SStefano Zampini   PetscCall(PetscSFComputeDegreeEnd(sf, &rootdegree));
1213633354d9SStefano Zampini   for (i = 0, mnroots = 0; i < nroots; i++) mnroots += rootdegree[i];
1214633354d9SStefano Zampini   PetscCall(PetscMalloc3(2 * PetscMax(mnroots, nroots), &mrootdata, 2 * nleaves, &leafdata, nleaves, &leafrd));
1215633354d9SStefano Zampini   for (i = 0, m = 0; i < nroots; i++) {
1216633354d9SStefano Zampini     mrootdata[2 * i + 0] = rootdegree[i];
1217633354d9SStefano Zampini     mrootdata[2 * i + 1] = m;
1218633354d9SStefano Zampini     m += rootdegree[i];
1219633354d9SStefano Zampini   }
1220633354d9SStefano Zampini   PetscCall(PetscSFBcastBegin(sf, MPIU_2INT, mrootdata, leafdata, MPI_REPLACE));
1221633354d9SStefano Zampini   PetscCall(PetscSFBcastEnd(sf, MPIU_2INT, mrootdata, leafdata, MPI_REPLACE));
1222633354d9SStefano Zampini 
1223633354d9SStefano Zampini   /* allocate enough space to store ranks */
1224633354d9SStefano Zampini   for (i = 0, newnleaves = 0; i < nleaves; i++) {
1225633354d9SStefano Zampini     newnleaves += leafdata[2 * i];
1226633354d9SStefano Zampini     leafrd[i] = leafdata[2 * i];
122724cf384cSBarry Smith   }
122824cf384cSBarry Smith 
1229633354d9SStefano Zampini   /* create new SF nodes to collect multi-root data at leaves */
1230633354d9SStefano Zampini   PetscCall(PetscMalloc1(newnleaves, &newsfnode));
1231633354d9SStefano Zampini   for (i = 0, m = 0; i < nleaves; i++) {
1232633354d9SStefano Zampini     for (j = 0; j < leafrd[i]; j++) {
1233633354d9SStefano Zampini       newsfnode[m].rank  = sfnode[i].rank;
1234633354d9SStefano Zampini       newsfnode[m].index = leafdata[2 * i + 1] + j;
1235633354d9SStefano Zampini       m++;
1236bc8ff85bSBarry Smith     }
1237bc8ff85bSBarry Smith   }
1238bc8ff85bSBarry Smith 
1239633354d9SStefano Zampini   /* gather ranks at multi roots */
1240633354d9SStefano Zampini   for (i = 0; i < mnroots; i++) mrootdata[i] = -1;
1241633354d9SStefano Zampini   for (i = 0; i < nleaves; i++) leafdata[i] = (PetscInt)rank;
124230dcb7c9SBarry Smith 
1243633354d9SStefano Zampini   PetscCall(PetscSFGatherBegin(sf, MPIU_INT, leafdata, mrootdata));
1244633354d9SStefano Zampini   PetscCall(PetscSFGatherEnd(sf, MPIU_INT, leafdata, mrootdata));
12453677ff5aSBarry Smith 
1246633354d9SStefano Zampini   /* set new multi-leaves graph into the SF */
1247633354d9SStefano Zampini   PetscCall(PetscSFSetGraph(sf, mnroots, newnleaves, NULL, PETSC_OWN_POINTER, newsfnode, PETSC_OWN_POINTER));
1248633354d9SStefano Zampini   PetscCall(PetscSFSetUp(sf));
1249f6e5521dSKarl Rupp 
1250633354d9SStefano Zampini   /* broadcast multi-root data to multi-leaves */
1251633354d9SStefano Zampini   PetscCall(PetscMalloc1(newnleaves, &newleafdata));
1252633354d9SStefano Zampini   PetscCall(PetscSFBcastBegin(sf, MPIU_INT, mrootdata, newleafdata, MPI_REPLACE));
1253633354d9SStefano Zampini   PetscCall(PetscSFBcastEnd(sf, MPIU_INT, mrootdata, newleafdata, MPI_REPLACE));
1254f6e5521dSKarl Rupp 
1255633354d9SStefano Zampini   /* sort sharing ranks */
1256633354d9SStefano Zampini   for (i = 0, m = 0; i < nleaves; i++) {
1257633354d9SStefano Zampini     PetscCall(PetscSortInt(leafrd[i], newleafdata + m));
1258633354d9SStefano Zampini     m += leafrd[i];
125930dcb7c9SBarry Smith   }
126030dcb7c9SBarry Smith 
1261633354d9SStefano Zampini   /* Number of neighbors and their ranks */
1262633354d9SStefano Zampini   PetscCall(PetscHMapICreate(&neighs));
1263633354d9SStefano Zampini   for (i = 0; i < newnleaves; i++) PetscCall(PetscHMapIPut(neighs, newleafdata[i], &iter, &missing));
1264633354d9SStefano Zampini   PetscCall(PetscHMapIGetSize(neighs, &mapping->info_nproc));
1265633354d9SStefano Zampini   PetscCall(PetscMalloc1(mapping->info_nproc + 1, &mapping->info_procs));
1266633354d9SStefano Zampini   PetscCall(PetscHMapIGetKeys(neighs, (i = 0, &i), mapping->info_procs));
1267633354d9SStefano Zampini   for (i = 0; i < mapping->info_nproc; i++) { /* put info for self first */
1268633354d9SStefano Zampini     if (mapping->info_procs[i] == rank) {
1269633354d9SStefano Zampini       PetscInt newr = mapping->info_procs[0];
1270d44834fbSBarry Smith 
1271633354d9SStefano Zampini       mapping->info_procs[0] = rank;
1272633354d9SStefano Zampini       mapping->info_procs[i] = newr;
127324cf384cSBarry Smith       break;
127424cf384cSBarry Smith     }
127524cf384cSBarry Smith   }
1276633354d9SStefano Zampini   if (mapping->info_nproc) PetscCall(PetscSortInt(mapping->info_nproc - 1, mapping->info_procs + 1));
1277633354d9SStefano Zampini   PetscCall(PetscHMapIDestroy(&neighs));
1278268a049cSStefano Zampini 
1279633354d9SStefano Zampini   /* collect info data */
1280633354d9SStefano Zampini   PetscCall(PetscMalloc1(mapping->info_nproc + 1, &mapping->info_numprocs));
1281633354d9SStefano Zampini   PetscCall(PetscMalloc1(mapping->info_nproc + 1, &mapping->info_indices));
1282633354d9SStefano Zampini   for (i = 0; i < mapping->info_nproc + 1; i++) mapping->info_indices[i] = NULL;
1283633354d9SStefano Zampini 
1284633354d9SStefano Zampini   PetscCall(PetscMalloc1(nleaves, &mask));
1285633354d9SStefano Zampini   PetscCall(PetscMalloc1(nleaves, &tmpg));
1286633354d9SStefano Zampini   for (p = 0; p < mapping->info_nproc; p++) {
1287633354d9SStefano Zampini     PetscInt *tmp, trank = mapping->info_procs[p];
1288633354d9SStefano Zampini 
1289633354d9SStefano Zampini     PetscCall(PetscMemzero(mask, nleaves * sizeof(*mask)));
1290633354d9SStefano Zampini     for (i = 0, m = 0; i < nleaves; i++) {
1291633354d9SStefano Zampini       for (j = 0; j < leafrd[i]; j++) {
1292633354d9SStefano Zampini         if (newleafdata[m] == trank) mask[i]++;
1293633354d9SStefano Zampini         if (!p && newleafdata[m] != rank) mask[i]++;
1294633354d9SStefano Zampini         m++;
1295633354d9SStefano Zampini       }
1296633354d9SStefano Zampini     }
1297633354d9SStefano Zampini     for (i = 0, m = 0; i < nleaves; i++)
1298633354d9SStefano Zampini       if (mask[i] > (!p ? 1 : 0)) m++;
1299633354d9SStefano Zampini 
1300633354d9SStefano Zampini     PetscCall(PetscMalloc1(m, &tmp));
1301633354d9SStefano Zampini     for (i = 0, m = 0; i < nleaves; i++)
1302633354d9SStefano Zampini       if (mask[i] > (!p ? 1 : 0)) {
1303633354d9SStefano Zampini         tmp[m]  = i;
1304633354d9SStefano Zampini         tmpg[m] = gidxs[i];
1305633354d9SStefano Zampini         m++;
1306633354d9SStefano Zampini       }
1307633354d9SStefano Zampini     PetscCall(PetscSortIntWithArray(m, tmpg, tmp));
1308633354d9SStefano Zampini     mapping->info_indices[p]  = tmp;
1309633354d9SStefano Zampini     mapping->info_numprocs[p] = m;
1310633354d9SStefano Zampini   }
1311633354d9SStefano Zampini 
1312633354d9SStefano Zampini   /* Node info */
1313633354d9SStefano Zampini   PetscCall(PetscMalloc2(nleaves, &mapping->info_nodec, nleaves + 1, &mapping->info_nodei));
1314633354d9SStefano Zampini   PetscCall(PetscArraycpy(mapping->info_nodec, leafrd, nleaves));
1315633354d9SStefano Zampini   PetscCall(PetscMalloc1(newnleaves, &mapping->info_nodei[0]));
1316633354d9SStefano Zampini   for (i = 0; i < nleaves - 1; i++) mapping->info_nodei[i + 1] = mapping->info_nodei[i] + mapping->info_nodec[i];
1317633354d9SStefano Zampini   PetscCall(PetscArraycpy(mapping->info_nodei[0], newleafdata, newnleaves));
1318633354d9SStefano Zampini 
1319633354d9SStefano Zampini   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(mapping, &gidxs));
1320633354d9SStefano Zampini   PetscCall(PetscFree(tmpg));
1321633354d9SStefano Zampini   PetscCall(PetscFree(mask));
1322633354d9SStefano Zampini   PetscCall(PetscSFDestroy(&sf));
1323633354d9SStefano Zampini   PetscCall(PetscFree3(mrootdata, leafdata, leafrd));
1324633354d9SStefano Zampini   PetscCall(PetscFree(newleafdata));
13253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
132689d82c54SBarry Smith }
132789d82c54SBarry Smith 
13286a818285SBarry Smith /*@C
1329cab54364SBarry Smith   ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by `ISLocalToGlobalMappingGetBlockInfo()`
13306a818285SBarry Smith 
1331633354d9SStefano Zampini   Not Collective
13326a818285SBarry Smith 
1333633354d9SStefano Zampini   Input Parameters:
1334633354d9SStefano Zampini + mapping  - the mapping from local to global indexing
1335633354d9SStefano Zampini . nproc    - number of processes that are connected to the calling process
1336633354d9SStefano Zampini . procs    - neighboring processes
1337633354d9SStefano Zampini . numprocs - number of block indices for each process
1338633354d9SStefano Zampini - indices  - block indices (in local numbering) shared with neighbors (sorted by global numbering)
13396a818285SBarry Smith 
13406a818285SBarry Smith   Level: advanced
13416a818285SBarry Smith 
1342cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1343db781477SPatrick Sanan           `ISLocalToGlobalMappingGetInfo()`
13446a818285SBarry Smith @*/
1345d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[])
1346d71ae5a4SJacob Faibussowitsch {
13476a818285SBarry Smith   PetscFunctionBegin;
1348cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
1349633354d9SStefano Zampini   if (nproc) *nproc = 0;
1350633354d9SStefano Zampini   if (procs) *procs = NULL;
1351633354d9SStefano Zampini   if (numprocs) *numprocs = NULL;
1352633354d9SStefano Zampini   if (indices) *indices = NULL;
13533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13546a818285SBarry Smith }
13556a818285SBarry Smith 
13566a818285SBarry Smith /*@C
1357633354d9SStefano Zampini   ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each process
13586a818285SBarry Smith 
1359633354d9SStefano Zampini   Collective the first time it is called
13606a818285SBarry Smith 
1361f899ff85SJose E. Roman   Input Parameter:
13626a818285SBarry Smith . mapping - the mapping from local to global indexing
13636a818285SBarry Smith 
1364d8d19677SJose E. Roman   Output Parameters:
1365633354d9SStefano Zampini + nproc    - number of processes that are connected to the calling process
1366633354d9SStefano Zampini . procs    - neighboring processes
1367633354d9SStefano Zampini . numprocs - number of indices for each process
1368633354d9SStefano Zampini - indices  - indices (in local numbering) shared with neighbors (sorted by global numbering)
13696a818285SBarry Smith 
13706a818285SBarry Smith   Level: advanced
13716a818285SBarry Smith 
1372cab54364SBarry Smith   Note:
1373cab54364SBarry Smith   The user needs to call `ISLocalToGlobalMappingRestoreInfo()` when the data is no longer needed.
13741bd0b88eSStefano Zampini 
137538b5cf2dSJacob Faibussowitsch   Fortran Notes:
1376e33f79d8SJacob Faibussowitsch   There is no `ISLocalToGlobalMappingRestoreInfo()` in Fortran. You must make sure that
1377e33f79d8SJacob Faibussowitsch   `procs`[], `numprocs`[] and `indices`[][] are large enough arrays, either by allocating them
1378e33f79d8SJacob Faibussowitsch   dynamically or defining static ones large enough.
1379e33f79d8SJacob Faibussowitsch 
1380cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1381633354d9SStefano Zampini           `ISLocalToGlobalMappingRestoreInfo()`, `ISLocalToGlobalMappingGetNodeInfo()`
13826a818285SBarry Smith @*/
1383d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[])
1384d71ae5a4SJacob Faibussowitsch {
1385633354d9SStefano Zampini   PetscInt **bindices = NULL, *bnumprocs = NULL, bs, i, j, k, n, *bprocs;
13866a818285SBarry Smith 
13876a818285SBarry Smith   PetscFunctionBegin;
13886a818285SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
13898a1f772fSStefano Zampini   bs = mapping->bs;
1390633354d9SStefano Zampini   PetscCall(ISLocalToGlobalMappingGetBlockInfo(mapping, &n, &bprocs, &bnumprocs, &bindices));
1391268a049cSStefano Zampini   if (bs > 1) { /* we need to expand the cached info */
1392633354d9SStefano Zampini     if (indices) PetscCall(PetscCalloc1(n, indices));
1393633354d9SStefano Zampini     if (numprocs) PetscCall(PetscCalloc1(n, numprocs));
1394633354d9SStefano Zampini     if (indices || numprocs) {
1395633354d9SStefano Zampini       for (i = 0; i < n; i++) {
1396633354d9SStefano Zampini         if (indices) {
13979566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(bs * bnumprocs[i], &(*indices)[i]));
1398268a049cSStefano Zampini           for (j = 0; j < bnumprocs[i]; j++) {
1399ad540459SPierre Jolivet             for (k = 0; k < bs; k++) (*indices)[i][j * bs + k] = bs * bindices[i][j] + k;
14006a818285SBarry Smith           }
14016a818285SBarry Smith         }
1402633354d9SStefano Zampini         if (numprocs) (*numprocs)[i] = bnumprocs[i] * bs;
1403633354d9SStefano Zampini       }
1404633354d9SStefano Zampini     }
1405268a049cSStefano Zampini   } else {
1406633354d9SStefano Zampini     if (numprocs) *numprocs = bnumprocs;
1407633354d9SStefano Zampini     if (indices) *indices = bindices;
14086a818285SBarry Smith   }
1409633354d9SStefano Zampini   if (nproc) *nproc = n;
1410633354d9SStefano Zampini   if (procs) *procs = bprocs;
14113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14126a818285SBarry Smith }
14136a818285SBarry Smith 
141407b52d57SBarry Smith /*@C
1415cab54364SBarry Smith   ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by `ISLocalToGlobalMappingGetInfo()`
141689d82c54SBarry Smith 
1417633354d9SStefano Zampini   Not Collective
141807b52d57SBarry Smith 
1419633354d9SStefano Zampini   Input Parameters:
1420633354d9SStefano Zampini + mapping  - the mapping from local to global indexing
1421633354d9SStefano Zampini . nproc    - number of processes that are connected to the calling process
1422633354d9SStefano Zampini . procs    - neighboring processes
1423633354d9SStefano Zampini . numprocs - number of indices for each process
1424633354d9SStefano Zampini - indices  - indices (in local numbering) shared with neighbors (sorted by global numbering)
142507b52d57SBarry Smith 
142607b52d57SBarry Smith   Level: advanced
142707b52d57SBarry Smith 
1428cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1429db781477SPatrick Sanan           `ISLocalToGlobalMappingGetInfo()`
143007b52d57SBarry Smith @*/
1431d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[])
1432d71ae5a4SJacob Faibussowitsch {
143307b52d57SBarry Smith   PetscFunctionBegin;
1434633354d9SStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
1435633354d9SStefano Zampini   if (mapping->bs > 1) {
1436633354d9SStefano Zampini     if (numprocs) PetscCall(PetscFree(*numprocs));
1437633354d9SStefano Zampini     if (indices) {
1438633354d9SStefano Zampini       if (*indices)
1439633354d9SStefano Zampini         for (PetscInt i = 0; i < *nproc; i++) PetscCall(PetscFree((*indices)[i]));
1440633354d9SStefano Zampini       PetscCall(PetscFree(*indices));
1441633354d9SStefano Zampini     }
1442633354d9SStefano Zampini   }
14433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
144407b52d57SBarry Smith }
144586994e45SJed Brown 
144686994e45SJed Brown /*@C
1447633354d9SStefano Zampini   ISLocalToGlobalMappingGetNodeInfo - Gets the neighbor information of local nodes
14481bd0b88eSStefano Zampini 
1449633354d9SStefano Zampini   Collective the first time it is called
14501bd0b88eSStefano Zampini 
1451f899ff85SJose E. Roman   Input Parameter:
14521bd0b88eSStefano Zampini . mapping - the mapping from local to global indexing
14531bd0b88eSStefano Zampini 
1454d8d19677SJose E. Roman   Output Parameters:
1455633354d9SStefano Zampini + n       - number of local nodes
1456633354d9SStefano Zampini . n_procs - an array storing the number of processes for each local node (including self)
1457633354d9SStefano Zampini - procs   - the processes' rank for each local node (sorted, self is first)
14581bd0b88eSStefano Zampini 
14591bd0b88eSStefano Zampini   Level: advanced
14601bd0b88eSStefano Zampini 
1461cab54364SBarry Smith   Note:
1462633354d9SStefano Zampini   The user needs to call `ISLocalToGlobalMappingRestoreNodeInfo()` when the data is no longer needed.
14631bd0b88eSStefano Zampini 
1464cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1465633354d9SStefano Zampini           `ISLocalToGlobalMappingGetInfo()`, `ISLocalToGlobalMappingRestoreNodeInfo()`, `ISLocalToGlobalMappingGetBlockNodeInfo()`
14661bd0b88eSStefano Zampini @*/
1467633354d9SStefano Zampini PetscErrorCode ISLocalToGlobalMappingGetNodeInfo(ISLocalToGlobalMapping mapping, PetscInt *n, PetscInt *n_procs[], PetscInt **procs[])
1468d71ae5a4SJacob Faibussowitsch {
1469633354d9SStefano Zampini   PetscInt **bprocs = NULL, *bn_procs = NULL, bs, i, j, k, bn;
14701bd0b88eSStefano Zampini 
14711bd0b88eSStefano Zampini   PetscFunctionBegin;
14721bd0b88eSStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
1473633354d9SStefano Zampini   bs = mapping->bs;
1474633354d9SStefano Zampini   PetscCall(ISLocalToGlobalMappingGetBlockNodeInfo(mapping, &bn, &bn_procs, &bprocs));
1475633354d9SStefano Zampini   if (bs > 1) { /* we need to expand the cached info */
1476633354d9SStefano Zampini     PetscInt *tn_procs;
1477633354d9SStefano Zampini     PetscInt  c;
14781bd0b88eSStefano Zampini 
1479633354d9SStefano Zampini     PetscCall(PetscMalloc1(bn * bs, &tn_procs));
1480633354d9SStefano Zampini     for (i = 0, c = 0; i < bn; i++) {
1481633354d9SStefano Zampini       for (k = 0; k < bs; k++) tn_procs[i * bs + k] = bn_procs[i];
1482633354d9SStefano Zampini       c += bs * bn_procs[i];
1483633354d9SStefano Zampini     }
1484633354d9SStefano Zampini     if (n) *n = bn * bs;
1485633354d9SStefano Zampini     if (procs) {
1486633354d9SStefano Zampini       PetscInt **tprocs;
1487633354d9SStefano Zampini       PetscInt   tn = bn * bs;
14881bd0b88eSStefano Zampini 
1489633354d9SStefano Zampini       PetscCall(PetscMalloc1(tn, &tprocs));
1490633354d9SStefano Zampini       if (tn) PetscCall(PetscMalloc1(c, &tprocs[0]));
1491633354d9SStefano Zampini       for (i = 0; i < tn - 1; i++) tprocs[i + 1] = tprocs[i] + tn_procs[i];
1492633354d9SStefano Zampini       for (i = 0; i < bn; i++) {
1493633354d9SStefano Zampini         for (k = 0; k < bs; k++) {
1494633354d9SStefano Zampini           for (j = 0; j < bn_procs[i]; j++) tprocs[i * bs + k][j] = bprocs[i][j];
14951bd0b88eSStefano Zampini         }
14961bd0b88eSStefano Zampini       }
1497633354d9SStefano Zampini       *procs = tprocs;
14981bd0b88eSStefano Zampini     }
1499633354d9SStefano Zampini     if (n_procs) *n_procs = tn_procs;
1500633354d9SStefano Zampini     else PetscCall(PetscFree(tn_procs));
1501633354d9SStefano Zampini   } else {
1502633354d9SStefano Zampini     if (n) *n = bn;
1503633354d9SStefano Zampini     if (n_procs) *n_procs = bn_procs;
1504633354d9SStefano Zampini     if (procs) *procs = bprocs;
1505633354d9SStefano Zampini   }
15063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15071bd0b88eSStefano Zampini }
15081bd0b88eSStefano Zampini 
15091bd0b88eSStefano Zampini /*@C
1510cab54364SBarry Smith   ISLocalToGlobalMappingRestoreNodeInfo - Frees the memory allocated by `ISLocalToGlobalMappingGetNodeInfo()`
15111bd0b88eSStefano Zampini 
1512633354d9SStefano Zampini   Not Collective
15131bd0b88eSStefano Zampini 
1514633354d9SStefano Zampini   Input Parameters:
1515633354d9SStefano Zampini + mapping - the mapping from local to global indexing
1516633354d9SStefano Zampini . n       - number of local nodes
1517633354d9SStefano Zampini . n_procs - an array storing the number of processes for each local node (including self)
1518633354d9SStefano Zampini - procs   - the processes' rank for each local node (sorted, self is first)
15191bd0b88eSStefano Zampini 
15201bd0b88eSStefano Zampini   Level: advanced
15211bd0b88eSStefano Zampini 
1522cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1523db781477SPatrick Sanan           `ISLocalToGlobalMappingGetInfo()`
15241bd0b88eSStefano Zampini @*/
1525633354d9SStefano Zampini PetscErrorCode ISLocalToGlobalMappingRestoreNodeInfo(ISLocalToGlobalMapping mapping, PetscInt *n, PetscInt *n_procs[], PetscInt **procs[])
1526d71ae5a4SJacob Faibussowitsch {
15271bd0b88eSStefano Zampini   PetscFunctionBegin;
15281bd0b88eSStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
1529633354d9SStefano Zampini   if (mapping->bs > 1) {
1530633354d9SStefano Zampini     if (n_procs) PetscCall(PetscFree(*n_procs));
1531633354d9SStefano Zampini     if (procs) {
1532633354d9SStefano Zampini       if (*procs) PetscCall(PetscFree((*procs)[0]));
1533633354d9SStefano Zampini       PetscCall(PetscFree(*procs));
1534633354d9SStefano Zampini     }
1535633354d9SStefano Zampini   }
1536633354d9SStefano Zampini   PetscCall(ISLocalToGlobalMappingRestoreBlockNodeInfo(mapping, n, n_procs, procs));
15373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15381bd0b88eSStefano Zampini }
15391bd0b88eSStefano Zampini 
15406ce40d10SBarry Smith /*MC
15416ce40d10SBarry Smith     ISLocalToGlobalMappingGetIndicesF90 - Get global indices for every local point that is mapped
15426ce40d10SBarry Smith 
15436ce40d10SBarry Smith     Synopsis:
15446ce40d10SBarry Smith     ISLocalToGlobalMappingGetIndicesF90(ISLocalToGlobalMapping ltog, PetscInt, pointer :: array(:)}, integer ierr)
15456ce40d10SBarry Smith 
15466ce40d10SBarry Smith     Not Collective
15476ce40d10SBarry Smith 
15486ce40d10SBarry Smith     Input Parameter:
15496ce40d10SBarry Smith .   A - the matrix
15506ce40d10SBarry Smith 
15512fe279fdSBarry Smith     Output Parameter:
15526ce40d10SBarry Smith .   array - array of indices, the length of this array may be obtained with `ISLocalToGlobalMappingGetSize()`
15536ce40d10SBarry Smith 
15546ce40d10SBarry Smith     Level: advanced
15556ce40d10SBarry Smith 
15566ce40d10SBarry Smith     Note:
15576ce40d10SBarry Smith     Use  `ISLocalToGlobalMappingGetIndicesF90()` when you no longer need access to the data
15586ce40d10SBarry Smith 
155920662ed9SBarry Smith .seealso: [](sec_fortranarrays), [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingGetIndices()`, `ISLocalToGlobalMappingRestoreIndices()`,
156020662ed9SBarry Smith           `ISLocalToGlobalMappingRestoreIndicesF90()`
15616ce40d10SBarry Smith M*/
15626ce40d10SBarry Smith 
15636ce40d10SBarry Smith /*MC
15646ce40d10SBarry Smith     ISLocalToGlobalMappingRestoreIndicesF90 - restores the global indices for every local point that is mapped obtained with `ISLocalToGlobalMappingGetIndicesF90()`
15656ce40d10SBarry Smith 
15666ce40d10SBarry Smith     Synopsis:
15676ce40d10SBarry Smith     ISLocalToGlobalMappingRestoreIndicesF90(ISLocalToGlobalMapping ltog, PetscInt, pointer :: array(:)}, integer ierr)
15686ce40d10SBarry Smith 
15696ce40d10SBarry Smith     Not Collective
15706ce40d10SBarry Smith 
15716ce40d10SBarry Smith     Input Parameters:
15726ce40d10SBarry Smith +   A - the matrix
15736ce40d10SBarry Smith -   array - array of indices, the length of this array may be obtained with `ISLocalToGlobalMappingGetSize()`
15746ce40d10SBarry Smith 
15756ce40d10SBarry Smith     Level: advanced
15766ce40d10SBarry Smith 
157720662ed9SBarry Smith .seealso: [](sec_fortranarrays), [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingGetIndices()`, `ISLocalToGlobalMappingRestoreIndices()`,
157820662ed9SBarry Smith           `ISLocalToGlobalMappingGetIndicesF90()`
15796ce40d10SBarry Smith M*/
15806ce40d10SBarry Smith 
15811bd0b88eSStefano Zampini /*@C
1582107e9a97SBarry Smith   ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped
158386994e45SJed Brown 
158486994e45SJed Brown   Not Collective
158586994e45SJed Brown 
15864165533cSJose E. Roman   Input Parameter:
158786994e45SJed Brown . ltog - local to global mapping
158886994e45SJed Brown 
15894165533cSJose E. Roman   Output Parameter:
1590cab54364SBarry Smith . array - array of indices, the length of this array may be obtained with `ISLocalToGlobalMappingGetSize()`
159186994e45SJed Brown 
159286994e45SJed Brown   Level: advanced
159386994e45SJed Brown 
1594cab54364SBarry Smith   Note:
1595cab54364SBarry Smith   `ISLocalToGlobalMappingGetSize()` returns the length the this array
1596107e9a97SBarry Smith 
159720662ed9SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingRestoreIndices()`,
159820662ed9SBarry Smith           `ISLocalToGlobalMappingGetBlockIndices()`, `ISLocalToGlobalMappingRestoreBlockIndices()`
159986994e45SJed Brown @*/
1600d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog, const PetscInt **array)
1601d71ae5a4SJacob Faibussowitsch {
160286994e45SJed Brown   PetscFunctionBegin;
160386994e45SJed Brown   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
16044f572ea9SToby Isaac   PetscAssertPointer(array, 2);
160545b6f7e9SBarry Smith   if (ltog->bs == 1) {
160686994e45SJed Brown     *array = ltog->indices;
160745b6f7e9SBarry Smith   } else {
160845b6f7e9SBarry Smith     PetscInt       *jj, k, i, j, n = ltog->n, bs = ltog->bs;
160945b6f7e9SBarry Smith     const PetscInt *ii;
161045b6f7e9SBarry Smith 
16119566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(bs * n, &jj));
161245b6f7e9SBarry Smith     *array = jj;
161345b6f7e9SBarry Smith     k      = 0;
161445b6f7e9SBarry Smith     ii     = ltog->indices;
161545b6f7e9SBarry Smith     for (i = 0; i < n; i++)
16169371c9d4SSatish Balay       for (j = 0; j < bs; j++) jj[k++] = bs * ii[i] + j;
161745b6f7e9SBarry Smith   }
16183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
161986994e45SJed Brown }
162086994e45SJed Brown 
162186994e45SJed Brown /*@C
1622cab54364SBarry Smith   ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with `ISLocalToGlobalMappingGetIndices()`
162386994e45SJed Brown 
162486994e45SJed Brown   Not Collective
162586994e45SJed Brown 
16264165533cSJose E. Roman   Input Parameters:
162786994e45SJed Brown + ltog  - local to global mapping
162886994e45SJed Brown - array - array of indices
162986994e45SJed Brown 
163086994e45SJed Brown   Level: advanced
163186994e45SJed Brown 
1632cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingGetIndices()`
163386994e45SJed Brown @*/
1634d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog, const PetscInt **array)
1635d71ae5a4SJacob Faibussowitsch {
163686994e45SJed Brown   PetscFunctionBegin;
163786994e45SJed Brown   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
16384f572ea9SToby Isaac   PetscAssertPointer(array, 2);
1639c9cc58a2SBarry Smith   PetscCheck(ltog->bs != 1 || *array == ltog->indices, PETSC_COMM_SELF, PETSC_ERR_ARG_BADPTR, "Trying to return mismatched pointer");
164048a46eb9SPierre Jolivet   if (ltog->bs > 1) PetscCall(PetscFree(*(void **)array));
16413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
164245b6f7e9SBarry Smith }
164345b6f7e9SBarry Smith 
164445b6f7e9SBarry Smith /*@C
164545b6f7e9SBarry Smith   ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block
164645b6f7e9SBarry Smith 
164745b6f7e9SBarry Smith   Not Collective
164845b6f7e9SBarry Smith 
16494165533cSJose E. Roman   Input Parameter:
165045b6f7e9SBarry Smith . ltog - local to global mapping
165145b6f7e9SBarry Smith 
16524165533cSJose E. Roman   Output Parameter:
165345b6f7e9SBarry Smith . array - array of indices
165445b6f7e9SBarry Smith 
165545b6f7e9SBarry Smith   Level: advanced
165645b6f7e9SBarry Smith 
1657cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingRestoreBlockIndices()`
165845b6f7e9SBarry Smith @*/
1659d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog, const PetscInt **array)
1660d71ae5a4SJacob Faibussowitsch {
166145b6f7e9SBarry Smith   PetscFunctionBegin;
166245b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
16634f572ea9SToby Isaac   PetscAssertPointer(array, 2);
166445b6f7e9SBarry Smith   *array = ltog->indices;
16653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
166645b6f7e9SBarry Smith }
166745b6f7e9SBarry Smith 
166845b6f7e9SBarry Smith /*@C
1669cab54364SBarry Smith   ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with `ISLocalToGlobalMappingGetBlockIndices()`
167045b6f7e9SBarry Smith 
167145b6f7e9SBarry Smith   Not Collective
167245b6f7e9SBarry Smith 
16734165533cSJose E. Roman   Input Parameters:
167445b6f7e9SBarry Smith + ltog  - local to global mapping
167545b6f7e9SBarry Smith - array - array of indices
167645b6f7e9SBarry Smith 
167745b6f7e9SBarry Smith   Level: advanced
167845b6f7e9SBarry Smith 
1679cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingGetIndices()`
168045b6f7e9SBarry Smith @*/
1681d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog, const PetscInt **array)
1682d71ae5a4SJacob Faibussowitsch {
168345b6f7e9SBarry Smith   PetscFunctionBegin;
168445b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
16854f572ea9SToby Isaac   PetscAssertPointer(array, 2);
168608401ef6SPierre Jolivet   PetscCheck(*array == ltog->indices, PETSC_COMM_SELF, PETSC_ERR_ARG_BADPTR, "Trying to return mismatched pointer");
16870298fd71SBarry Smith   *array = NULL;
16883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
168986994e45SJed Brown }
1690f7efa3c7SJed Brown 
1691f7efa3c7SJed Brown /*@C
1692f7efa3c7SJed Brown   ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings
1693f7efa3c7SJed Brown 
1694f7efa3c7SJed Brown   Not Collective
1695f7efa3c7SJed Brown 
16964165533cSJose E. Roman   Input Parameters:
1697f7efa3c7SJed Brown + comm  - communicator for the new mapping, must contain the communicator of every mapping to concatenate
1698f7efa3c7SJed Brown . n     - number of mappings to concatenate
1699f7efa3c7SJed Brown - ltogs - local to global mappings
1700f7efa3c7SJed Brown 
17014165533cSJose E. Roman   Output Parameter:
1702f7efa3c7SJed Brown . ltogcat - new mapping
1703f7efa3c7SJed Brown 
1704f7efa3c7SJed Brown   Level: advanced
1705f7efa3c7SJed Brown 
1706cab54364SBarry Smith   Note:
1707cab54364SBarry Smith   This currently always returns a mapping with block size of 1
1708cab54364SBarry Smith 
170938b5cf2dSJacob Faibussowitsch   Developer Notes:
1710cab54364SBarry Smith   If all the input mapping have the same block size we could easily handle that as a special case
1711cab54364SBarry Smith 
1712cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingCreate()`
1713f7efa3c7SJed Brown @*/
1714d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm, PetscInt n, const ISLocalToGlobalMapping ltogs[], ISLocalToGlobalMapping *ltogcat)
1715d71ae5a4SJacob Faibussowitsch {
1716f7efa3c7SJed Brown   PetscInt i, cnt, m, *idx;
1717f7efa3c7SJed Brown 
1718f7efa3c7SJed Brown   PetscFunctionBegin;
171908401ef6SPierre Jolivet   PetscCheck(n >= 0, comm, PETSC_ERR_ARG_OUTOFRANGE, "Must have a non-negative number of mappings, given %" PetscInt_FMT, n);
17204f572ea9SToby Isaac   if (n > 0) PetscAssertPointer(ltogs, 3);
1721f7efa3c7SJed Brown   for (i = 0; i < n; i++) PetscValidHeaderSpecific(ltogs[i], IS_LTOGM_CLASSID, 3);
17224f572ea9SToby Isaac   PetscAssertPointer(ltogcat, 4);
1723f7efa3c7SJed Brown   for (cnt = 0, i = 0; i < n; i++) {
17249566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(ltogs[i], &m));
1725f7efa3c7SJed Brown     cnt += m;
1726f7efa3c7SJed Brown   }
17279566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(cnt, &idx));
1728f7efa3c7SJed Brown   for (cnt = 0, i = 0; i < n; i++) {
1729f7efa3c7SJed Brown     const PetscInt *subidx;
17309566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(ltogs[i], &m));
17319566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetIndices(ltogs[i], &subidx));
17329566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(&idx[cnt], subidx, m));
17339566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreIndices(ltogs[i], &subidx));
1734f7efa3c7SJed Brown     cnt += m;
1735f7efa3c7SJed Brown   }
17369566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(comm, 1, cnt, idx, PETSC_OWN_POINTER, ltogcat));
17373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1738f7efa3c7SJed Brown }
173904a59952SBarry Smith 
1740413f72f0SBarry Smith /*MC
1741cab54364SBarry Smith       ISLOCALTOGLOBALMAPPINGBASIC - basic implementation of the `ISLocalToGlobalMapping` object. When `ISGlobalToLocalMappingApply()` is
1742413f72f0SBarry Smith                                     used this is good for only small and moderate size problems.
1743413f72f0SBarry Smith 
174420662ed9SBarry Smith    Options Database Key:
1745a2b725a8SWilliam Gropp .   -islocaltoglobalmapping_type basic - select this method
1746413f72f0SBarry Smith 
1747413f72f0SBarry Smith    Level: beginner
1748413f72f0SBarry Smith 
1749cab54364SBarry Smith    Developer Note:
1750cab54364SBarry Smith    This stores all the mapping information on each MPI rank.
1751cab54364SBarry Smith 
1752cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`, `ISLOCALTOGLOBALMAPPINGHASH`
1753413f72f0SBarry Smith M*/
1754d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Basic(ISLocalToGlobalMapping ltog)
1755d71ae5a4SJacob Faibussowitsch {
1756413f72f0SBarry Smith   PetscFunctionBegin;
1757413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Basic;
1758413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Basic;
1759413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Basic;
1760413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Basic;
17613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1762413f72f0SBarry Smith }
1763413f72f0SBarry Smith 
1764413f72f0SBarry Smith /*MC
1765cab54364SBarry Smith       ISLOCALTOGLOBALMAPPINGHASH - hash implementation of the `ISLocalToGlobalMapping` object. When `ISGlobalToLocalMappingApply()` is
1766ed56e8eaSBarry Smith                                     used this is good for large memory problems.
1767413f72f0SBarry Smith 
176820662ed9SBarry Smith    Options Database Key:
1769a2b725a8SWilliam Gropp .   -islocaltoglobalmapping_type hash - select this method
1770413f72f0SBarry Smith 
1771413f72f0SBarry Smith    Level: beginner
1772413f72f0SBarry Smith 
1773cab54364SBarry Smith    Note:
1774cab54364SBarry Smith     This is selected automatically for large problems if the user does not set the type.
1775cab54364SBarry Smith 
1776cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`, `ISLOCALTOGLOBALMAPPINGBASIC`
1777413f72f0SBarry Smith M*/
1778d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Hash(ISLocalToGlobalMapping ltog)
1779d71ae5a4SJacob Faibussowitsch {
1780413f72f0SBarry Smith   PetscFunctionBegin;
1781413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Hash;
1782413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Hash;
1783413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Hash;
1784413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Hash;
17853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1786413f72f0SBarry Smith }
1787413f72f0SBarry Smith 
1788413f72f0SBarry Smith /*@C
1789cab54364SBarry Smith   ISLocalToGlobalMappingRegister -  Registers a method for applying a global to local mapping with an `ISLocalToGlobalMapping`
1790413f72f0SBarry Smith 
1791413f72f0SBarry Smith   Not Collective
1792413f72f0SBarry Smith 
1793413f72f0SBarry Smith   Input Parameters:
1794413f72f0SBarry Smith + sname    - name of a new method
17952fe279fdSBarry Smith - function - routine to create method context
1796413f72f0SBarry Smith 
179738b5cf2dSJacob Faibussowitsch   Example Usage:
1798413f72f0SBarry Smith .vb
1799413f72f0SBarry Smith    ISLocalToGlobalMappingRegister("my_mapper", MyCreate);
1800413f72f0SBarry Smith .ve
1801413f72f0SBarry Smith 
1802ed56e8eaSBarry Smith   Then, your mapping can be chosen with the procedural interface via
1803413f72f0SBarry Smith $     ISLocalToGlobalMappingSetType(ltog, "my_mapper")
1804413f72f0SBarry Smith   or at runtime via the option
1805ed56e8eaSBarry Smith $     -islocaltoglobalmapping_type my_mapper
1806413f72f0SBarry Smith 
1807413f72f0SBarry Smith   Level: advanced
1808413f72f0SBarry Smith 
1809cab54364SBarry Smith   Note:
1810cab54364SBarry Smith   `ISLocalToGlobalMappingRegister()` may be called multiple times to add several user-defined mappings.
1811413f72f0SBarry Smith 
1812cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingRegisterAll()`, `ISLocalToGlobalMappingRegisterDestroy()`, `ISLOCALTOGLOBALMAPPINGBASIC`,
1813cab54364SBarry Smith           `ISLOCALTOGLOBALMAPPINGHASH`, `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApply()`
1814413f72f0SBarry Smith @*/
1815d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRegister(const char sname[], PetscErrorCode (*function)(ISLocalToGlobalMapping))
1816d71ae5a4SJacob Faibussowitsch {
1817413f72f0SBarry Smith   PetscFunctionBegin;
18189566063dSJacob Faibussowitsch   PetscCall(ISInitializePackage());
18199566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&ISLocalToGlobalMappingList, sname, function));
18203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1821413f72f0SBarry Smith }
1822413f72f0SBarry Smith 
1823413f72f0SBarry Smith /*@C
1824cab54364SBarry Smith   ISLocalToGlobalMappingSetType - Sets the implementation type `ISLocalToGlobalMapping` will use
1825413f72f0SBarry Smith 
1826c3339decSBarry Smith   Logically Collective
1827413f72f0SBarry Smith 
1828413f72f0SBarry Smith   Input Parameters:
1829cab54364SBarry Smith + ltog - the `ISLocalToGlobalMapping` object
1830413f72f0SBarry Smith - type - a known method
1831413f72f0SBarry Smith 
1832413f72f0SBarry Smith   Options Database Key:
1833cab54364SBarry Smith . -islocaltoglobalmapping_type  <method> - Sets the method; use -help for a list of available methods (for instance, basic or hash)
1834cab54364SBarry Smith 
1835cab54364SBarry Smith   Level: intermediate
1836413f72f0SBarry Smith 
1837413f72f0SBarry Smith   Notes:
183820662ed9SBarry Smith   See `ISLocalToGlobalMappingType` for available methods
1839413f72f0SBarry Smith 
1840cab54364SBarry Smith   Normally, it is best to use the `ISLocalToGlobalMappingSetFromOptions()` command and
1841cab54364SBarry Smith   then set the `ISLocalToGlobalMappingType` from the options database rather than by using
1842413f72f0SBarry Smith   this routine.
1843413f72f0SBarry Smith 
184438b5cf2dSJacob Faibussowitsch   Developer Notes:
1845cab54364SBarry Smith   `ISLocalToGlobalMappingRegister()` is used to add new types to `ISLocalToGlobalMappingList` from which they
1846cab54364SBarry Smith   are accessed by `ISLocalToGlobalMappingSetType()`.
1847413f72f0SBarry Smith 
184838b5cf2dSJacob Faibussowitsch .seealso: [](sec_scatter), `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingRegister()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingGetType()`
1849413f72f0SBarry Smith @*/
1850d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType type)
1851d71ae5a4SJacob Faibussowitsch {
1852413f72f0SBarry Smith   PetscBool match;
18535f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(ISLocalToGlobalMapping) = NULL;
1854413f72f0SBarry Smith 
1855413f72f0SBarry Smith   PetscFunctionBegin;
1856413f72f0SBarry Smith   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
18574f572ea9SToby Isaac   if (type) PetscAssertPointer(type, 2);
1858413f72f0SBarry Smith 
18599566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ltog, type, &match));
18603ba16761SJacob Faibussowitsch   if (match) PetscFunctionReturn(PETSC_SUCCESS);
1861413f72f0SBarry Smith 
1862a0d79125SStefano Zampini   /* L2G maps defer type setup at globaltolocal calls, allow passing NULL here */
1863a0d79125SStefano Zampini   if (type) {
18649566063dSJacob Faibussowitsch     PetscCall(PetscFunctionListFind(ISLocalToGlobalMappingList, type, &r));
1865a0d79125SStefano Zampini     PetscCheck(r, PetscObjectComm((PetscObject)ltog), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested ISLocalToGlobalMapping type %s", type);
1866a0d79125SStefano Zampini   }
1867413f72f0SBarry Smith   /* Destroy the previous private LTOG context */
1868dbbe0bcdSBarry Smith   PetscTryTypeMethod(ltog, destroy);
1869413f72f0SBarry Smith   ltog->ops->destroy = NULL;
1870dbbe0bcdSBarry Smith 
18719566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)ltog, type));
18729566063dSJacob Faibussowitsch   if (r) PetscCall((*r)(ltog));
18733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1874a0d79125SStefano Zampini }
1875a0d79125SStefano Zampini 
1876a0d79125SStefano Zampini /*@C
1877cab54364SBarry Smith   ISLocalToGlobalMappingGetType - Get the type of the `ISLocalToGlobalMapping`
1878a0d79125SStefano Zampini 
1879a0d79125SStefano Zampini   Not Collective
1880a0d79125SStefano Zampini 
1881a0d79125SStefano Zampini   Input Parameter:
1882cab54364SBarry Smith . ltog - the `ISLocalToGlobalMapping` object
1883a0d79125SStefano Zampini 
1884a0d79125SStefano Zampini   Output Parameter:
1885a0d79125SStefano Zampini . type - the type
1886a0d79125SStefano Zampini 
188749762cbcSSatish Balay   Level: intermediate
188849762cbcSSatish Balay 
188938b5cf2dSJacob Faibussowitsch .seealso: [](sec_scatter), `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingRegister()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`
1890a0d79125SStefano Zampini @*/
1891d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType *type)
1892d71ae5a4SJacob Faibussowitsch {
1893a0d79125SStefano Zampini   PetscFunctionBegin;
1894a0d79125SStefano Zampini   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
18954f572ea9SToby Isaac   PetscAssertPointer(type, 2);
1896a0d79125SStefano Zampini   *type = ((PetscObject)ltog)->type_name;
18973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1898413f72f0SBarry Smith }
1899413f72f0SBarry Smith 
1900413f72f0SBarry Smith PetscBool ISLocalToGlobalMappingRegisterAllCalled = PETSC_FALSE;
1901413f72f0SBarry Smith 
1902413f72f0SBarry Smith /*@C
1903cab54364SBarry Smith   ISLocalToGlobalMappingRegisterAll - Registers all of the local to global mapping components in the `IS` package.
1904413f72f0SBarry Smith 
1905413f72f0SBarry Smith   Not Collective
1906413f72f0SBarry Smith 
1907413f72f0SBarry Smith   Level: advanced
1908413f72f0SBarry Smith 
1909cab54364SBarry Smith .seealso: [](sec_scatter), `ISRegister()`, `ISLocalToGlobalRegister()`
1910413f72f0SBarry Smith @*/
1911d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRegisterAll(void)
1912d71ae5a4SJacob Faibussowitsch {
1913413f72f0SBarry Smith   PetscFunctionBegin;
19143ba16761SJacob Faibussowitsch   if (ISLocalToGlobalMappingRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
1915413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAllCalled = PETSC_TRUE;
19169566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGBASIC, ISLocalToGlobalMappingCreate_Basic));
19179566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGHASH, ISLocalToGlobalMappingCreate_Hash));
19183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1919413f72f0SBarry Smith }
1920