xref: /petsc/src/vec/is/utils/isltog.c (revision f1957bc3d7b66ed8a0ce45d1acff8920f4c07d7b)
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 @*/
48cc4c1da9SBarry Smith 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 @*/
79cc4c1da9SBarry Smith 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 
91cc4c1da9SBarry Smith /*@
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 @*/
112cc4c1da9SBarry Smith 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;
1391690c2aeSBarry Smith   start = PETSC_INT_MAX;
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 
2217de13914SStefano Zampini static PetscErrorCode ISLocalToGlobalMappingResetBlockInfo_Private(ISLocalToGlobalMapping mapping)
2227de13914SStefano Zampini {
2237de13914SStefano Zampini   PetscFunctionBegin;
2247de13914SStefano Zampini   PetscCall(PetscFree(mapping->info_procs));
2257de13914SStefano Zampini   PetscCall(PetscFree(mapping->info_numprocs));
2267de13914SStefano Zampini   if (mapping->info_indices) {
2277de13914SStefano Zampini     for (PetscInt i = 0; i < mapping->info_nproc; i++) PetscCall(PetscFree(mapping->info_indices[i]));
2287de13914SStefano Zampini     PetscCall(PetscFree(mapping->info_indices));
2297de13914SStefano Zampini   }
2307de13914SStefano Zampini   if (mapping->info_nodei) PetscCall(PetscFree(mapping->info_nodei[0]));
2317de13914SStefano Zampini   PetscCall(PetscFree2(mapping->info_nodec, mapping->info_nodei));
232d4df40f3SStefano Zampini   PetscCall(PetscSFDestroy(&mapping->multileaves_sf));
2337de13914SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
2347de13914SStefano Zampini }
2357de13914SStefano Zampini 
236413f72f0SBarry Smith #define GTOLTYPE _Basic
237413f72f0SBarry Smith #define GTOLNAME _Basic
238541bf97eSAdrian Croucher #define GTOLBS   mapping->bs
2399371c9d4SSatish Balay #define GTOL(g, local) \
2409371c9d4SSatish Balay   do { \
241413f72f0SBarry Smith     local = map->globals[g / bs - start]; \
2420040bde1SJunchao Zhang     if (local >= 0) local = bs * local + (g % bs); \
243413f72f0SBarry Smith   } while (0)
244541bf97eSAdrian Croucher 
245413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
246413f72f0SBarry Smith 
247413f72f0SBarry Smith #define GTOLTYPE _Basic
248413f72f0SBarry Smith #define GTOLNAME Block_Basic
249541bf97eSAdrian Croucher #define GTOLBS   1
2509371c9d4SSatish Balay #define GTOL(g, local) \
251d71ae5a4SJacob Faibussowitsch   do { \
252d71ae5a4SJacob Faibussowitsch     local = map->globals[g - start]; \
253d71ae5a4SJacob Faibussowitsch   } while (0)
254413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
255413f72f0SBarry Smith 
256413f72f0SBarry Smith #define GTOLTYPE _Hash
257413f72f0SBarry Smith #define GTOLNAME _Hash
258541bf97eSAdrian Croucher #define GTOLBS   mapping->bs
2599371c9d4SSatish Balay #define GTOL(g, local) \
2609371c9d4SSatish Balay   do { \
261e8f14785SLisandro Dalcin     (void)PetscHMapIGet(map->globalht, g / bs, &local); \
2620040bde1SJunchao Zhang     if (local >= 0) local = bs * local + (g % bs); \
263413f72f0SBarry Smith   } while (0)
264413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
265413f72f0SBarry Smith 
266413f72f0SBarry Smith #define GTOLTYPE _Hash
267413f72f0SBarry Smith #define GTOLNAME Block_Hash
268541bf97eSAdrian Croucher #define GTOLBS   1
2699371c9d4SSatish Balay #define GTOL(g, local) \
270d71ae5a4SJacob Faibussowitsch   do { \
271d71ae5a4SJacob Faibussowitsch     (void)PetscHMapIGet(map->globalht, g, &local); \
272d71ae5a4SJacob Faibussowitsch   } while (0)
273413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
274413f72f0SBarry Smith 
2756658fb44Sstefano_zampini /*@
2766658fb44Sstefano_zampini   ISLocalToGlobalMappingDuplicate - Duplicates the local to global mapping object
2776658fb44Sstefano_zampini 
2786658fb44Sstefano_zampini   Not Collective
2796658fb44Sstefano_zampini 
2806658fb44Sstefano_zampini   Input Parameter:
2816658fb44Sstefano_zampini . ltog - local to global mapping
2826658fb44Sstefano_zampini 
2836658fb44Sstefano_zampini   Output Parameter:
2846658fb44Sstefano_zampini . nltog - the duplicated local to global mapping
2856658fb44Sstefano_zampini 
2866658fb44Sstefano_zampini   Level: advanced
2876658fb44Sstefano_zampini 
288cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`
2896658fb44Sstefano_zampini @*/
290d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingDuplicate(ISLocalToGlobalMapping ltog, ISLocalToGlobalMapping *nltog)
291d71ae5a4SJacob Faibussowitsch {
292a0d79125SStefano Zampini   ISLocalToGlobalMappingType l2gtype;
2936658fb44Sstefano_zampini 
2946658fb44Sstefano_zampini   PetscFunctionBegin;
2956658fb44Sstefano_zampini   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
2969566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)ltog), ltog->bs, ltog->n, ltog->indices, PETSC_COPY_VALUES, nltog));
2979566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetType(ltog, &l2gtype));
2989566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingSetType(*nltog, l2gtype));
2993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3006658fb44Sstefano_zampini }
3016658fb44Sstefano_zampini 
302565245c5SBarry Smith /*@
303107e9a97SBarry Smith   ISLocalToGlobalMappingGetSize - Gets the local size of a local to global mapping
3043b9aefa3SBarry Smith 
3053b9aefa3SBarry Smith   Not Collective
3063b9aefa3SBarry Smith 
3073b9aefa3SBarry Smith   Input Parameter:
30838b5cf2dSJacob Faibussowitsch . mapping - local to global mapping
3093b9aefa3SBarry Smith 
3103b9aefa3SBarry Smith   Output Parameter:
311cab54364SBarry Smith . n - the number of entries in the local mapping, `ISLocalToGlobalMappingGetIndices()` returns an array of this length
3123b9aefa3SBarry Smith 
3133b9aefa3SBarry Smith   Level: advanced
3143b9aefa3SBarry Smith 
315cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`
3163b9aefa3SBarry Smith @*/
317d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping mapping, PetscInt *n)
318d71ae5a4SJacob Faibussowitsch {
3193b9aefa3SBarry Smith   PetscFunctionBegin;
3200700a824SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
3214f572ea9SToby Isaac   PetscAssertPointer(n, 2);
322107e9a97SBarry Smith   *n = mapping->bs * mapping->n;
3233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3243b9aefa3SBarry Smith }
3253b9aefa3SBarry Smith 
326ffeef943SBarry Smith /*@
327cab54364SBarry Smith   ISLocalToGlobalMappingViewFromOptions - View an `ISLocalToGlobalMapping` based on values in the options database
328fe2efc57SMark 
329c3339decSBarry Smith   Collective
330fe2efc57SMark 
331fe2efc57SMark   Input Parameters:
332fe2efc57SMark + A    - the local to global mapping object
33320662ed9SBarry Smith . obj  - Optional object that provides the options prefix used for the options database query
334736c3998SJose E. Roman - name - command line option
335fe2efc57SMark 
336fe2efc57SMark   Level: intermediate
337cab54364SBarry Smith 
33820662ed9SBarry Smith   Note:
33920662ed9SBarry Smith   See `PetscObjectViewFromOptions()` for the available `PetscViewer` and `PetscViewerFormat`
34020662ed9SBarry Smith 
341a94f484eSPierre Jolivet .seealso: [](sec_scatter), `PetscViewer`, `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingView`, `PetscObjectViewFromOptions()`, `ISLocalToGlobalMappingCreate()`
342fe2efc57SMark @*/
343d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingViewFromOptions(ISLocalToGlobalMapping A, PetscObject obj, const char name[])
344d71ae5a4SJacob Faibussowitsch {
345fe2efc57SMark   PetscFunctionBegin;
346fe2efc57SMark   PetscValidHeaderSpecific(A, IS_LTOGM_CLASSID, 1);
3479566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
3483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
349fe2efc57SMark }
350fe2efc57SMark 
351ffeef943SBarry Smith /*@
3525a5d4f66SBarry Smith   ISLocalToGlobalMappingView - View a local to global mapping
3535a5d4f66SBarry Smith 
3547de13914SStefano Zampini   Collective on viewer
355b9cd556bSLois Curfman McInnes 
3565a5d4f66SBarry Smith   Input Parameters:
35738b5cf2dSJacob Faibussowitsch + mapping - local to global mapping
3583b9aefa3SBarry Smith - viewer  - viewer
3595a5d4f66SBarry Smith 
3607de13914SStefano Zampini   Level: intermediate
361a997ad1aSLois Curfman McInnes 
36220662ed9SBarry Smith .seealso: [](sec_scatter), `PetscViewer`, `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`
3635a5d4f66SBarry Smith @*/
364d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingView(ISLocalToGlobalMapping mapping, PetscViewer viewer)
365d71ae5a4SJacob Faibussowitsch {
3669f196a02SMartin Diehl   PetscBool         isascii, isbinary;
3677de13914SStefano Zampini   PetscViewerFormat format;
3685a5d4f66SBarry Smith 
3695a5d4f66SBarry Smith   PetscFunctionBegin;
3700700a824SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
37148a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mapping), &viewer));
3720700a824SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
3735a5d4f66SBarry Smith 
3749f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
3757de13914SStefano Zampini   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
3767de13914SStefano Zampini   PetscCall(PetscViewerGetFormat(viewer, &format));
3779f196a02SMartin Diehl   if (isascii) {
3787de13914SStefano Zampini     if (format == PETSC_VIEWER_ASCII_MATLAB) {
3797de13914SStefano Zampini       const PetscInt *idxs;
3807de13914SStefano Zampini       IS              is;
3817de13914SStefano Zampini       const char     *name = ((PetscObject)mapping)->name;
3827de13914SStefano Zampini       char            iname[PETSC_MAX_PATH_LEN];
3837de13914SStefano Zampini 
3847de13914SStefano Zampini       PetscCall(PetscSNPrintf(iname, sizeof(iname), "%sl2g", name ? name : ""));
3857de13914SStefano Zampini       PetscCall(ISLocalToGlobalMappingGetIndices(mapping, &idxs));
3867de13914SStefano Zampini       PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)viewer), mapping->n * mapping->bs, idxs, PETSC_USE_POINTER, &is));
3877de13914SStefano Zampini       PetscCall(PetscObjectSetName((PetscObject)is, iname));
3887de13914SStefano Zampini       PetscCall(ISView(is, viewer));
3897de13914SStefano Zampini       PetscCall(ISLocalToGlobalMappingRestoreIndices(mapping, &idxs));
3907de13914SStefano Zampini       PetscCall(ISDestroy(&is));
3917de13914SStefano Zampini     } else {
3927de13914SStefano Zampini       PetscMPIInt rank;
3937de13914SStefano Zampini 
3947de13914SStefano Zampini       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mapping), &rank));
3959566063dSJacob Faibussowitsch       PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)mapping, viewer));
3969566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
3977de13914SStefano Zampini       for (PetscInt i = 0; i < mapping->n; i++) {
398f2c6b1a2SJed Brown         PetscInt bs = mapping->bs, g = mapping->indices[i];
399f2c6b1a2SJed Brown         if (bs == 1) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %" PetscInt_FMT " %" PetscInt_FMT "\n", rank, i, g));
400f2c6b1a2SJed 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));
401f2c6b1a2SJed Brown       }
4029566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
4039566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
4041575c14dSBarry Smith     }
4057de13914SStefano Zampini   } else if (isbinary) {
4067de13914SStefano Zampini     PetscBool skipHeader;
4077de13914SStefano Zampini 
4087de13914SStefano Zampini     PetscCall(PetscViewerSetUp(viewer));
4097de13914SStefano Zampini     PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &skipHeader));
4107de13914SStefano Zampini     if (!skipHeader) {
4117de13914SStefano Zampini       PetscMPIInt size;
4127de13914SStefano Zampini       PetscInt    tr[3];
4137de13914SStefano Zampini 
4147de13914SStefano Zampini       PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)viewer), &size));
4157de13914SStefano Zampini       tr[0] = IS_LTOGM_FILE_CLASSID;
4167de13914SStefano Zampini       tr[1] = mapping->bs;
4177de13914SStefano Zampini       tr[2] = size;
4187de13914SStefano Zampini       PetscCall(PetscViewerBinaryWrite(viewer, tr, 3, PETSC_INT));
4197de13914SStefano Zampini       PetscCall(PetscViewerBinaryWriteAll(viewer, &mapping->n, 1, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_INT));
4207de13914SStefano Zampini     }
4217de13914SStefano Zampini     /* write block indices */
4227de13914SStefano Zampini     PetscCall(PetscViewerBinaryWriteAll(viewer, mapping->indices, mapping->n, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_INT));
4237de13914SStefano Zampini   }
4247de13914SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
4257de13914SStefano Zampini }
4267de13914SStefano Zampini 
4277de13914SStefano Zampini /*@
4287de13914SStefano Zampini   ISLocalToGlobalMappingLoad - Loads a local-to-global mapping that has been stored in binary format.
4297de13914SStefano Zampini 
4307de13914SStefano Zampini   Collective on viewer
4317de13914SStefano Zampini 
4327de13914SStefano Zampini   Input Parameters:
4337de13914SStefano Zampini + mapping - the newly loaded map, this needs to have been created with `ISLocalToGlobalMappingCreate()` or some related function before a call to `ISLocalToGlobalMappingLoad()`
4347de13914SStefano Zampini - viewer  - binary file viewer, obtained from `PetscViewerBinaryOpen()`
4357de13914SStefano Zampini 
4367de13914SStefano Zampini   Level: intermediate
4377de13914SStefano Zampini 
4387de13914SStefano Zampini .seealso: [](sec_scatter), `PetscViewer`, `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingView()`, `ISLocalToGlobalMappingCreate()`
4397de13914SStefano Zampini @*/
4407de13914SStefano Zampini PetscErrorCode ISLocalToGlobalMappingLoad(ISLocalToGlobalMapping mapping, PetscViewer viewer)
4417de13914SStefano Zampini {
4427de13914SStefano Zampini   PetscBool isbinary, skipHeader;
4437de13914SStefano Zampini 
4447de13914SStefano Zampini   PetscFunctionBegin;
4457de13914SStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
4467de13914SStefano Zampini   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
4477de13914SStefano Zampini   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
4487de13914SStefano Zampini   PetscCheck(isbinary, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Invalid viewer of type %s", ((PetscObject)viewer)->type_name);
4497de13914SStefano Zampini 
4507de13914SStefano Zampini   /* reset previous data */
4517de13914SStefano Zampini   PetscCall(ISLocalToGlobalMappingResetBlockInfo_Private(mapping));
4527de13914SStefano Zampini 
4537de13914SStefano Zampini   PetscCall(PetscViewerSetUp(viewer));
4547de13914SStefano Zampini   PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &skipHeader));
4557de13914SStefano Zampini 
4567de13914SStefano Zampini   /* When skipping header, it assumes bs and n have been already set */
4577de13914SStefano Zampini   if (!skipHeader) {
4587de13914SStefano Zampini     MPI_Comm comm = PetscObjectComm((PetscObject)viewer);
4597de13914SStefano Zampini     PetscInt tr[3], nold = mapping->n, *sizes, nmaps = PETSC_DECIDE, st = 0;
4607de13914SStefano Zampini 
4617de13914SStefano Zampini     PetscCall(PetscViewerBinaryRead(viewer, tr, 3, NULL, PETSC_INT));
4627de13914SStefano Zampini     PetscCheck(tr[0] == IS_LTOGM_FILE_CLASSID, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Not a local-to-global map next in file");
4637de13914SStefano Zampini 
4647de13914SStefano Zampini     mapping->bs = tr[1];
4657de13914SStefano Zampini     PetscCall(PetscMalloc1(tr[2], &sizes));
4667de13914SStefano Zampini     PetscCall(PetscViewerBinaryRead(viewer, sizes, tr[2], NULL, PETSC_INT));
4677de13914SStefano Zampini 
4687de13914SStefano Zampini     /* consume the input, read multiple maps per process if needed */
4697de13914SStefano Zampini     PetscCall(PetscSplitOwnership(comm, &nmaps, &tr[2]));
4707de13914SStefano Zampini     PetscCallMPI(MPI_Exscan(&nmaps, &st, 1, MPIU_INT, MPI_SUM, comm));
4717de13914SStefano Zampini     mapping->n = 0;
4727de13914SStefano Zampini     for (PetscInt i = st; i < st + nmaps; i++) mapping->n += sizes[i];
4737de13914SStefano Zampini     PetscCall(PetscFree(sizes));
4747de13914SStefano Zampini 
4757de13914SStefano Zampini     if (nold != mapping->n) {
4767de13914SStefano Zampini       if (mapping->dealloc_indices) PetscCall(PetscFree(mapping->indices));
4777de13914SStefano Zampini       mapping->indices = NULL;
4787de13914SStefano Zampini     }
4797de13914SStefano Zampini   }
4807de13914SStefano Zampini 
4817de13914SStefano Zampini   /* read indices */
4827de13914SStefano Zampini   if (mapping->n && !mapping->indices) {
4837de13914SStefano Zampini     PetscCall(PetscMalloc1(mapping->n, &mapping->indices));
4847de13914SStefano Zampini     mapping->dealloc_indices = PETSC_TRUE;
4857de13914SStefano Zampini   }
4867de13914SStefano Zampini   PetscCall(PetscViewerBinaryReadAll(viewer, mapping->indices, mapping->n, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_INT));
4873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4885a5d4f66SBarry Smith }
4895a5d4f66SBarry Smith 
4901f428162SBarry Smith /*@
4912bdab257SBarry Smith   ISLocalToGlobalMappingCreateIS - Creates a mapping between a local (0 to n)
4922bdab257SBarry Smith   ordering and a global parallel ordering.
4932bdab257SBarry Smith 
49420662ed9SBarry Smith   Not Collective
495b9cd556bSLois Curfman McInnes 
496a997ad1aSLois Curfman McInnes   Input Parameter:
4978c03b21aSDmitry Karpeev . is - index set containing the global numbers for each local number
4982bdab257SBarry Smith 
499a997ad1aSLois Curfman McInnes   Output Parameter:
5002bdab257SBarry Smith . mapping - new mapping data structure
5012bdab257SBarry Smith 
502a997ad1aSLois Curfman McInnes   Level: advanced
503a997ad1aSLois Curfman McInnes 
504cab54364SBarry Smith   Note:
505cab54364SBarry Smith   the block size of the `IS` determines the block size of the mapping
506cab54364SBarry Smith 
507cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetFromOptions()`
5082bdab257SBarry Smith @*/
509d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingCreateIS(IS is, ISLocalToGlobalMapping *mapping)
510d71ae5a4SJacob Faibussowitsch {
5113bbf0e92SBarry Smith   PetscInt        n, bs;
5125d0c19d7SBarry Smith   const PetscInt *indices;
5132bdab257SBarry Smith   MPI_Comm        comm;
5143bbf0e92SBarry Smith   PetscBool       isblock;
5153a40ed3dSBarry Smith 
5163a40ed3dSBarry Smith   PetscFunctionBegin;
5170700a824SBarry Smith   PetscValidHeaderSpecific(is, IS_CLASSID, 1);
5184f572ea9SToby Isaac   PetscAssertPointer(mapping, 2);
5192bdab257SBarry Smith 
5209566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)is, &comm));
5219566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is, &n));
5229566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)is, ISBLOCK, &isblock));
5236006e8d2SBarry Smith   if (!isblock) {
5249566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(is, &indices));
5259566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(comm, 1, n, indices, PETSC_COPY_VALUES, mapping));
5269566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(is, &indices));
5276006e8d2SBarry Smith   } else {
5289566063dSJacob Faibussowitsch     PetscCall(ISGetBlockSize(is, &bs));
5299566063dSJacob Faibussowitsch     PetscCall(ISBlockGetIndices(is, &indices));
5309566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(comm, bs, n / bs, indices, PETSC_COPY_VALUES, mapping));
5319566063dSJacob Faibussowitsch     PetscCall(ISBlockRestoreIndices(is, &indices));
5326006e8d2SBarry Smith   }
5333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5342bdab257SBarry Smith }
5355a5d4f66SBarry Smith 
536cc4c1da9SBarry Smith /*@
537d4df40f3SStefano Zampini   ISLocalToGlobalMappingCreateSF - Creates a mapping between a local (0 to n) ordering and a global parallel ordering induced by a star forest.
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);
75857508eceSPierre Jolivet   PetscCall(PetscOptionsFList("-islocaltoglobalmapping_type", "ISLocalToGlobalMapping method", "ISLocalToGlobalMappingSetType", ISLocalToGlobalMappingList, ((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:
810586b16e1SPierre Jolivet   The output `IS` will have the same communicator as the input `IS` as well as the same block size.
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 {
817586b16e1SPierre Jolivet   PetscInt        n, *idxout, bs;
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));
826586b16e1SPierre Jolivet   PetscCall(ISGetBlockSize(is, &bs));
8279566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is, &idxin));
8289566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &idxout));
8299566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(mapping, n, idxin, idxout));
8309566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is, &idxin));
8319566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)is), n, idxout, PETSC_OWN_POINTER, newis));
832586b16e1SPierre Jolivet   PetscCall(ISSetBlockSize(*newis, bs));
8333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
83490f02eecSBarry Smith }
83590f02eecSBarry Smith 
836b89cb25eSSatish Balay /*@
8373acfe500SLois Curfman McInnes   ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering
8383acfe500SLois Curfman McInnes   and converts them to the global numbering.
83990f02eecSBarry Smith 
84020662ed9SBarry Smith   Not Collective
841b9cd556bSLois Curfman McInnes 
842bb25748dSBarry Smith   Input Parameters:
843b9cd556bSLois Curfman McInnes + mapping - the local to global mapping context
844bb25748dSBarry Smith . N       - number of integers
845b9cd556bSLois Curfman McInnes - in      - input indices in local numbering
846bb25748dSBarry Smith 
847bb25748dSBarry Smith   Output Parameter:
848bb25748dSBarry Smith . out - indices in global numbering
849bb25748dSBarry Smith 
850a997ad1aSLois Curfman McInnes   Level: advanced
851a997ad1aSLois Curfman McInnes 
852cab54364SBarry Smith   Note:
85320662ed9SBarry Smith   The `in` and `out` array parameters may be identical.
854cab54364SBarry Smith 
855cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApplyBlock()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingDestroy()`,
856c2e3fba1SPatrick Sanan           `ISLocalToGlobalMappingApplyIS()`, `AOCreateBasic()`, `AOApplicationToPetsc()`,
857db781477SPatrick Sanan           `AOPetscToApplication()`, `ISGlobalToLocalMappingApply()`
858afcb2eb5SJed Brown @*/
859d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping, PetscInt N, const PetscInt in[], PetscInt out[])
860d71ae5a4SJacob Faibussowitsch {
861cbc1caf0SMatthew G. Knepley   PetscInt i, bs, Nmax;
86245b6f7e9SBarry Smith 
86345b6f7e9SBarry Smith   PetscFunctionBegin;
864cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
865cbc1caf0SMatthew G. Knepley   bs   = mapping->bs;
866cbc1caf0SMatthew G. Knepley   Nmax = bs * mapping->n;
86745b6f7e9SBarry Smith   if (bs == 1) {
868cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
86945b6f7e9SBarry Smith     for (i = 0; i < N; i++) {
87045b6f7e9SBarry Smith       if (in[i] < 0) {
87145b6f7e9SBarry Smith         out[i] = in[i];
87245b6f7e9SBarry Smith         continue;
87345b6f7e9SBarry Smith       }
87408401ef6SPierre 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);
87545b6f7e9SBarry Smith       out[i] = idx[in[i]];
87645b6f7e9SBarry Smith     }
87745b6f7e9SBarry Smith   } else {
878cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
87945b6f7e9SBarry Smith     for (i = 0; i < N; i++) {
88045b6f7e9SBarry Smith       if (in[i] < 0) {
88145b6f7e9SBarry Smith         out[i] = in[i];
88245b6f7e9SBarry Smith         continue;
88345b6f7e9SBarry Smith       }
88408401ef6SPierre 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);
88545b6f7e9SBarry Smith       out[i] = idx[in[i] / bs] * bs + (in[i] % bs);
88645b6f7e9SBarry Smith     }
88745b6f7e9SBarry Smith   }
8883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
88945b6f7e9SBarry Smith }
89045b6f7e9SBarry Smith 
89145b6f7e9SBarry Smith /*@
8926006e8d2SBarry Smith   ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local block numbering and converts them to the global block numbering
89345b6f7e9SBarry Smith 
89420662ed9SBarry Smith   Not Collective
89545b6f7e9SBarry Smith 
89645b6f7e9SBarry Smith   Input Parameters:
89745b6f7e9SBarry Smith + mapping - the local to global mapping context
89845b6f7e9SBarry Smith . N       - number of integers
8996006e8d2SBarry Smith - in      - input indices in local block numbering
90045b6f7e9SBarry Smith 
90145b6f7e9SBarry Smith   Output Parameter:
9026006e8d2SBarry Smith . out - indices in global block numbering
90345b6f7e9SBarry Smith 
9046006e8d2SBarry Smith   Example:
905cab54364SBarry 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
9066006e8d2SBarry Smith   (the first block) would produce 0 and the mapping applied to 1 (the second block) would produce 3.
9076006e8d2SBarry Smith 
90820662ed9SBarry Smith   Level: advanced
90920662ed9SBarry Smith 
91020662ed9SBarry Smith   Note:
91120662ed9SBarry Smith   The `in` and `out` array parameters may be identical.
91220662ed9SBarry Smith 
913cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingDestroy()`,
914c2e3fba1SPatrick Sanan           `ISLocalToGlobalMappingApplyIS()`, `AOCreateBasic()`, `AOApplicationToPetsc()`,
915db781477SPatrick Sanan           `AOPetscToApplication()`, `ISGlobalToLocalMappingApply()`
91645b6f7e9SBarry Smith @*/
917d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping, PetscInt N, const PetscInt in[], PetscInt out[])
918d71ae5a4SJacob Faibussowitsch {
9198a1f772fSStefano Zampini   PetscInt        i, Nmax;
9208a1f772fSStefano Zampini   const PetscInt *idx;
921d4bb536fSBarry Smith 
922a0d79125SStefano Zampini   PetscFunctionBegin;
923a0d79125SStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
9248a1f772fSStefano Zampini   Nmax = mapping->n;
9258a1f772fSStefano Zampini   idx  = mapping->indices;
926afcb2eb5SJed Brown   for (i = 0; i < N; i++) {
927afcb2eb5SJed Brown     if (in[i] < 0) {
928afcb2eb5SJed Brown       out[i] = in[i];
929afcb2eb5SJed Brown       continue;
930afcb2eb5SJed Brown     }
93108401ef6SPierre 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);
932afcb2eb5SJed Brown     out[i] = idx[in[i]];
933afcb2eb5SJed Brown   }
9343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
935afcb2eb5SJed Brown }
936d4bb536fSBarry Smith 
9377e99dc12SLawrence Mitchell /*@
938a997ad1aSLois Curfman McInnes   ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers
939a997ad1aSLois Curfman McInnes   specified with a global numbering.
940d4bb536fSBarry Smith 
94120662ed9SBarry Smith   Not Collective
942b9cd556bSLois Curfman McInnes 
943d4bb536fSBarry Smith   Input Parameters:
944b9cd556bSLois Curfman McInnes + mapping - mapping between local and global numbering
945cab54364SBarry Smith . type    - `IS_GTOLM_MASK` - maps global indices with no local value to -1 in the output list (i.e., mask them)
946cab54364SBarry Smith            `IS_GTOLM_DROP` - drops the indices with no local value from the output list
947d4bb536fSBarry Smith . n       - number of global indices to map
948b9cd556bSLois Curfman McInnes - idx     - global indices to map
949d4bb536fSBarry Smith 
950d4bb536fSBarry Smith   Output Parameters:
951cab54364SBarry Smith + nout   - number of indices in output array (if type == `IS_GTOLM_MASK` then nout = n)
952b9cd556bSLois Curfman McInnes - idxout - local index of each global index, one must pass in an array long enough
953cab54364SBarry Smith              to hold all the indices. You can call `ISGlobalToLocalMappingApply()` with
9540298fd71SBarry Smith              idxout == NULL to determine the required length (returned in nout)
955cab54364SBarry Smith              and then allocate the required space and call `ISGlobalToLocalMappingApply()`
956e182c471SBarry Smith              a second time to set the values.
957d4bb536fSBarry Smith 
958cab54364SBarry Smith   Level: advanced
959cab54364SBarry Smith 
960b9cd556bSLois Curfman McInnes   Notes:
96120662ed9SBarry Smith   Either `nout` or `idxout` may be `NULL`. `idx` and `idxout` may be identical.
962d4bb536fSBarry Smith 
96320662ed9SBarry Smith   For "small" problems when using `ISGlobalToLocalMappingApply()` and `ISGlobalToLocalMappingApplyBlock()`, the `ISLocalToGlobalMappingType` of
96420662ed9SBarry Smith   `ISLOCALTOGLOBALMAPPINGBASIC` will be used;
965cab54364SBarry 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.
966cab54364SBarry Smith   Use `ISLocalToGlobalMappingSetType()` or call `ISLocalToGlobalMappingSetFromOptions()` with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
9670f5bd95cSBarry Smith 
96838b5cf2dSJacob Faibussowitsch   Developer Notes:
96920662ed9SBarry Smith   The manual page states that `idx` and `idxout` may be identical but the calling
97020662ed9SBarry Smith   sequence declares `idx` as const so it cannot be the same as `idxout`.
97132fd6b96SBarry Smith 
972cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApply()`, `ISGlobalToLocalMappingApplyBlock()`, `ISLocalToGlobalMappingCreate()`,
973db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`
974d4bb536fSBarry Smith @*/
975d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, PetscInt n, const PetscInt idx[], PetscInt *nout, PetscInt idxout[])
976d71ae5a4SJacob Faibussowitsch {
9779d90f715SBarry Smith   PetscFunctionBegin;
9789d90f715SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
97948a46eb9SPierre Jolivet   if (!mapping->data) PetscCall(ISGlobalToLocalMappingSetUp(mapping));
980dbbe0bcdSBarry Smith   PetscUseTypeMethod(mapping, globaltolocalmappingapply, type, n, idx, nout, idxout);
9813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9829d90f715SBarry Smith }
9839d90f715SBarry Smith 
984d4fe737eSStefano Zampini /*@
985cab54364SBarry Smith   ISGlobalToLocalMappingApplyIS - Creates from an `IS` in the global numbering
986cab54364SBarry Smith   a new index set using the local numbering defined in an `ISLocalToGlobalMapping`
987d4fe737eSStefano Zampini   context.
988d4fe737eSStefano Zampini 
98920662ed9SBarry Smith   Not Collective
990d4fe737eSStefano Zampini 
991d4fe737eSStefano Zampini   Input Parameters:
992d4fe737eSStefano Zampini + mapping - mapping between local and global numbering
993cab54364SBarry Smith . type    - `IS_GTOLM_MASK` - maps global indices with no local value to -1 in the output list (i.e., mask them)
994cab54364SBarry Smith            `IS_GTOLM_DROP` - drops the indices with no local value from the output list
995d4fe737eSStefano Zampini - is      - index set in global numbering
996d4fe737eSStefano Zampini 
9972fe279fdSBarry Smith   Output Parameter:
998d4fe737eSStefano Zampini . newis - index set in local numbering
999d4fe737eSStefano Zampini 
1000d4fe737eSStefano Zampini   Level: advanced
1001d4fe737eSStefano Zampini 
1002586b16e1SPierre Jolivet   Notes:
1003cab54364SBarry Smith   The output `IS` will be sequential, as it encodes a purely local operation
1004cab54364SBarry Smith 
1005586b16e1SPierre Jolivet   If `type` is `IS_GTOLM_MASK`, `newis` will have the same block size as `is`
1006586b16e1SPierre Jolivet 
1007cab54364SBarry Smith .seealso: [](sec_scatter), `ISGlobalToLocalMapping`, `ISGlobalToLocalMappingApply()`, `ISLocalToGlobalMappingCreate()`,
1008db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`
1009d4fe737eSStefano Zampini @*/
1010d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, IS is, IS *newis)
1011d71ae5a4SJacob Faibussowitsch {
1012586b16e1SPierre Jolivet   PetscInt        n, nout, *idxout, bs;
1013d4fe737eSStefano Zampini   const PetscInt *idxin;
1014d4fe737eSStefano Zampini 
1015d4fe737eSStefano Zampini   PetscFunctionBegin;
1016d4fe737eSStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
1017d4fe737eSStefano Zampini   PetscValidHeaderSpecific(is, IS_CLASSID, 3);
10184f572ea9SToby Isaac   PetscAssertPointer(newis, 4);
1019d4fe737eSStefano Zampini 
10209566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is, &n));
10219566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is, &idxin));
1022d4fe737eSStefano Zampini   if (type == IS_GTOLM_MASK) {
10239566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &idxout));
1024d4fe737eSStefano Zampini   } else {
10259566063dSJacob Faibussowitsch     PetscCall(ISGlobalToLocalMappingApply(mapping, type, n, idxin, &nout, NULL));
10269566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nout, &idxout));
1027d4fe737eSStefano Zampini   }
10289566063dSJacob Faibussowitsch   PetscCall(ISGlobalToLocalMappingApply(mapping, type, n, idxin, &nout, idxout));
10299566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is, &idxin));
10309566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nout, idxout, PETSC_OWN_POINTER, newis));
1031586b16e1SPierre Jolivet   if (type == IS_GTOLM_MASK) {
1032586b16e1SPierre Jolivet     PetscCall(ISGetBlockSize(is, &bs));
1033586b16e1SPierre Jolivet     PetscCall(ISSetBlockSize(*newis, bs));
1034586b16e1SPierre Jolivet   }
10353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1036d4fe737eSStefano Zampini }
1037d4fe737eSStefano Zampini 
10389d90f715SBarry Smith /*@
10399d90f715SBarry Smith   ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers
10409d90f715SBarry Smith   specified with a block global numbering.
10419d90f715SBarry Smith 
104220662ed9SBarry Smith   Not Collective
10439d90f715SBarry Smith 
10449d90f715SBarry Smith   Input Parameters:
10459d90f715SBarry Smith + mapping - mapping between local and global numbering
1046cab54364SBarry Smith . type    - `IS_GTOLM_MASK` - maps global indices with no local value to -1 in the output list (i.e., mask them)
1047cab54364SBarry Smith            `IS_GTOLM_DROP` - drops the indices with no local value from the output list
10489d90f715SBarry Smith . n       - number of global indices to map
10499d90f715SBarry Smith - idx     - global indices to map
10509d90f715SBarry Smith 
10519d90f715SBarry Smith   Output Parameters:
1052cab54364SBarry Smith + nout   - number of indices in output array (if type == `IS_GTOLM_MASK` then nout = n)
10539d90f715SBarry Smith - idxout - local index of each global index, one must pass in an array long enough
1054cab54364SBarry Smith              to hold all the indices. You can call `ISGlobalToLocalMappingApplyBlock()` with
10559d90f715SBarry Smith              idxout == NULL to determine the required length (returned in nout)
1056cab54364SBarry Smith              and then allocate the required space and call `ISGlobalToLocalMappingApplyBlock()`
10579d90f715SBarry Smith              a second time to set the values.
10589d90f715SBarry Smith 
1059cab54364SBarry Smith   Level: advanced
1060cab54364SBarry Smith 
10619d90f715SBarry Smith   Notes:
106220662ed9SBarry Smith   Either `nout` or `idxout` may be `NULL`. `idx` and `idxout` may be identical.
10639d90f715SBarry Smith 
106420662ed9SBarry Smith   For "small" problems when using `ISGlobalToLocalMappingApply()` and `ISGlobalToLocalMappingApplyBlock()`, the `ISLocalToGlobalMappingType` of
106520662ed9SBarry Smith   `ISLOCALTOGLOBALMAPPINGBASIC` will be used;
1066cab54364SBarry 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.
1067cab54364SBarry Smith   Use `ISLocalToGlobalMappingSetType()` or call `ISLocalToGlobalMappingSetFromOptions()` with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
10689d90f715SBarry Smith 
106938b5cf2dSJacob Faibussowitsch   Developer Notes:
107020662ed9SBarry Smith   The manual page states that `idx` and `idxout` may be identical but the calling
107120662ed9SBarry Smith   sequence declares `idx` as const so it cannot be the same as `idxout`.
10729d90f715SBarry Smith 
1073cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApply()`, `ISGlobalToLocalMappingApply()`, `ISLocalToGlobalMappingCreate()`,
1074db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`
10759d90f715SBarry Smith @*/
1076d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, PetscInt n, const PetscInt idx[], PetscInt *nout, PetscInt idxout[])
1077d71ae5a4SJacob Faibussowitsch {
10783a40ed3dSBarry Smith   PetscFunctionBegin;
10790700a824SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
108048a46eb9SPierre Jolivet   if (!mapping->data) PetscCall(ISGlobalToLocalMappingSetUp(mapping));
1081dbbe0bcdSBarry Smith   PetscUseTypeMethod(mapping, globaltolocalmappingapplyblock, type, n, idx, nout, idxout);
10823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1083d4bb536fSBarry Smith }
108490f02eecSBarry Smith 
108589d82c54SBarry Smith /*@C
1086633354d9SStefano Zampini   ISLocalToGlobalMappingGetBlockInfo - Gets the neighbor information
108789d82c54SBarry Smith 
1088633354d9SStefano Zampini   Collective the first time it is called
108989d82c54SBarry Smith 
1090f899ff85SJose E. Roman   Input Parameter:
109189d82c54SBarry Smith . mapping - the mapping from local to global indexing
109289d82c54SBarry Smith 
1093d8d19677SJose E. Roman   Output Parameters:
1094633354d9SStefano Zampini + nproc    - number of processes that are connected to the calling process
1095633354d9SStefano Zampini . procs    - neighboring processes
1096633354d9SStefano Zampini . numprocs - number of block indices for each process
1097633354d9SStefano Zampini - indices  - block indices (in local numbering) shared with neighbors (sorted by global numbering)
109889d82c54SBarry Smith 
109989d82c54SBarry Smith   Level: advanced
110089d82c54SBarry Smith 
1101cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1102d4df40f3SStefano Zampini           `ISLocalToGlobalMappingRestoreBlockInfo()`, `ISLocalToGlobalMappingGetBlockMultiLeavesSF()`
110389d82c54SBarry Smith @*/
1104d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[])
1105d71ae5a4SJacob Faibussowitsch {
1106268a049cSStefano Zampini   PetscFunctionBegin;
1107268a049cSStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
1108633354d9SStefano Zampini   PetscCall(ISLocalToGlobalMappingSetUpBlockInfo_Private(mapping));
1109633354d9SStefano Zampini   if (nproc) *nproc = mapping->info_nproc;
1110633354d9SStefano Zampini   if (procs) *procs = mapping->info_procs;
1111633354d9SStefano Zampini   if (numprocs) *numprocs = mapping->info_numprocs;
1112633354d9SStefano Zampini   if (indices) *indices = mapping->info_indices;
11133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1114268a049cSStefano Zampini }
1115268a049cSStefano Zampini 
1116633354d9SStefano Zampini /*@C
1117633354d9SStefano Zampini   ISLocalToGlobalMappingGetBlockNodeInfo - Gets the neighbor information for each local block index
1118633354d9SStefano Zampini 
1119633354d9SStefano Zampini   Collective the first time it is called
1120633354d9SStefano Zampini 
1121633354d9SStefano Zampini   Input Parameter:
1122633354d9SStefano Zampini . mapping - the mapping from local to global indexing
1123633354d9SStefano Zampini 
11249cde84edSJose E. Roman   Output Parameters:
1125633354d9SStefano Zampini + n       - number of local block nodes
1126633354d9SStefano Zampini . n_procs - an array storing the number of processes for each local block node (including self)
1127633354d9SStefano Zampini - procs   - the processes' rank for each local block node (sorted, self is first)
1128633354d9SStefano Zampini 
1129633354d9SStefano Zampini   Level: advanced
1130633354d9SStefano Zampini 
1131633354d9SStefano Zampini   Notes:
1132633354d9SStefano Zampini   The user needs to call `ISLocalToGlobalMappingRestoreBlockNodeInfo()` when the data is no longer needed.
1133633354d9SStefano Zampini   The information returned by this function complements that of `ISLocalToGlobalMappingGetBlockInfo()`.
1134633354d9SStefano Zampini   The latter only provides local information, and the neighboring information
1135633354d9SStefano Zampini   cannot be inferred in the general case, unless the mapping is locally one-to-one on each process.
1136633354d9SStefano Zampini 
1137633354d9SStefano Zampini .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1138633354d9SStefano Zampini           `ISLocalToGlobalMappingGetBlockInfo()`, `ISLocalToGlobalMappingRestoreBlockNodeInfo()`, `ISLocalToGlobalMappingGetNodeInfo()`
1139633354d9SStefano Zampini @*/
1140633354d9SStefano Zampini PetscErrorCode ISLocalToGlobalMappingGetBlockNodeInfo(ISLocalToGlobalMapping mapping, PetscInt *n, PetscInt *n_procs[], PetscInt **procs[])
1141d71ae5a4SJacob Faibussowitsch {
1142633354d9SStefano Zampini   PetscFunctionBegin;
1143633354d9SStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
1144633354d9SStefano Zampini   PetscCall(ISLocalToGlobalMappingSetUpBlockInfo_Private(mapping));
1145633354d9SStefano Zampini   if (n) *n = mapping->n;
1146633354d9SStefano Zampini   if (n_procs) *n_procs = mapping->info_nodec;
1147633354d9SStefano Zampini   if (procs) *procs = mapping->info_nodei;
1148633354d9SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1149633354d9SStefano Zampini }
1150633354d9SStefano Zampini 
1151633354d9SStefano Zampini /*@C
1152633354d9SStefano Zampini   ISLocalToGlobalMappingRestoreBlockNodeInfo - Frees the memory allocated by `ISLocalToGlobalMappingGetBlockNodeInfo()`
1153633354d9SStefano Zampini 
1154633354d9SStefano Zampini   Not Collective
1155633354d9SStefano Zampini 
1156633354d9SStefano Zampini   Input Parameters:
1157633354d9SStefano Zampini + mapping - the mapping from local to global indexing
1158633354d9SStefano Zampini . n       - number of local block nodes
1159633354d9SStefano Zampini . n_procs - an array storing the number of processes for each local block nodes (including self)
1160633354d9SStefano Zampini - procs   - the processes' rank for each local block node (sorted, self is first)
1161633354d9SStefano Zampini 
1162633354d9SStefano Zampini   Level: advanced
1163633354d9SStefano Zampini 
1164633354d9SStefano Zampini .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1165633354d9SStefano Zampini           `ISLocalToGlobalMappingGetBlockNodeInfo()`
1166633354d9SStefano Zampini @*/
1167633354d9SStefano Zampini PetscErrorCode ISLocalToGlobalMappingRestoreBlockNodeInfo(ISLocalToGlobalMapping mapping, PetscInt *n, PetscInt *n_procs[], PetscInt **procs[])
1168633354d9SStefano Zampini {
1169633354d9SStefano Zampini   PetscFunctionBegin;
1170633354d9SStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
1171633354d9SStefano Zampini   if (n) *n = 0;
1172633354d9SStefano Zampini   if (n_procs) *n_procs = NULL;
1173633354d9SStefano Zampini   if (procs) *procs = NULL;
1174633354d9SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1175633354d9SStefano Zampini }
1176633354d9SStefano Zampini 
1177d4df40f3SStefano Zampini /*@C
1178d4df40f3SStefano Zampini   ISLocalToGlobalMappingGetBlockMultiLeavesSF - Get the star-forest to communicate multi-leaf block data
1179d4df40f3SStefano Zampini 
1180d4df40f3SStefano Zampini   Collective the first time it is called
1181d4df40f3SStefano Zampini 
1182d4df40f3SStefano Zampini   Input Parameter:
1183d4df40f3SStefano Zampini . mapping - the mapping from local to global indexing
1184d4df40f3SStefano Zampini 
1185d4df40f3SStefano Zampini   Output Parameter:
1186d4df40f3SStefano Zampini . mlsf - the `PetscSF`
1187d4df40f3SStefano Zampini 
1188d4df40f3SStefano Zampini   Level: advanced
1189d4df40f3SStefano Zampini 
1190d4df40f3SStefano Zampini   Notes:
1191d4df40f3SStefano Zampini   The returned star forest is suitable to exchange local information with other processes sharing the same global block index.
1192d4df40f3SStefano Zampini   For example, suppose a mapping with two processes has been created with
1193d4df40f3SStefano Zampini .vb
1194d4df40f3SStefano Zampini     rank 0 global block indices: [0, 1, 2]
1195d4df40f3SStefano Zampini     rank 1 global block indices: [2, 3, 4]
1196d4df40f3SStefano Zampini .ve
1197d4df40f3SStefano Zampini   and we want to share the local information
1198d4df40f3SStefano Zampini .vb
1199d4df40f3SStefano Zampini     rank 0 data: [-1, -2, -3]
1200d4df40f3SStefano Zampini     rank 1 data: [1, 2, 3]
1201d4df40f3SStefano Zampini .ve
1202d4df40f3SStefano Zampini   then, the broadcasting action of `mlsf` will allow to collect
1203d4df40f3SStefano Zampini .vb
1204d4df40f3SStefano Zampini     rank 0 mlleafdata: [-1, -2, -3, 3]
1205d4df40f3SStefano Zampini     rank 1 mlleafdata: [-3, 3, 1, 2]
1206d4df40f3SStefano Zampini .ve
1207d4df40f3SStefano Zampini   Use ``ISLocalToGlobalMappingGetBlockNodeInfo()`` to index into the multi-leaf data.
1208d4df40f3SStefano Zampini 
1209d4df40f3SStefano Zampini .seealso: [](sec_scatter), `ISLocalToGlobalMappingGetBlockNodeInfo()`, `PetscSF`
1210d4df40f3SStefano Zampini @*/
1211d4df40f3SStefano Zampini PetscErrorCode ISLocalToGlobalMappingGetBlockMultiLeavesSF(ISLocalToGlobalMapping mapping, PetscSF *mlsf)
1212d4df40f3SStefano Zampini {
1213d4df40f3SStefano Zampini   PetscFunctionBegin;
1214d4df40f3SStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
1215d4df40f3SStefano Zampini   PetscAssertPointer(mlsf, 2);
1216d4df40f3SStefano Zampini   PetscCall(ISLocalToGlobalMappingSetUpBlockInfo_Private(mapping));
1217d4df40f3SStefano Zampini   *mlsf = mapping->multileaves_sf;
1218d4df40f3SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1219d4df40f3SStefano Zampini }
1220d4df40f3SStefano Zampini 
1221633354d9SStefano Zampini static PetscErrorCode ISLocalToGlobalMappingSetUpBlockInfo_Private(ISLocalToGlobalMapping mapping)
1222633354d9SStefano Zampini {
1223d4df40f3SStefano Zampini   PetscSF            sf, sf2, imsf, msf;
1224ce94432eSBarry Smith   MPI_Comm           comm;
1225633354d9SStefano Zampini   const PetscSFNode *sfnode;
1226633354d9SStefano Zampini   PetscSFNode       *newsfnode;
1227633354d9SStefano Zampini   PetscLayout        layout;
1228633354d9SStefano Zampini   PetscHMapI         neighs;
1229633354d9SStefano Zampini   PetscHashIter      iter;
1230633354d9SStefano Zampini   PetscBool          missing;
1231633354d9SStefano Zampini   const PetscInt    *gidxs, *rootdegree;
1232633354d9SStefano Zampini   PetscInt          *mask, *mrootdata, *leafdata, *newleafdata, *leafrd, *tmpg;
1233633354d9SStefano Zampini   PetscInt           nroots, nleaves, newnleaves, bs, i, j, m, mnroots, p;
1234633354d9SStefano Zampini   PetscMPIInt        rank, size;
123589d82c54SBarry Smith 
123689d82c54SBarry Smith   PetscFunctionBegin;
1237d4df40f3SStefano Zampini   if (mapping->multileaves_sf) PetscFunctionReturn(PETSC_SUCCESS);
12389566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mapping, &comm));
12399566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
12409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1241633354d9SStefano Zampini 
1242633354d9SStefano Zampini   /* Get mapping indices */
1243633354d9SStefano Zampini   PetscCall(ISLocalToGlobalMappingGetBlockSize(mapping, &bs));
1244633354d9SStefano Zampini   PetscCall(ISLocalToGlobalMappingGetBlockIndices(mapping, &gidxs));
1245633354d9SStefano Zampini   PetscCall(ISLocalToGlobalMappingGetSize(mapping, &nleaves));
1246633354d9SStefano Zampini   nleaves /= bs;
1247633354d9SStefano Zampini 
1248633354d9SStefano Zampini   /* Create layout for global indices */
1249633354d9SStefano Zampini   for (i = 0, m = 0; i < nleaves; i++) m = PetscMax(m, gidxs[i]);
1250462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &m, 1, MPIU_INT, MPI_MAX, comm));
1251633354d9SStefano Zampini   PetscCall(PetscLayoutCreate(comm, &layout));
1252633354d9SStefano Zampini   PetscCall(PetscLayoutSetSize(layout, m + 1));
1253633354d9SStefano Zampini   PetscCall(PetscLayoutSetUp(layout));
1254633354d9SStefano Zampini 
1255633354d9SStefano Zampini   /* Create SF to share global indices */
1256633354d9SStefano Zampini   PetscCall(PetscSFCreate(comm, &sf));
1257633354d9SStefano Zampini   PetscCall(PetscSFSetGraphLayout(sf, layout, nleaves, NULL, PETSC_OWN_POINTER, gidxs));
1258633354d9SStefano Zampini   PetscCall(PetscSFSetUp(sf));
1259633354d9SStefano Zampini   PetscCall(PetscLayoutDestroy(&layout));
1260633354d9SStefano Zampini 
1261633354d9SStefano Zampini   /* communicate root degree to leaves */
1262633354d9SStefano Zampini   PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, &sfnode));
1263633354d9SStefano Zampini   PetscCall(PetscSFComputeDegreeBegin(sf, &rootdegree));
1264633354d9SStefano Zampini   PetscCall(PetscSFComputeDegreeEnd(sf, &rootdegree));
1265633354d9SStefano Zampini   for (i = 0, mnroots = 0; i < nroots; i++) mnroots += rootdegree[i];
1266633354d9SStefano Zampini   PetscCall(PetscMalloc3(2 * PetscMax(mnroots, nroots), &mrootdata, 2 * nleaves, &leafdata, nleaves, &leafrd));
1267633354d9SStefano Zampini   for (i = 0, m = 0; i < nroots; i++) {
1268633354d9SStefano Zampini     mrootdata[2 * i + 0] = rootdegree[i];
1269633354d9SStefano Zampini     mrootdata[2 * i + 1] = m;
1270633354d9SStefano Zampini     m += rootdegree[i];
1271633354d9SStefano Zampini   }
1272633354d9SStefano Zampini   PetscCall(PetscSFBcastBegin(sf, MPIU_2INT, mrootdata, leafdata, MPI_REPLACE));
1273633354d9SStefano Zampini   PetscCall(PetscSFBcastEnd(sf, MPIU_2INT, mrootdata, leafdata, MPI_REPLACE));
1274633354d9SStefano Zampini 
1275633354d9SStefano Zampini   /* allocate enough space to store ranks */
1276633354d9SStefano Zampini   for (i = 0, newnleaves = 0; i < nleaves; i++) {
1277633354d9SStefano Zampini     newnleaves += leafdata[2 * i];
1278633354d9SStefano Zampini     leafrd[i] = leafdata[2 * i];
127924cf384cSBarry Smith   }
128024cf384cSBarry Smith 
1281633354d9SStefano Zampini   /* create new SF nodes to collect multi-root data at leaves */
1282633354d9SStefano Zampini   PetscCall(PetscMalloc1(newnleaves, &newsfnode));
1283633354d9SStefano Zampini   for (i = 0, m = 0; i < nleaves; i++) {
1284633354d9SStefano Zampini     for (j = 0; j < leafrd[i]; j++) {
1285633354d9SStefano Zampini       newsfnode[m].rank  = sfnode[i].rank;
1286633354d9SStefano Zampini       newsfnode[m].index = leafdata[2 * i + 1] + j;
1287633354d9SStefano Zampini       m++;
1288bc8ff85bSBarry Smith     }
1289bc8ff85bSBarry Smith   }
1290bc8ff85bSBarry Smith 
1291633354d9SStefano Zampini   /* gather ranks at multi roots */
1292633354d9SStefano Zampini   for (i = 0; i < mnroots; i++) mrootdata[i] = -1;
1293835f2295SStefano Zampini   for (i = 0; i < nleaves; i++) leafdata[i] = rank;
129430dcb7c9SBarry Smith 
1295633354d9SStefano Zampini   PetscCall(PetscSFGatherBegin(sf, MPIU_INT, leafdata, mrootdata));
1296633354d9SStefano Zampini   PetscCall(PetscSFGatherEnd(sf, MPIU_INT, leafdata, mrootdata));
12973677ff5aSBarry Smith 
1298d4df40f3SStefano Zampini   /* from multi-roots to multi-leaves */
1299d4df40f3SStefano Zampini   PetscCall(PetscSFCreate(comm, &sf2));
1300d4df40f3SStefano Zampini   PetscCall(PetscSFSetGraph(sf2, mnroots, newnleaves, NULL, PETSC_OWN_POINTER, newsfnode, PETSC_OWN_POINTER));
1301d4df40f3SStefano Zampini   PetscCall(PetscSFSetUp(sf2));
1302f6e5521dSKarl Rupp 
1303633354d9SStefano Zampini   /* broadcast multi-root data to multi-leaves */
1304633354d9SStefano Zampini   PetscCall(PetscMalloc1(newnleaves, &newleafdata));
1305d4df40f3SStefano Zampini   PetscCall(PetscSFBcastBegin(sf2, MPIU_INT, mrootdata, newleafdata, MPI_REPLACE));
1306d4df40f3SStefano Zampini   PetscCall(PetscSFBcastEnd(sf2, MPIU_INT, mrootdata, newleafdata, MPI_REPLACE));
1307f6e5521dSKarl Rupp 
1308633354d9SStefano Zampini   /* sort sharing ranks */
1309633354d9SStefano Zampini   for (i = 0, m = 0; i < nleaves; i++) {
1310633354d9SStefano Zampini     PetscCall(PetscSortInt(leafrd[i], newleafdata + m));
1311633354d9SStefano Zampini     m += leafrd[i];
131230dcb7c9SBarry Smith   }
131330dcb7c9SBarry Smith 
1314633354d9SStefano Zampini   /* Number of neighbors and their ranks */
1315633354d9SStefano Zampini   PetscCall(PetscHMapICreate(&neighs));
1316633354d9SStefano Zampini   for (i = 0; i < newnleaves; i++) PetscCall(PetscHMapIPut(neighs, newleafdata[i], &iter, &missing));
1317633354d9SStefano Zampini   PetscCall(PetscHMapIGetSize(neighs, &mapping->info_nproc));
1318*f1957bc3SPierre Jolivet   PetscCall(PetscMalloc1(mapping->info_nproc, &mapping->info_procs));
1319633354d9SStefano Zampini   PetscCall(PetscHMapIGetKeys(neighs, (i = 0, &i), mapping->info_procs));
1320633354d9SStefano Zampini   for (i = 0; i < mapping->info_nproc; i++) { /* put info for self first */
1321633354d9SStefano Zampini     if (mapping->info_procs[i] == rank) {
1322633354d9SStefano Zampini       PetscInt newr = mapping->info_procs[0];
1323d44834fbSBarry Smith 
1324633354d9SStefano Zampini       mapping->info_procs[0] = rank;
1325633354d9SStefano Zampini       mapping->info_procs[i] = newr;
132624cf384cSBarry Smith       break;
132724cf384cSBarry Smith     }
132824cf384cSBarry Smith   }
1329633354d9SStefano Zampini   if (mapping->info_nproc) PetscCall(PetscSortInt(mapping->info_nproc - 1, mapping->info_procs + 1));
1330633354d9SStefano Zampini   PetscCall(PetscHMapIDestroy(&neighs));
1331268a049cSStefano Zampini 
1332633354d9SStefano Zampini   /* collect info data */
1333d4df40f3SStefano Zampini   PetscCall(PetscMalloc1(mapping->info_nproc, &mapping->info_numprocs));
1334d4df40f3SStefano Zampini   PetscCall(PetscMalloc1(mapping->info_nproc, &mapping->info_indices));
1335d4df40f3SStefano Zampini   for (i = 0; i < mapping->info_nproc; i++) mapping->info_indices[i] = NULL;
1336633354d9SStefano Zampini 
1337633354d9SStefano Zampini   PetscCall(PetscMalloc1(nleaves, &mask));
1338633354d9SStefano Zampini   PetscCall(PetscMalloc1(nleaves, &tmpg));
1339633354d9SStefano Zampini   for (p = 0; p < mapping->info_nproc; p++) {
1340633354d9SStefano Zampini     PetscInt *tmp, trank = mapping->info_procs[p];
1341633354d9SStefano Zampini 
1342633354d9SStefano Zampini     PetscCall(PetscMemzero(mask, nleaves * sizeof(*mask)));
1343633354d9SStefano Zampini     for (i = 0, m = 0; i < nleaves; i++) {
1344633354d9SStefano Zampini       for (j = 0; j < leafrd[i]; j++) {
1345633354d9SStefano Zampini         if (newleafdata[m] == trank) mask[i]++;
1346633354d9SStefano Zampini         if (!p && newleafdata[m] != rank) mask[i]++;
1347633354d9SStefano Zampini         m++;
1348633354d9SStefano Zampini       }
1349633354d9SStefano Zampini     }
1350633354d9SStefano Zampini     for (i = 0, m = 0; i < nleaves; i++)
1351633354d9SStefano Zampini       if (mask[i] > (!p ? 1 : 0)) m++;
1352633354d9SStefano Zampini 
1353633354d9SStefano Zampini     PetscCall(PetscMalloc1(m, &tmp));
1354633354d9SStefano Zampini     for (i = 0, m = 0; i < nleaves; i++)
1355633354d9SStefano Zampini       if (mask[i] > (!p ? 1 : 0)) {
1356633354d9SStefano Zampini         tmp[m]  = i;
1357633354d9SStefano Zampini         tmpg[m] = gidxs[i];
1358633354d9SStefano Zampini         m++;
1359633354d9SStefano Zampini       }
1360633354d9SStefano Zampini     PetscCall(PetscSortIntWithArray(m, tmpg, tmp));
1361633354d9SStefano Zampini     mapping->info_indices[p]  = tmp;
1362633354d9SStefano Zampini     mapping->info_numprocs[p] = m;
1363633354d9SStefano Zampini   }
1364633354d9SStefano Zampini 
1365633354d9SStefano Zampini   /* Node info */
1366633354d9SStefano Zampini   PetscCall(PetscMalloc2(nleaves, &mapping->info_nodec, nleaves + 1, &mapping->info_nodei));
1367633354d9SStefano Zampini   PetscCall(PetscArraycpy(mapping->info_nodec, leafrd, nleaves));
1368633354d9SStefano Zampini   PetscCall(PetscMalloc1(newnleaves, &mapping->info_nodei[0]));
1369633354d9SStefano Zampini   for (i = 0; i < nleaves - 1; i++) mapping->info_nodei[i + 1] = mapping->info_nodei[i] + mapping->info_nodec[i];
1370633354d9SStefano Zampini   PetscCall(PetscArraycpy(mapping->info_nodei[0], newleafdata, newnleaves));
1371633354d9SStefano Zampini 
1372d4df40f3SStefano Zampini   /* Create SF from leaves to multi-leaves */
1373d4df40f3SStefano Zampini   PetscCall(PetscSFGetMultiSF(sf, &msf));
1374d4df40f3SStefano Zampini   PetscCall(PetscSFCreateInverseSF(msf, &imsf));
1375d4df40f3SStefano Zampini   PetscCall(PetscSFCompose(imsf, sf2, &mapping->multileaves_sf));
1376d4df40f3SStefano Zampini   PetscCall(PetscSFDestroy(&imsf));
1377d4df40f3SStefano Zampini   PetscCall(PetscSFDestroy(&sf));
1378d4df40f3SStefano Zampini   PetscCall(PetscSFDestroy(&sf2));
1379d4df40f3SStefano Zampini 
1380633354d9SStefano Zampini   PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(mapping, &gidxs));
1381633354d9SStefano Zampini   PetscCall(PetscFree(tmpg));
1382633354d9SStefano Zampini   PetscCall(PetscFree(mask));
1383633354d9SStefano Zampini   PetscCall(PetscFree3(mrootdata, leafdata, leafrd));
1384633354d9SStefano Zampini   PetscCall(PetscFree(newleafdata));
13853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
138689d82c54SBarry Smith }
138789d82c54SBarry Smith 
13886a818285SBarry Smith /*@C
1389cab54364SBarry Smith   ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by `ISLocalToGlobalMappingGetBlockInfo()`
13906a818285SBarry Smith 
1391633354d9SStefano Zampini   Not Collective
13926a818285SBarry Smith 
1393633354d9SStefano Zampini   Input Parameters:
1394633354d9SStefano Zampini + mapping  - the mapping from local to global indexing
1395633354d9SStefano Zampini . nproc    - number of processes that are connected to the calling process
1396633354d9SStefano Zampini . procs    - neighboring processes
1397633354d9SStefano Zampini . numprocs - number of block indices for each process
1398633354d9SStefano Zampini - indices  - block indices (in local numbering) shared with neighbors (sorted by global numbering)
13996a818285SBarry Smith 
14006a818285SBarry Smith   Level: advanced
14016a818285SBarry Smith 
1402cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1403db781477SPatrick Sanan           `ISLocalToGlobalMappingGetInfo()`
14046a818285SBarry Smith @*/
1405d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[])
1406d71ae5a4SJacob Faibussowitsch {
14076a818285SBarry Smith   PetscFunctionBegin;
1408cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
1409633354d9SStefano Zampini   if (nproc) *nproc = 0;
1410633354d9SStefano Zampini   if (procs) *procs = NULL;
1411633354d9SStefano Zampini   if (numprocs) *numprocs = NULL;
1412633354d9SStefano Zampini   if (indices) *indices = NULL;
14133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14146a818285SBarry Smith }
14156a818285SBarry Smith 
14166a818285SBarry Smith /*@C
1417633354d9SStefano Zampini   ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each process
14186a818285SBarry Smith 
1419633354d9SStefano Zampini   Collective the first time it is called
14206a818285SBarry Smith 
1421f899ff85SJose E. Roman   Input Parameter:
14226a818285SBarry Smith . mapping - the mapping from local to global indexing
14236a818285SBarry Smith 
1424d8d19677SJose E. Roman   Output Parameters:
1425633354d9SStefano Zampini + nproc    - number of processes that are connected to the calling process
1426633354d9SStefano Zampini . procs    - neighboring processes
1427633354d9SStefano Zampini . numprocs - number of indices for each process
1428633354d9SStefano Zampini - indices  - indices (in local numbering) shared with neighbors (sorted by global numbering)
14296a818285SBarry Smith 
14306a818285SBarry Smith   Level: advanced
14316a818285SBarry Smith 
1432cab54364SBarry Smith   Note:
1433cab54364SBarry Smith   The user needs to call `ISLocalToGlobalMappingRestoreInfo()` when the data is no longer needed.
14341bd0b88eSStefano Zampini 
143538b5cf2dSJacob Faibussowitsch   Fortran Notes:
1436e33f79d8SJacob Faibussowitsch   There is no `ISLocalToGlobalMappingRestoreInfo()` in Fortran. You must make sure that
1437e33f79d8SJacob Faibussowitsch   `procs`[], `numprocs`[] and `indices`[][] are large enough arrays, either by allocating them
1438e33f79d8SJacob Faibussowitsch   dynamically or defining static ones large enough.
1439e33f79d8SJacob Faibussowitsch 
1440cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1441633354d9SStefano Zampini           `ISLocalToGlobalMappingRestoreInfo()`, `ISLocalToGlobalMappingGetNodeInfo()`
14426a818285SBarry Smith @*/
1443d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[])
1444d71ae5a4SJacob Faibussowitsch {
1445633354d9SStefano Zampini   PetscInt **bindices = NULL, *bnumprocs = NULL, bs, i, j, k, n, *bprocs;
14466a818285SBarry Smith 
14476a818285SBarry Smith   PetscFunctionBegin;
14486a818285SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
14498a1f772fSStefano Zampini   bs = mapping->bs;
1450633354d9SStefano Zampini   PetscCall(ISLocalToGlobalMappingGetBlockInfo(mapping, &n, &bprocs, &bnumprocs, &bindices));
1451268a049cSStefano Zampini   if (bs > 1) { /* we need to expand the cached info */
1452633354d9SStefano Zampini     if (indices) PetscCall(PetscCalloc1(n, indices));
1453633354d9SStefano Zampini     if (numprocs) PetscCall(PetscCalloc1(n, numprocs));
1454633354d9SStefano Zampini     if (indices || numprocs) {
1455633354d9SStefano Zampini       for (i = 0; i < n; i++) {
1456633354d9SStefano Zampini         if (indices) {
14579566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(bs * bnumprocs[i], &(*indices)[i]));
1458268a049cSStefano Zampini           for (j = 0; j < bnumprocs[i]; j++) {
1459ad540459SPierre Jolivet             for (k = 0; k < bs; k++) (*indices)[i][j * bs + k] = bs * bindices[i][j] + k;
14606a818285SBarry Smith           }
14616a818285SBarry Smith         }
1462633354d9SStefano Zampini         if (numprocs) (*numprocs)[i] = bnumprocs[i] * bs;
1463633354d9SStefano Zampini       }
1464633354d9SStefano Zampini     }
1465268a049cSStefano Zampini   } else {
1466633354d9SStefano Zampini     if (numprocs) *numprocs = bnumprocs;
1467633354d9SStefano Zampini     if (indices) *indices = bindices;
14686a818285SBarry Smith   }
1469633354d9SStefano Zampini   if (nproc) *nproc = n;
1470633354d9SStefano Zampini   if (procs) *procs = bprocs;
14713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14726a818285SBarry Smith }
14736a818285SBarry Smith 
147407b52d57SBarry Smith /*@C
1475cab54364SBarry Smith   ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by `ISLocalToGlobalMappingGetInfo()`
147689d82c54SBarry Smith 
1477633354d9SStefano Zampini   Not Collective
147807b52d57SBarry Smith 
1479633354d9SStefano Zampini   Input Parameters:
1480633354d9SStefano Zampini + mapping  - the mapping from local to global indexing
1481633354d9SStefano Zampini . nproc    - number of processes that are connected to the calling process
1482633354d9SStefano Zampini . procs    - neighboring processes
1483633354d9SStefano Zampini . numprocs - number of indices for each process
1484633354d9SStefano Zampini - indices  - indices (in local numbering) shared with neighbors (sorted by global numbering)
148507b52d57SBarry Smith 
148607b52d57SBarry Smith   Level: advanced
148707b52d57SBarry Smith 
1488cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1489db781477SPatrick Sanan           `ISLocalToGlobalMappingGetInfo()`
149007b52d57SBarry Smith @*/
1491d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[])
1492d71ae5a4SJacob Faibussowitsch {
149307b52d57SBarry Smith   PetscFunctionBegin;
1494633354d9SStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
1495633354d9SStefano Zampini   if (mapping->bs > 1) {
1496633354d9SStefano Zampini     if (numprocs) PetscCall(PetscFree(*numprocs));
1497633354d9SStefano Zampini     if (indices) {
1498633354d9SStefano Zampini       if (*indices)
1499633354d9SStefano Zampini         for (PetscInt i = 0; i < *nproc; i++) PetscCall(PetscFree((*indices)[i]));
1500633354d9SStefano Zampini       PetscCall(PetscFree(*indices));
1501633354d9SStefano Zampini     }
1502633354d9SStefano Zampini   }
15033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
150407b52d57SBarry Smith }
150586994e45SJed Brown 
150686994e45SJed Brown /*@C
1507633354d9SStefano Zampini   ISLocalToGlobalMappingGetNodeInfo - Gets the neighbor information of local nodes
15081bd0b88eSStefano Zampini 
1509633354d9SStefano Zampini   Collective the first time it is called
15101bd0b88eSStefano Zampini 
1511f899ff85SJose E. Roman   Input Parameter:
15121bd0b88eSStefano Zampini . mapping - the mapping from local to global indexing
15131bd0b88eSStefano Zampini 
1514d8d19677SJose E. Roman   Output Parameters:
1515633354d9SStefano Zampini + n       - number of local nodes
1516633354d9SStefano Zampini . n_procs - an array storing the number of processes for each local node (including self)
1517633354d9SStefano Zampini - procs   - the processes' rank for each local node (sorted, self is first)
15181bd0b88eSStefano Zampini 
15191bd0b88eSStefano Zampini   Level: advanced
15201bd0b88eSStefano Zampini 
1521cab54364SBarry Smith   Note:
1522633354d9SStefano Zampini   The user needs to call `ISLocalToGlobalMappingRestoreNodeInfo()` when the data is no longer needed.
15231bd0b88eSStefano Zampini 
1524cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1525633354d9SStefano Zampini           `ISLocalToGlobalMappingGetInfo()`, `ISLocalToGlobalMappingRestoreNodeInfo()`, `ISLocalToGlobalMappingGetBlockNodeInfo()`
15261bd0b88eSStefano Zampini @*/
1527633354d9SStefano Zampini PetscErrorCode ISLocalToGlobalMappingGetNodeInfo(ISLocalToGlobalMapping mapping, PetscInt *n, PetscInt *n_procs[], PetscInt **procs[])
1528d71ae5a4SJacob Faibussowitsch {
1529633354d9SStefano Zampini   PetscInt **bprocs = NULL, *bn_procs = NULL, bs, i, j, k, bn;
15301bd0b88eSStefano Zampini 
15311bd0b88eSStefano Zampini   PetscFunctionBegin;
15321bd0b88eSStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
1533633354d9SStefano Zampini   bs = mapping->bs;
1534633354d9SStefano Zampini   PetscCall(ISLocalToGlobalMappingGetBlockNodeInfo(mapping, &bn, &bn_procs, &bprocs));
1535633354d9SStefano Zampini   if (bs > 1) { /* we need to expand the cached info */
1536633354d9SStefano Zampini     PetscInt *tn_procs;
1537633354d9SStefano Zampini     PetscInt  c;
15381bd0b88eSStefano Zampini 
1539633354d9SStefano Zampini     PetscCall(PetscMalloc1(bn * bs, &tn_procs));
1540633354d9SStefano Zampini     for (i = 0, c = 0; i < bn; i++) {
1541633354d9SStefano Zampini       for (k = 0; k < bs; k++) tn_procs[i * bs + k] = bn_procs[i];
1542633354d9SStefano Zampini       c += bs * bn_procs[i];
1543633354d9SStefano Zampini     }
1544633354d9SStefano Zampini     if (n) *n = bn * bs;
1545633354d9SStefano Zampini     if (procs) {
1546633354d9SStefano Zampini       PetscInt **tprocs;
1547633354d9SStefano Zampini       PetscInt   tn = bn * bs;
15481bd0b88eSStefano Zampini 
1549633354d9SStefano Zampini       PetscCall(PetscMalloc1(tn, &tprocs));
1550633354d9SStefano Zampini       if (tn) PetscCall(PetscMalloc1(c, &tprocs[0]));
1551633354d9SStefano Zampini       for (i = 0; i < tn - 1; i++) tprocs[i + 1] = tprocs[i] + tn_procs[i];
1552633354d9SStefano Zampini       for (i = 0; i < bn; i++) {
1553633354d9SStefano Zampini         for (k = 0; k < bs; k++) {
1554633354d9SStefano Zampini           for (j = 0; j < bn_procs[i]; j++) tprocs[i * bs + k][j] = bprocs[i][j];
15551bd0b88eSStefano Zampini         }
15561bd0b88eSStefano Zampini       }
1557633354d9SStefano Zampini       *procs = tprocs;
15581bd0b88eSStefano Zampini     }
1559633354d9SStefano Zampini     if (n_procs) *n_procs = tn_procs;
1560633354d9SStefano Zampini     else PetscCall(PetscFree(tn_procs));
1561633354d9SStefano Zampini   } else {
1562633354d9SStefano Zampini     if (n) *n = bn;
1563633354d9SStefano Zampini     if (n_procs) *n_procs = bn_procs;
1564633354d9SStefano Zampini     if (procs) *procs = bprocs;
1565633354d9SStefano Zampini   }
15663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15671bd0b88eSStefano Zampini }
15681bd0b88eSStefano Zampini 
15691bd0b88eSStefano Zampini /*@C
1570cab54364SBarry Smith   ISLocalToGlobalMappingRestoreNodeInfo - Frees the memory allocated by `ISLocalToGlobalMappingGetNodeInfo()`
15711bd0b88eSStefano Zampini 
1572633354d9SStefano Zampini   Not Collective
15731bd0b88eSStefano Zampini 
1574633354d9SStefano Zampini   Input Parameters:
1575633354d9SStefano Zampini + mapping - the mapping from local to global indexing
1576633354d9SStefano Zampini . n       - number of local nodes
1577633354d9SStefano Zampini . n_procs - an array storing the number of processes for each local node (including self)
1578633354d9SStefano Zampini - procs   - the processes' rank for each local node (sorted, self is first)
15791bd0b88eSStefano Zampini 
15801bd0b88eSStefano Zampini   Level: advanced
15811bd0b88eSStefano Zampini 
1582cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1583db781477SPatrick Sanan           `ISLocalToGlobalMappingGetInfo()`
15841bd0b88eSStefano Zampini @*/
1585633354d9SStefano Zampini PetscErrorCode ISLocalToGlobalMappingRestoreNodeInfo(ISLocalToGlobalMapping mapping, PetscInt *n, PetscInt *n_procs[], PetscInt **procs[])
1586d71ae5a4SJacob Faibussowitsch {
15871bd0b88eSStefano Zampini   PetscFunctionBegin;
15881bd0b88eSStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
1589633354d9SStefano Zampini   if (mapping->bs > 1) {
1590633354d9SStefano Zampini     if (n_procs) PetscCall(PetscFree(*n_procs));
1591633354d9SStefano Zampini     if (procs) {
1592633354d9SStefano Zampini       if (*procs) PetscCall(PetscFree((*procs)[0]));
1593633354d9SStefano Zampini       PetscCall(PetscFree(*procs));
1594633354d9SStefano Zampini     }
1595633354d9SStefano Zampini   }
1596633354d9SStefano Zampini   PetscCall(ISLocalToGlobalMappingRestoreBlockNodeInfo(mapping, n, n_procs, procs));
15973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15981bd0b88eSStefano Zampini }
15991bd0b88eSStefano Zampini 
16001bd0b88eSStefano Zampini /*@C
1601107e9a97SBarry Smith   ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped
160286994e45SJed Brown 
160386994e45SJed Brown   Not Collective
160486994e45SJed Brown 
16054165533cSJose E. Roman   Input Parameter:
160686994e45SJed Brown . ltog - local to global mapping
160786994e45SJed Brown 
16084165533cSJose E. Roman   Output Parameter:
1609cab54364SBarry Smith . array - array of indices, the length of this array may be obtained with `ISLocalToGlobalMappingGetSize()`
161086994e45SJed Brown 
161186994e45SJed Brown   Level: advanced
161286994e45SJed Brown 
1613cab54364SBarry Smith   Note:
1614cab54364SBarry Smith   `ISLocalToGlobalMappingGetSize()` returns the length the this array
1615107e9a97SBarry Smith 
161620662ed9SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingRestoreIndices()`,
161720662ed9SBarry Smith           `ISLocalToGlobalMappingGetBlockIndices()`, `ISLocalToGlobalMappingRestoreBlockIndices()`
161886994e45SJed Brown @*/
1619cc4c1da9SBarry Smith PetscErrorCode ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog, const PetscInt *array[])
1620d71ae5a4SJacob Faibussowitsch {
162186994e45SJed Brown   PetscFunctionBegin;
162286994e45SJed Brown   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
16234f572ea9SToby Isaac   PetscAssertPointer(array, 2);
162445b6f7e9SBarry Smith   if (ltog->bs == 1) {
162586994e45SJed Brown     *array = ltog->indices;
162645b6f7e9SBarry Smith   } else {
162745b6f7e9SBarry Smith     PetscInt       *jj, k, i, j, n = ltog->n, bs = ltog->bs;
162845b6f7e9SBarry Smith     const PetscInt *ii;
162945b6f7e9SBarry Smith 
16309566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(bs * n, &jj));
163145b6f7e9SBarry Smith     *array = jj;
163245b6f7e9SBarry Smith     k      = 0;
163345b6f7e9SBarry Smith     ii     = ltog->indices;
163445b6f7e9SBarry Smith     for (i = 0; i < n; i++)
16359371c9d4SSatish Balay       for (j = 0; j < bs; j++) jj[k++] = bs * ii[i] + j;
163645b6f7e9SBarry Smith   }
16373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
163886994e45SJed Brown }
163986994e45SJed Brown 
164086994e45SJed Brown /*@C
1641cab54364SBarry Smith   ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with `ISLocalToGlobalMappingGetIndices()`
164286994e45SJed Brown 
164386994e45SJed Brown   Not Collective
164486994e45SJed Brown 
16454165533cSJose E. Roman   Input Parameters:
164686994e45SJed Brown + ltog  - local to global mapping
164786994e45SJed Brown - array - array of indices
164886994e45SJed Brown 
164986994e45SJed Brown   Level: advanced
165086994e45SJed Brown 
1651cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingGetIndices()`
165286994e45SJed Brown @*/
1653cc4c1da9SBarry Smith PetscErrorCode ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog, const PetscInt *array[])
1654d71ae5a4SJacob Faibussowitsch {
165586994e45SJed Brown   PetscFunctionBegin;
165686994e45SJed Brown   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
16574f572ea9SToby Isaac   PetscAssertPointer(array, 2);
1658c9cc58a2SBarry Smith   PetscCheck(ltog->bs != 1 || *array == ltog->indices, PETSC_COMM_SELF, PETSC_ERR_ARG_BADPTR, "Trying to return mismatched pointer");
165948a46eb9SPierre Jolivet   if (ltog->bs > 1) PetscCall(PetscFree(*(void **)array));
16603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
166145b6f7e9SBarry Smith }
166245b6f7e9SBarry Smith 
166345b6f7e9SBarry Smith /*@C
1664ce78bad3SBarry Smith   ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block in a `ISLocalToGlobalMapping`
166545b6f7e9SBarry Smith 
166645b6f7e9SBarry Smith   Not Collective
166745b6f7e9SBarry Smith 
16684165533cSJose E. Roman   Input Parameter:
166945b6f7e9SBarry Smith . ltog - local to global mapping
167045b6f7e9SBarry Smith 
16714165533cSJose E. Roman   Output Parameter:
167245b6f7e9SBarry Smith . array - array of indices
167345b6f7e9SBarry Smith 
167445b6f7e9SBarry Smith   Level: advanced
167545b6f7e9SBarry Smith 
1676ce78bad3SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`,
1677ce78bad3SBarry Smith           `ISLocalToGlobalMappingRestoreBlockIndices()`
167845b6f7e9SBarry Smith @*/
1679cc4c1da9SBarry Smith PetscErrorCode ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog, const PetscInt *array[])
1680d71ae5a4SJacob Faibussowitsch {
168145b6f7e9SBarry Smith   PetscFunctionBegin;
168245b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
16834f572ea9SToby Isaac   PetscAssertPointer(array, 2);
168445b6f7e9SBarry Smith   *array = ltog->indices;
16853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
168645b6f7e9SBarry Smith }
168745b6f7e9SBarry Smith 
168845b6f7e9SBarry Smith /*@C
1689cab54364SBarry Smith   ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with `ISLocalToGlobalMappingGetBlockIndices()`
169045b6f7e9SBarry Smith 
169145b6f7e9SBarry Smith   Not Collective
169245b6f7e9SBarry Smith 
16934165533cSJose E. Roman   Input Parameters:
169445b6f7e9SBarry Smith + ltog  - local to global mapping
169545b6f7e9SBarry Smith - array - array of indices
169645b6f7e9SBarry Smith 
169745b6f7e9SBarry Smith   Level: advanced
169845b6f7e9SBarry Smith 
1699cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingGetIndices()`
170045b6f7e9SBarry Smith @*/
1701cc4c1da9SBarry Smith PetscErrorCode ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog, const PetscInt *array[])
1702d71ae5a4SJacob Faibussowitsch {
170345b6f7e9SBarry Smith   PetscFunctionBegin;
170445b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
17054f572ea9SToby Isaac   PetscAssertPointer(array, 2);
170608401ef6SPierre Jolivet   PetscCheck(*array == ltog->indices, PETSC_COMM_SELF, PETSC_ERR_ARG_BADPTR, "Trying to return mismatched pointer");
17070298fd71SBarry Smith   *array = NULL;
17083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
170986994e45SJed Brown }
1710f7efa3c7SJed Brown 
1711cc4c1da9SBarry Smith /*@
1712f7efa3c7SJed Brown   ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings
1713f7efa3c7SJed Brown 
1714f7efa3c7SJed Brown   Not Collective
1715f7efa3c7SJed Brown 
17164165533cSJose E. Roman   Input Parameters:
1717f7efa3c7SJed Brown + comm  - communicator for the new mapping, must contain the communicator of every mapping to concatenate
1718f7efa3c7SJed Brown . n     - number of mappings to concatenate
1719f7efa3c7SJed Brown - ltogs - local to global mappings
1720f7efa3c7SJed Brown 
17214165533cSJose E. Roman   Output Parameter:
1722f7efa3c7SJed Brown . ltogcat - new mapping
1723f7efa3c7SJed Brown 
1724f7efa3c7SJed Brown   Level: advanced
1725f7efa3c7SJed Brown 
1726cab54364SBarry Smith   Note:
1727cab54364SBarry Smith   This currently always returns a mapping with block size of 1
1728cab54364SBarry Smith 
172938b5cf2dSJacob Faibussowitsch   Developer Notes:
1730cab54364SBarry Smith   If all the input mapping have the same block size we could easily handle that as a special case
1731cab54364SBarry Smith 
1732cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingCreate()`
1733f7efa3c7SJed Brown @*/
1734d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm, PetscInt n, const ISLocalToGlobalMapping ltogs[], ISLocalToGlobalMapping *ltogcat)
1735d71ae5a4SJacob Faibussowitsch {
1736f7efa3c7SJed Brown   PetscInt i, cnt, m, *idx;
1737f7efa3c7SJed Brown 
1738f7efa3c7SJed Brown   PetscFunctionBegin;
173908401ef6SPierre Jolivet   PetscCheck(n >= 0, comm, PETSC_ERR_ARG_OUTOFRANGE, "Must have a non-negative number of mappings, given %" PetscInt_FMT, n);
17404f572ea9SToby Isaac   if (n > 0) PetscAssertPointer(ltogs, 3);
1741f7efa3c7SJed Brown   for (i = 0; i < n; i++) PetscValidHeaderSpecific(ltogs[i], IS_LTOGM_CLASSID, 3);
17424f572ea9SToby Isaac   PetscAssertPointer(ltogcat, 4);
1743f7efa3c7SJed Brown   for (cnt = 0, i = 0; i < n; i++) {
17449566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(ltogs[i], &m));
1745f7efa3c7SJed Brown     cnt += m;
1746f7efa3c7SJed Brown   }
17479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(cnt, &idx));
1748f7efa3c7SJed Brown   for (cnt = 0, i = 0; i < n; i++) {
1749f7efa3c7SJed Brown     const PetscInt *subidx;
17509566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(ltogs[i], &m));
17519566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetIndices(ltogs[i], &subidx));
17529566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(&idx[cnt], subidx, m));
17539566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreIndices(ltogs[i], &subidx));
1754f7efa3c7SJed Brown     cnt += m;
1755f7efa3c7SJed Brown   }
17569566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(comm, 1, cnt, idx, PETSC_OWN_POINTER, ltogcat));
17573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1758f7efa3c7SJed Brown }
175904a59952SBarry Smith 
1760413f72f0SBarry Smith /*MC
1761cab54364SBarry Smith       ISLOCALTOGLOBALMAPPINGBASIC - basic implementation of the `ISLocalToGlobalMapping` object. When `ISGlobalToLocalMappingApply()` is
1762413f72f0SBarry Smith                                     used this is good for only small and moderate size problems.
1763413f72f0SBarry Smith 
176420662ed9SBarry Smith    Options Database Key:
1765a2b725a8SWilliam Gropp .   -islocaltoglobalmapping_type basic - select this method
1766413f72f0SBarry Smith 
1767413f72f0SBarry Smith    Level: beginner
1768413f72f0SBarry Smith 
1769cab54364SBarry Smith    Developer Note:
1770cab54364SBarry Smith    This stores all the mapping information on each MPI rank.
1771cab54364SBarry Smith 
1772cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`, `ISLOCALTOGLOBALMAPPINGHASH`
1773413f72f0SBarry Smith M*/
1774d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Basic(ISLocalToGlobalMapping ltog)
1775d71ae5a4SJacob Faibussowitsch {
1776413f72f0SBarry Smith   PetscFunctionBegin;
1777413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Basic;
1778413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Basic;
1779413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Basic;
1780413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Basic;
17813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1782413f72f0SBarry Smith }
1783413f72f0SBarry Smith 
1784413f72f0SBarry Smith /*MC
1785cab54364SBarry Smith       ISLOCALTOGLOBALMAPPINGHASH - hash implementation of the `ISLocalToGlobalMapping` object. When `ISGlobalToLocalMappingApply()` is
1786ed56e8eaSBarry Smith                                     used this is good for large memory problems.
1787413f72f0SBarry Smith 
178820662ed9SBarry Smith    Options Database Key:
1789a2b725a8SWilliam Gropp .   -islocaltoglobalmapping_type hash - select this method
1790413f72f0SBarry Smith 
1791413f72f0SBarry Smith    Level: beginner
1792413f72f0SBarry Smith 
1793cab54364SBarry Smith    Note:
1794cab54364SBarry Smith     This is selected automatically for large problems if the user does not set the type.
1795cab54364SBarry Smith 
1796cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`, `ISLOCALTOGLOBALMAPPINGBASIC`
1797413f72f0SBarry Smith M*/
1798d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Hash(ISLocalToGlobalMapping ltog)
1799d71ae5a4SJacob Faibussowitsch {
1800413f72f0SBarry Smith   PetscFunctionBegin;
1801413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Hash;
1802413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Hash;
1803413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Hash;
1804413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Hash;
18053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1806413f72f0SBarry Smith }
1807413f72f0SBarry Smith 
1808413f72f0SBarry Smith /*@C
1809cab54364SBarry Smith   ISLocalToGlobalMappingRegister -  Registers a method for applying a global to local mapping with an `ISLocalToGlobalMapping`
1810413f72f0SBarry Smith 
1811cc4c1da9SBarry Smith   Not Collective, No Fortran Support
1812413f72f0SBarry Smith 
1813413f72f0SBarry Smith   Input Parameters:
1814413f72f0SBarry Smith + sname    - name of a new method
18152fe279fdSBarry Smith - function - routine to create method context
1816413f72f0SBarry Smith 
181738b5cf2dSJacob Faibussowitsch   Example Usage:
1818413f72f0SBarry Smith .vb
1819413f72f0SBarry Smith    ISLocalToGlobalMappingRegister("my_mapper", MyCreate);
1820413f72f0SBarry Smith .ve
1821413f72f0SBarry Smith 
1822ed56e8eaSBarry Smith   Then, your mapping can be chosen with the procedural interface via
1823b44f4de4SBarry Smith .vb
1824b44f4de4SBarry Smith   ISLocalToGlobalMappingSetType(ltog, "my_mapper")
1825b44f4de4SBarry Smith .ve
1826413f72f0SBarry Smith   or at runtime via the option
1827b44f4de4SBarry Smith .vb
1828b44f4de4SBarry Smith   -islocaltoglobalmapping_type my_mapper
1829b44f4de4SBarry Smith .ve
1830413f72f0SBarry Smith 
1831413f72f0SBarry Smith   Level: advanced
1832413f72f0SBarry Smith 
1833cab54364SBarry Smith   Note:
1834cab54364SBarry Smith   `ISLocalToGlobalMappingRegister()` may be called multiple times to add several user-defined mappings.
1835413f72f0SBarry Smith 
1836cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingRegisterAll()`, `ISLocalToGlobalMappingRegisterDestroy()`, `ISLOCALTOGLOBALMAPPINGBASIC`,
1837cab54364SBarry Smith           `ISLOCALTOGLOBALMAPPINGHASH`, `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApply()`
1838413f72f0SBarry Smith @*/
1839d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRegister(const char sname[], PetscErrorCode (*function)(ISLocalToGlobalMapping))
1840d71ae5a4SJacob Faibussowitsch {
1841413f72f0SBarry Smith   PetscFunctionBegin;
18429566063dSJacob Faibussowitsch   PetscCall(ISInitializePackage());
18439566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&ISLocalToGlobalMappingList, sname, function));
18443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1845413f72f0SBarry Smith }
1846413f72f0SBarry Smith 
1847cc4c1da9SBarry Smith /*@
1848cab54364SBarry Smith   ISLocalToGlobalMappingSetType - Sets the implementation type `ISLocalToGlobalMapping` will use
1849413f72f0SBarry Smith 
1850c3339decSBarry Smith   Logically Collective
1851413f72f0SBarry Smith 
1852413f72f0SBarry Smith   Input Parameters:
1853cab54364SBarry Smith + ltog - the `ISLocalToGlobalMapping` object
1854413f72f0SBarry Smith - type - a known method
1855413f72f0SBarry Smith 
1856413f72f0SBarry Smith   Options Database Key:
1857cab54364SBarry Smith . -islocaltoglobalmapping_type  <method> - Sets the method; use -help for a list of available methods (for instance, basic or hash)
1858cab54364SBarry Smith 
1859cab54364SBarry Smith   Level: intermediate
1860413f72f0SBarry Smith 
1861413f72f0SBarry Smith   Notes:
186220662ed9SBarry Smith   See `ISLocalToGlobalMappingType` for available methods
1863413f72f0SBarry Smith 
1864cab54364SBarry Smith   Normally, it is best to use the `ISLocalToGlobalMappingSetFromOptions()` command and
1865cab54364SBarry Smith   then set the `ISLocalToGlobalMappingType` from the options database rather than by using
1866413f72f0SBarry Smith   this routine.
1867413f72f0SBarry Smith 
186838b5cf2dSJacob Faibussowitsch   Developer Notes:
1869cab54364SBarry Smith   `ISLocalToGlobalMappingRegister()` is used to add new types to `ISLocalToGlobalMappingList` from which they
1870cab54364SBarry Smith   are accessed by `ISLocalToGlobalMappingSetType()`.
1871413f72f0SBarry Smith 
187238b5cf2dSJacob Faibussowitsch .seealso: [](sec_scatter), `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingRegister()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingGetType()`
1873413f72f0SBarry Smith @*/
1874d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType type)
1875d71ae5a4SJacob Faibussowitsch {
1876413f72f0SBarry Smith   PetscBool match;
18775f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(ISLocalToGlobalMapping) = NULL;
1878413f72f0SBarry Smith 
1879413f72f0SBarry Smith   PetscFunctionBegin;
1880413f72f0SBarry Smith   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
18814f572ea9SToby Isaac   if (type) PetscAssertPointer(type, 2);
1882413f72f0SBarry Smith 
18839566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ltog, type, &match));
18843ba16761SJacob Faibussowitsch   if (match) PetscFunctionReturn(PETSC_SUCCESS);
1885413f72f0SBarry Smith 
1886a0d79125SStefano Zampini   /* L2G maps defer type setup at globaltolocal calls, allow passing NULL here */
1887a0d79125SStefano Zampini   if (type) {
18889566063dSJacob Faibussowitsch     PetscCall(PetscFunctionListFind(ISLocalToGlobalMappingList, type, &r));
1889a0d79125SStefano Zampini     PetscCheck(r, PetscObjectComm((PetscObject)ltog), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested ISLocalToGlobalMapping type %s", type);
1890a0d79125SStefano Zampini   }
1891413f72f0SBarry Smith   /* Destroy the previous private LTOG context */
1892dbbe0bcdSBarry Smith   PetscTryTypeMethod(ltog, destroy);
1893413f72f0SBarry Smith   ltog->ops->destroy = NULL;
1894dbbe0bcdSBarry Smith 
18959566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)ltog, type));
18969566063dSJacob Faibussowitsch   if (r) PetscCall((*r)(ltog));
18973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1898a0d79125SStefano Zampini }
1899a0d79125SStefano Zampini 
1900cc4c1da9SBarry Smith /*@
1901cab54364SBarry Smith   ISLocalToGlobalMappingGetType - Get the type of the `ISLocalToGlobalMapping`
1902a0d79125SStefano Zampini 
1903a0d79125SStefano Zampini   Not Collective
1904a0d79125SStefano Zampini 
1905a0d79125SStefano Zampini   Input Parameter:
1906cab54364SBarry Smith . ltog - the `ISLocalToGlobalMapping` object
1907a0d79125SStefano Zampini 
1908a0d79125SStefano Zampini   Output Parameter:
1909a0d79125SStefano Zampini . type - the type
1910a0d79125SStefano Zampini 
191149762cbcSSatish Balay   Level: intermediate
191249762cbcSSatish Balay 
191338b5cf2dSJacob Faibussowitsch .seealso: [](sec_scatter), `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingRegister()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`
1914a0d79125SStefano Zampini @*/
1915d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType *type)
1916d71ae5a4SJacob Faibussowitsch {
1917a0d79125SStefano Zampini   PetscFunctionBegin;
1918a0d79125SStefano Zampini   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
19194f572ea9SToby Isaac   PetscAssertPointer(type, 2);
1920a0d79125SStefano Zampini   *type = ((PetscObject)ltog)->type_name;
19213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1922413f72f0SBarry Smith }
1923413f72f0SBarry Smith 
1924413f72f0SBarry Smith PetscBool ISLocalToGlobalMappingRegisterAllCalled = PETSC_FALSE;
1925413f72f0SBarry Smith 
1926413f72f0SBarry Smith /*@C
1927cab54364SBarry Smith   ISLocalToGlobalMappingRegisterAll - Registers all of the local to global mapping components in the `IS` package.
1928413f72f0SBarry Smith 
1929413f72f0SBarry Smith   Not Collective
1930413f72f0SBarry Smith 
1931413f72f0SBarry Smith   Level: advanced
1932413f72f0SBarry Smith 
1933cab54364SBarry Smith .seealso: [](sec_scatter), `ISRegister()`, `ISLocalToGlobalRegister()`
1934413f72f0SBarry Smith @*/
1935d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRegisterAll(void)
1936d71ae5a4SJacob Faibussowitsch {
1937413f72f0SBarry Smith   PetscFunctionBegin;
19383ba16761SJacob Faibussowitsch   if (ISLocalToGlobalMappingRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
1939413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAllCalled = PETSC_TRUE;
19409566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGBASIC, ISLocalToGlobalMappingCreate_Basic));
19419566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGHASH, ISLocalToGlobalMappingCreate_Hash));
19423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1943413f72f0SBarry Smith }
1944