xref: /petsc/src/vec/is/utils/isltog.c (revision f2c6b1a247e0aba1e6cff92019aae48a2a13617a)
12362add9SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/isimpl.h> /*I "petscis.h"  I*/
3e8f14785SLisandro Dalcin #include <petsc/private/hashmapi.h>
40c312b8eSJed Brown #include <petscsf.h>
5665c2dedSJed Brown #include <petscviewer.h>
62362add9SBarry Smith 
77087cfbeSBarry Smith PetscClassId          IS_LTOGM_CLASSID;
8268a049cSStefano Zampini static PetscErrorCode ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping, PetscInt *, PetscInt **, PetscInt **, PetscInt ***);
98e58c17dSMatthew Knepley 
10413f72f0SBarry Smith typedef struct {
11413f72f0SBarry Smith   PetscInt *globals;
12413f72f0SBarry Smith } ISLocalToGlobalMapping_Basic;
13413f72f0SBarry Smith 
14413f72f0SBarry Smith typedef struct {
15e8f14785SLisandro Dalcin   PetscHMapI globalht;
16413f72f0SBarry Smith } ISLocalToGlobalMapping_Hash;
17413f72f0SBarry Smith 
186528b96dSMatthew G. Knepley /*@C
19cab54364SBarry Smith   ISGetPointRange - Returns a description of the points in an `IS` suitable for traversal
20413f72f0SBarry Smith 
216528b96dSMatthew G. Knepley   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:
34cab54364SBarry Smith   If the `IS` contains contiguous indices in an `ISSTRIDE`, then the indices are contained in [pStart, pEnd) and points = NULL. Otherwise, pStart = 0, pEnd = numIndices, and points is an array of the indices. This supports the following pattern
35cab54364SBarry Smith .vb
36cab54364SBarry Smith   ISGetPointRange(is, &pStart, &pEnd, &points);
37cab54364SBarry Smith   for (p = pStart; p < pEnd; ++p) {
38cab54364SBarry Smith     const PetscInt point = points ? points[p] : p;
39cab54364SBarry Smith   }
40cab54364SBarry Smith   ISRestorePointRange(is, &pstart, &pEnd, &points);
41cab54364SBarry Smith .ve
42cab54364SBarry Smith 
43cab54364SBarry Smith .seealso: [](sec_scatter), `IS`, `ISRestorePointRange()`, `ISGetPointSubrange()`, `ISGetIndices()`, `ISCreateStride()`
446528b96dSMatthew G. Knepley @*/
45d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGetPointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points)
46d71ae5a4SJacob Faibussowitsch {
479305a4c7SMatthew G. Knepley   PetscInt  numCells, step = 1;
489305a4c7SMatthew G. Knepley   PetscBool isStride;
499305a4c7SMatthew G. Knepley 
509305a4c7SMatthew G. Knepley   PetscFunctionBeginHot;
519305a4c7SMatthew G. Knepley   *pStart = 0;
529305a4c7SMatthew G. Knepley   *points = NULL;
539566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(pointIS, &numCells));
549566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pointIS, ISSTRIDE, &isStride));
559566063dSJacob Faibussowitsch   if (isStride) PetscCall(ISStrideGetInfo(pointIS, pStart, &step));
569305a4c7SMatthew G. Knepley   *pEnd = *pStart + numCells;
579566063dSJacob Faibussowitsch   if (!isStride || step != 1) PetscCall(ISGetIndices(pointIS, points));
589305a4c7SMatthew G. Knepley   PetscFunctionReturn(0);
599305a4c7SMatthew G. Knepley }
609305a4c7SMatthew G. Knepley 
616528b96dSMatthew G. Knepley /*@C
62cab54364SBarry Smith   ISRestorePointRange - Destroys the traversal description created with `ISGetPointRange()`
636528b96dSMatthew G. Knepley 
646528b96dSMatthew G. Knepley   Not collective
656528b96dSMatthew G. Knepley 
666528b96dSMatthew G. Knepley   Input Parameters:
67cab54364SBarry Smith + pointIS - The `IS` object
68cab54364SBarry Smith . pStart  - The first index, from `ISGetPointRange()`
69cab54364SBarry Smith . pEnd    - One past the last index, from `ISGetPointRange()`
70cab54364SBarry Smith - points  - The indices, from `ISGetPointRange()`
716528b96dSMatthew G. Knepley 
726528b96dSMatthew G. Knepley   Level: intermediate
736528b96dSMatthew G. Knepley 
74cab54364SBarry Smith .seealso: [](sec_scatter), `IS`, `ISGetPointRange()`, `ISGetPointSubrange()`, `ISGetIndices()`, `ISCreateStride()`
756528b96dSMatthew G. Knepley @*/
76d71ae5a4SJacob Faibussowitsch PetscErrorCode ISRestorePointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points)
77d71ae5a4SJacob Faibussowitsch {
789305a4c7SMatthew G. Knepley   PetscInt  step = 1;
799305a4c7SMatthew G. Knepley   PetscBool isStride;
809305a4c7SMatthew G. Knepley 
819305a4c7SMatthew G. Knepley   PetscFunctionBeginHot;
829566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pointIS, ISSTRIDE, &isStride));
839566063dSJacob Faibussowitsch   if (isStride) PetscCall(ISStrideGetInfo(pointIS, pStart, &step));
849566063dSJacob Faibussowitsch   if (!isStride || step != 1) PetscCall(ISGetIndices(pointIS, points));
859305a4c7SMatthew G. Knepley   PetscFunctionReturn(0);
869305a4c7SMatthew G. Knepley }
879305a4c7SMatthew G. Knepley 
886528b96dSMatthew G. Knepley /*@C
89cab54364SBarry Smith   ISGetPointSubrange - Configures the input `IS` to be a subrange for the traversal information given
906528b96dSMatthew G. Knepley 
916528b96dSMatthew G. Knepley   Not collective
926528b96dSMatthew G. Knepley 
936528b96dSMatthew G. Knepley   Input Parameters:
94cab54364SBarry Smith + subpointIS - The `IS` object to be configured
956528b96dSMatthew G. Knepley . pStar   t  - The first index of the subrange
966528b96dSMatthew G. Knepley . pEnd       - One past the last index for the subrange
97cab54364SBarry Smith - points     - The indices for the entire range, from `ISGetPointRange()`
986528b96dSMatthew G. Knepley 
996528b96dSMatthew G. Knepley   Output Parameters:
100cab54364SBarry Smith . subpointIS - The `IS` object now configured to be a subrange
1016528b96dSMatthew G. Knepley 
1026528b96dSMatthew G. Knepley   Level: intermediate
1036528b96dSMatthew G. Knepley 
104cab54364SBarry Smith   Note:
105cab54364SBarry Smith   The input `IS` will now respond properly to calls to `ISGetPointRange()` and return the subrange.
106cab54364SBarry Smith 
107cab54364SBarry Smith .seealso: [](sec_scatter), `IS`, `ISGetPointRange()`, `ISRestorePointRange()`, `ISGetIndices()`, `ISCreateStride()`
1086528b96dSMatthew G. Knepley @*/
109d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGetPointSubrange(IS subpointIS, PetscInt pStart, PetscInt pEnd, const PetscInt *points)
110d71ae5a4SJacob Faibussowitsch {
1119305a4c7SMatthew G. Knepley   PetscFunctionBeginHot;
1129305a4c7SMatthew G. Knepley   if (points) {
1139566063dSJacob Faibussowitsch     PetscCall(ISSetType(subpointIS, ISGENERAL));
1149566063dSJacob Faibussowitsch     PetscCall(ISGeneralSetIndices(subpointIS, pEnd - pStart, &points[pStart], PETSC_USE_POINTER));
1159305a4c7SMatthew G. Knepley   } else {
1169566063dSJacob Faibussowitsch     PetscCall(ISSetType(subpointIS, ISSTRIDE));
1179566063dSJacob Faibussowitsch     PetscCall(ISStrideSetStride(subpointIS, pEnd - pStart, pStart, 1));
1189305a4c7SMatthew G. Knepley   }
1199305a4c7SMatthew G. Knepley   PetscFunctionReturn(0);
1209305a4c7SMatthew G. Knepley }
1219305a4c7SMatthew G. Knepley 
122413f72f0SBarry Smith /* -----------------------------------------------------------------------------------------*/
123413f72f0SBarry Smith 
124413f72f0SBarry Smith /*
125413f72f0SBarry Smith     Creates the global mapping information in the ISLocalToGlobalMapping structure
126413f72f0SBarry Smith 
127413f72f0SBarry Smith     If the user has not selected how to handle the global to local mapping then use HASH for "large" problems
128413f72f0SBarry Smith */
129d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISGlobalToLocalMappingSetUp(ISLocalToGlobalMapping mapping)
130d71ae5a4SJacob Faibussowitsch {
131413f72f0SBarry Smith   PetscInt i, *idx = mapping->indices, n = mapping->n, end, start;
132413f72f0SBarry Smith 
133413f72f0SBarry Smith   PetscFunctionBegin;
134413f72f0SBarry Smith   if (mapping->data) PetscFunctionReturn(0);
135413f72f0SBarry Smith   end   = 0;
136413f72f0SBarry Smith   start = PETSC_MAX_INT;
137413f72f0SBarry Smith 
138413f72f0SBarry Smith   for (i = 0; i < n; i++) {
139413f72f0SBarry Smith     if (idx[i] < 0) continue;
140413f72f0SBarry Smith     if (idx[i] < start) start = idx[i];
141413f72f0SBarry Smith     if (idx[i] > end) end = idx[i];
142413f72f0SBarry Smith   }
1439371c9d4SSatish Balay   if (start > end) {
1449371c9d4SSatish Balay     start = 0;
1459371c9d4SSatish Balay     end   = -1;
1469371c9d4SSatish Balay   }
147413f72f0SBarry Smith   mapping->globalstart = start;
148413f72f0SBarry Smith   mapping->globalend   = end;
149413f72f0SBarry Smith   if (!((PetscObject)mapping)->type_name) {
150413f72f0SBarry Smith     if ((end - start) > PetscMax(4 * n, 1000000)) {
1519566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingSetType(mapping, ISLOCALTOGLOBALMAPPINGHASH));
152413f72f0SBarry Smith     } else {
1539566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingSetType(mapping, ISLOCALTOGLOBALMAPPINGBASIC));
154413f72f0SBarry Smith     }
155413f72f0SBarry Smith   }
156dbbe0bcdSBarry Smith   PetscTryTypeMethod(mapping, globaltolocalmappingsetup);
157413f72f0SBarry Smith   PetscFunctionReturn(0);
158413f72f0SBarry Smith }
159413f72f0SBarry Smith 
160d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISGlobalToLocalMappingSetUp_Basic(ISLocalToGlobalMapping mapping)
161d71ae5a4SJacob Faibussowitsch {
162413f72f0SBarry Smith   PetscInt                      i, *idx = mapping->indices, n = mapping->n, end, start, *globals;
163413f72f0SBarry Smith   ISLocalToGlobalMapping_Basic *map;
164413f72f0SBarry Smith 
165413f72f0SBarry Smith   PetscFunctionBegin;
166413f72f0SBarry Smith   start = mapping->globalstart;
167413f72f0SBarry Smith   end   = mapping->globalend;
1689566063dSJacob Faibussowitsch   PetscCall(PetscNew(&map));
1699566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(end - start + 2, &globals));
170413f72f0SBarry Smith   map->globals = globals;
171413f72f0SBarry Smith   for (i = 0; i < end - start + 1; i++) globals[i] = -1;
172413f72f0SBarry Smith   for (i = 0; i < n; i++) {
173413f72f0SBarry Smith     if (idx[i] < 0) continue;
174413f72f0SBarry Smith     globals[idx[i] - start] = i;
175413f72f0SBarry Smith   }
176413f72f0SBarry Smith   mapping->data = (void *)map;
177413f72f0SBarry Smith   PetscFunctionReturn(0);
178413f72f0SBarry Smith }
179413f72f0SBarry Smith 
180d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISGlobalToLocalMappingSetUp_Hash(ISLocalToGlobalMapping mapping)
181d71ae5a4SJacob Faibussowitsch {
182413f72f0SBarry Smith   PetscInt                     i, *idx = mapping->indices, n = mapping->n;
183413f72f0SBarry Smith   ISLocalToGlobalMapping_Hash *map;
184413f72f0SBarry Smith 
185413f72f0SBarry Smith   PetscFunctionBegin;
1869566063dSJacob Faibussowitsch   PetscCall(PetscNew(&map));
1879566063dSJacob Faibussowitsch   PetscCall(PetscHMapICreate(&map->globalht));
188413f72f0SBarry Smith   for (i = 0; i < n; i++) {
189413f72f0SBarry Smith     if (idx[i] < 0) continue;
1909566063dSJacob Faibussowitsch     PetscCall(PetscHMapISet(map->globalht, idx[i], i));
191413f72f0SBarry Smith   }
192413f72f0SBarry Smith   mapping->data = (void *)map;
193413f72f0SBarry Smith   PetscFunctionReturn(0);
194413f72f0SBarry Smith }
195413f72f0SBarry Smith 
196d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISLocalToGlobalMappingDestroy_Basic(ISLocalToGlobalMapping mapping)
197d71ae5a4SJacob Faibussowitsch {
198413f72f0SBarry Smith   ISLocalToGlobalMapping_Basic *map = (ISLocalToGlobalMapping_Basic *)mapping->data;
199413f72f0SBarry Smith 
200413f72f0SBarry Smith   PetscFunctionBegin;
201413f72f0SBarry Smith   if (!map) PetscFunctionReturn(0);
2029566063dSJacob Faibussowitsch   PetscCall(PetscFree(map->globals));
2039566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping->data));
204413f72f0SBarry Smith   PetscFunctionReturn(0);
205413f72f0SBarry Smith }
206413f72f0SBarry Smith 
207d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISLocalToGlobalMappingDestroy_Hash(ISLocalToGlobalMapping mapping)
208d71ae5a4SJacob Faibussowitsch {
209413f72f0SBarry Smith   ISLocalToGlobalMapping_Hash *map = (ISLocalToGlobalMapping_Hash *)mapping->data;
210413f72f0SBarry Smith 
211413f72f0SBarry Smith   PetscFunctionBegin;
212413f72f0SBarry Smith   if (!map) PetscFunctionReturn(0);
2139566063dSJacob Faibussowitsch   PetscCall(PetscHMapIDestroy(&map->globalht));
2149566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping->data));
215413f72f0SBarry Smith   PetscFunctionReturn(0);
216413f72f0SBarry Smith }
217413f72f0SBarry Smith 
218413f72f0SBarry Smith #define GTOLTYPE _Basic
219413f72f0SBarry Smith #define GTOLNAME _Basic
220541bf97eSAdrian Croucher #define GTOLBS   mapping->bs
2219371c9d4SSatish Balay #define GTOL(g, local) \
2229371c9d4SSatish Balay   do { \
223413f72f0SBarry Smith     local = map->globals[g / bs - start]; \
2240040bde1SJunchao Zhang     if (local >= 0) local = bs * local + (g % bs); \
225413f72f0SBarry Smith   } while (0)
226541bf97eSAdrian Croucher 
227413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
228413f72f0SBarry Smith 
229413f72f0SBarry Smith #define GTOLTYPE _Basic
230413f72f0SBarry Smith #define GTOLNAME Block_Basic
231541bf97eSAdrian Croucher #define GTOLBS   1
2329371c9d4SSatish Balay #define GTOL(g, local) \
233d71ae5a4SJacob Faibussowitsch   do { \
234d71ae5a4SJacob Faibussowitsch     local = map->globals[g - start]; \
235d71ae5a4SJacob Faibussowitsch   } while (0)
236413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
237413f72f0SBarry Smith 
238413f72f0SBarry Smith #define GTOLTYPE _Hash
239413f72f0SBarry Smith #define GTOLNAME _Hash
240541bf97eSAdrian Croucher #define GTOLBS   mapping->bs
2419371c9d4SSatish Balay #define GTOL(g, local) \
2429371c9d4SSatish Balay   do { \
243e8f14785SLisandro Dalcin     (void)PetscHMapIGet(map->globalht, g / bs, &local); \
2440040bde1SJunchao Zhang     if (local >= 0) local = bs * local + (g % bs); \
245413f72f0SBarry Smith   } while (0)
246413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
247413f72f0SBarry Smith 
248413f72f0SBarry Smith #define GTOLTYPE _Hash
249413f72f0SBarry Smith #define GTOLNAME Block_Hash
250541bf97eSAdrian Croucher #define GTOLBS   1
2519371c9d4SSatish Balay #define GTOL(g, local) \
252d71ae5a4SJacob Faibussowitsch   do { \
253d71ae5a4SJacob Faibussowitsch     (void)PetscHMapIGet(map->globalht, g, &local); \
254d71ae5a4SJacob Faibussowitsch   } while (0)
255413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
256413f72f0SBarry Smith 
2576658fb44Sstefano_zampini /*@
2586658fb44Sstefano_zampini     ISLocalToGlobalMappingDuplicate - Duplicates the local to global mapping object
2596658fb44Sstefano_zampini 
2606658fb44Sstefano_zampini     Not Collective
2616658fb44Sstefano_zampini 
2626658fb44Sstefano_zampini     Input Parameter:
2636658fb44Sstefano_zampini .   ltog - local to global mapping
2646658fb44Sstefano_zampini 
2656658fb44Sstefano_zampini     Output Parameter:
2666658fb44Sstefano_zampini .   nltog - the duplicated local to global mapping
2676658fb44Sstefano_zampini 
2686658fb44Sstefano_zampini     Level: advanced
2696658fb44Sstefano_zampini 
270cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`
2716658fb44Sstefano_zampini @*/
272d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingDuplicate(ISLocalToGlobalMapping ltog, ISLocalToGlobalMapping *nltog)
273d71ae5a4SJacob Faibussowitsch {
274a0d79125SStefano Zampini   ISLocalToGlobalMappingType l2gtype;
2756658fb44Sstefano_zampini 
2766658fb44Sstefano_zampini   PetscFunctionBegin;
2776658fb44Sstefano_zampini   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
2789566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)ltog), ltog->bs, ltog->n, ltog->indices, PETSC_COPY_VALUES, nltog));
2799566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetType(ltog, &l2gtype));
2809566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingSetType(*nltog, l2gtype));
2816658fb44Sstefano_zampini   PetscFunctionReturn(0);
2826658fb44Sstefano_zampini }
2836658fb44Sstefano_zampini 
284565245c5SBarry Smith /*@
285107e9a97SBarry Smith     ISLocalToGlobalMappingGetSize - Gets the local size of a local to global mapping
2863b9aefa3SBarry Smith 
2873b9aefa3SBarry Smith     Not Collective
2883b9aefa3SBarry Smith 
2893b9aefa3SBarry Smith     Input Parameter:
2903b9aefa3SBarry Smith .   ltog - local to global mapping
2913b9aefa3SBarry Smith 
2923b9aefa3SBarry Smith     Output Parameter:
293cab54364SBarry Smith .   n - the number of entries in the local mapping, `ISLocalToGlobalMappingGetIndices()` returns an array of this length
2943b9aefa3SBarry Smith 
2953b9aefa3SBarry Smith     Level: advanced
2963b9aefa3SBarry Smith 
297cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`
2983b9aefa3SBarry Smith @*/
299d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping mapping, PetscInt *n)
300d71ae5a4SJacob Faibussowitsch {
3013b9aefa3SBarry Smith   PetscFunctionBegin;
3020700a824SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
3034482741eSBarry Smith   PetscValidIntPointer(n, 2);
304107e9a97SBarry Smith   *n = mapping->bs * mapping->n;
3053b9aefa3SBarry Smith   PetscFunctionReturn(0);
3063b9aefa3SBarry Smith }
3073b9aefa3SBarry Smith 
3085a5d4f66SBarry Smith /*@C
309cab54364SBarry Smith    ISLocalToGlobalMappingViewFromOptions - View an `ISLocalToGlobalMapping` based on values in the options database
310fe2efc57SMark 
311c3339decSBarry Smith    Collective
312fe2efc57SMark 
313fe2efc57SMark    Input Parameters:
314fe2efc57SMark +  A - the local to global mapping object
315736c3998SJose E. Roman .  obj - Optional object
316736c3998SJose E. Roman -  name - command line option
317fe2efc57SMark 
318fe2efc57SMark    Level: intermediate
319cab54364SBarry Smith 
320cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingView`, `PetscObjectViewFromOptions()`, `ISLocalToGlobalMappingCreate()`
321fe2efc57SMark @*/
322d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingViewFromOptions(ISLocalToGlobalMapping A, PetscObject obj, const char name[])
323d71ae5a4SJacob Faibussowitsch {
324fe2efc57SMark   PetscFunctionBegin;
325fe2efc57SMark   PetscValidHeaderSpecific(A, IS_LTOGM_CLASSID, 1);
3269566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
327fe2efc57SMark   PetscFunctionReturn(0);
328fe2efc57SMark }
329fe2efc57SMark 
330fe2efc57SMark /*@C
3315a5d4f66SBarry Smith     ISLocalToGlobalMappingView - View a local to global mapping
3325a5d4f66SBarry Smith 
333b9cd556bSLois Curfman McInnes     Not Collective
334b9cd556bSLois Curfman McInnes 
3355a5d4f66SBarry Smith     Input Parameters:
3363b9aefa3SBarry Smith +   ltog - local to global mapping
3373b9aefa3SBarry Smith -   viewer - viewer
3385a5d4f66SBarry Smith 
339a997ad1aSLois Curfman McInnes     Level: advanced
340a997ad1aSLois Curfman McInnes 
341cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`
3425a5d4f66SBarry Smith @*/
343d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingView(ISLocalToGlobalMapping mapping, PetscViewer viewer)
344d71ae5a4SJacob Faibussowitsch {
34532dcc486SBarry Smith   PetscInt    i;
34632dcc486SBarry Smith   PetscMPIInt rank;
347ace3abfcSBarry Smith   PetscBool   iascii;
3485a5d4f66SBarry Smith 
3495a5d4f66SBarry Smith   PetscFunctionBegin;
3500700a824SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
35148a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mapping), &viewer));
3520700a824SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
3535a5d4f66SBarry Smith 
3549566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mapping), &rank));
3559566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
35632077d6dSBarry Smith   if (iascii) {
3579566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)mapping, viewer));
3589566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
359*f2c6b1a2SJed Brown     for (i = 0; i < mapping->n; i++) {
360*f2c6b1a2SJed Brown       PetscInt bs = mapping->bs, g = mapping->indices[i];
361*f2c6b1a2SJed Brown       if (bs == 1) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %" PetscInt_FMT " %" PetscInt_FMT "\n", rank, i, g));
362*f2c6b1a2SJed 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));
363*f2c6b1a2SJed Brown     }
3649566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
3659566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
3661575c14dSBarry Smith   }
3675a5d4f66SBarry Smith   PetscFunctionReturn(0);
3685a5d4f66SBarry Smith }
3695a5d4f66SBarry Smith 
3701f428162SBarry Smith /*@
3712bdab257SBarry Smith     ISLocalToGlobalMappingCreateIS - Creates a mapping between a local (0 to n)
3722bdab257SBarry Smith     ordering and a global parallel ordering.
3732bdab257SBarry Smith 
3740f5bd95cSBarry Smith     Not collective
375b9cd556bSLois Curfman McInnes 
376a997ad1aSLois Curfman McInnes     Input Parameter:
3778c03b21aSDmitry Karpeev .   is - index set containing the global numbers for each local number
3782bdab257SBarry Smith 
379a997ad1aSLois Curfman McInnes     Output Parameter:
3802bdab257SBarry Smith .   mapping - new mapping data structure
3812bdab257SBarry Smith 
382a997ad1aSLois Curfman McInnes     Level: advanced
383a997ad1aSLois Curfman McInnes 
384cab54364SBarry Smith     Note:
385cab54364SBarry Smith     the block size of the `IS` determines the block size of the mapping
386cab54364SBarry Smith 
387cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetFromOptions()`
3882bdab257SBarry Smith @*/
389d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingCreateIS(IS is, ISLocalToGlobalMapping *mapping)
390d71ae5a4SJacob Faibussowitsch {
3913bbf0e92SBarry Smith   PetscInt        n, bs;
3925d0c19d7SBarry Smith   const PetscInt *indices;
3932bdab257SBarry Smith   MPI_Comm        comm;
3943bbf0e92SBarry Smith   PetscBool       isblock;
3953a40ed3dSBarry Smith 
3963a40ed3dSBarry Smith   PetscFunctionBegin;
3970700a824SBarry Smith   PetscValidHeaderSpecific(is, IS_CLASSID, 1);
3984482741eSBarry Smith   PetscValidPointer(mapping, 2);
3992bdab257SBarry Smith 
4009566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)is, &comm));
4019566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is, &n));
4029566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)is, ISBLOCK, &isblock));
4036006e8d2SBarry Smith   if (!isblock) {
4049566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(is, &indices));
4059566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(comm, 1, n, indices, PETSC_COPY_VALUES, mapping));
4069566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(is, &indices));
4076006e8d2SBarry Smith   } else {
4089566063dSJacob Faibussowitsch     PetscCall(ISGetBlockSize(is, &bs));
4099566063dSJacob Faibussowitsch     PetscCall(ISBlockGetIndices(is, &indices));
4109566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(comm, bs, n / bs, indices, PETSC_COPY_VALUES, mapping));
4119566063dSJacob Faibussowitsch     PetscCall(ISBlockRestoreIndices(is, &indices));
4126006e8d2SBarry Smith   }
4133a40ed3dSBarry Smith   PetscFunctionReturn(0);
4142bdab257SBarry Smith }
4155a5d4f66SBarry Smith 
416a4d96a55SJed Brown /*@C
417a4d96a55SJed Brown     ISLocalToGlobalMappingCreateSF - Creates a mapping between a local (0 to n)
418a4d96a55SJed Brown     ordering and a global parallel ordering.
419a4d96a55SJed Brown 
420a4d96a55SJed Brown     Collective
421a4d96a55SJed Brown 
422d8d19677SJose E. Roman     Input Parameters:
423a4d96a55SJed Brown +   sf - star forest mapping contiguous local indices to (rank, offset)
424cab54364SBarry Smith -   start - first global index on this process, or `PETSC_DECIDE` to compute contiguous global numbering automatically
425a4d96a55SJed Brown 
426a4d96a55SJed Brown     Output Parameter:
427a4d96a55SJed Brown .   mapping - new mapping data structure
428a4d96a55SJed Brown 
429a4d96a55SJed Brown     Level: advanced
430a4d96a55SJed Brown 
4319a535bafSVaclav Hapla     Notes:
432cab54364SBarry Smith     If any processor calls this with start = `PETSC_DECIDE` then all processors must, otherwise the program will hang.
4339a535bafSVaclav Hapla 
434cab54364SBarry Smith .seealso: [](sec_scatter), `PetscSF`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingSetFromOptions()`
435a4d96a55SJed Brown @*/
436d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF sf, PetscInt start, ISLocalToGlobalMapping *mapping)
437d71ae5a4SJacob Faibussowitsch {
438a4d96a55SJed Brown   PetscInt i, maxlocal, nroots, nleaves, *globals, *ltog;
439a4d96a55SJed Brown   MPI_Comm comm;
440a4d96a55SJed Brown 
441a4d96a55SJed Brown   PetscFunctionBegin;
442a4d96a55SJed Brown   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
443a4d96a55SJed Brown   PetscValidPointer(mapping, 3);
4449566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
44541f4c31fSVaclav Hapla   PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL));
4469a535bafSVaclav Hapla   if (start == PETSC_DECIDE) {
4479a535bafSVaclav Hapla     start = 0;
4489566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Exscan(&nroots, &start, 1, MPIU_INT, MPI_SUM, comm));
44941f4c31fSVaclav Hapla   } else PetscCheck(start >= 0, comm, PETSC_ERR_ARG_OUTOFRANGE, "start must be nonnegative or PETSC_DECIDE");
45041f4c31fSVaclav Hapla   PetscCall(PetscSFGetLeafRange(sf, NULL, &maxlocal));
45141f4c31fSVaclav Hapla   ++maxlocal;
4529566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nroots, &globals));
4539566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxlocal, &ltog));
454a4d96a55SJed Brown   for (i = 0; i < nroots; i++) globals[i] = start + i;
455a4d96a55SJed Brown   for (i = 0; i < maxlocal; i++) ltog[i] = -1;
4569566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sf, MPIU_INT, globals, ltog, MPI_REPLACE));
4579566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sf, MPIU_INT, globals, ltog, MPI_REPLACE));
4589566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(comm, 1, maxlocal, ltog, PETSC_OWN_POINTER, mapping));
4599566063dSJacob Faibussowitsch   PetscCall(PetscFree(globals));
460a4d96a55SJed Brown   PetscFunctionReturn(0);
461a4d96a55SJed Brown }
462b46b645bSBarry Smith 
46363fa5c83Sstefano_zampini /*@
46463fa5c83Sstefano_zampini     ISLocalToGlobalMappingSetBlockSize - Sets the blocksize of the mapping
46563fa5c83Sstefano_zampini 
46663fa5c83Sstefano_zampini     Not collective
46763fa5c83Sstefano_zampini 
46863fa5c83Sstefano_zampini     Input Parameters:
469a2b725a8SWilliam Gropp +   mapping - mapping data structure
470a2b725a8SWilliam Gropp -   bs - the blocksize
47163fa5c83Sstefano_zampini 
47263fa5c83Sstefano_zampini     Level: advanced
47363fa5c83Sstefano_zampini 
474cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`
47563fa5c83Sstefano_zampini @*/
476d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingSetBlockSize(ISLocalToGlobalMapping mapping, PetscInt bs)
477d71ae5a4SJacob Faibussowitsch {
478a59f3c4dSstefano_zampini   PetscInt       *nid;
479a59f3c4dSstefano_zampini   const PetscInt *oid;
480a59f3c4dSstefano_zampini   PetscInt        i, cn, on, obs, nn;
48163fa5c83Sstefano_zampini 
48263fa5c83Sstefano_zampini   PetscFunctionBegin;
48363fa5c83Sstefano_zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
48408401ef6SPierre Jolivet   PetscCheck(bs >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid block size %" PetscInt_FMT, bs);
48563fa5c83Sstefano_zampini   if (bs == mapping->bs) PetscFunctionReturn(0);
48663fa5c83Sstefano_zampini   on  = mapping->n;
48763fa5c83Sstefano_zampini   obs = mapping->bs;
48863fa5c83Sstefano_zampini   oid = mapping->indices;
48963fa5c83Sstefano_zampini   nn  = (on * obs) / bs;
49008401ef6SPierre 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);
491a59f3c4dSstefano_zampini 
4929566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nn, &nid));
4939566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetIndices(mapping, &oid));
494a59f3c4dSstefano_zampini   for (i = 0; i < nn; i++) {
495a59f3c4dSstefano_zampini     PetscInt j;
496a59f3c4dSstefano_zampini     for (j = 0, cn = 0; j < bs - 1; j++) {
4979371c9d4SSatish Balay       if (oid[i * bs + j] < 0) {
4989371c9d4SSatish Balay         cn++;
4999371c9d4SSatish Balay         continue;
5009371c9d4SSatish Balay       }
50108401ef6SPierre 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]);
502a59f3c4dSstefano_zampini     }
503a59f3c4dSstefano_zampini     if (oid[i * bs + j] < 0) cn++;
5048b7cb0e6Sstefano_zampini     if (cn) {
50508401ef6SPierre 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);
506a59f3c4dSstefano_zampini       nid[i] = -1;
5078b7cb0e6Sstefano_zampini     } else {
508a59f3c4dSstefano_zampini       nid[i] = oid[i * bs] / bs;
50963fa5c83Sstefano_zampini     }
51063fa5c83Sstefano_zampini   }
5119566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreIndices(mapping, &oid));
512a59f3c4dSstefano_zampini 
51363fa5c83Sstefano_zampini   mapping->n  = nn;
51463fa5c83Sstefano_zampini   mapping->bs = bs;
5159566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping->indices));
51663fa5c83Sstefano_zampini   mapping->indices     = nid;
517c9345713Sstefano_zampini   mapping->globalstart = 0;
518c9345713Sstefano_zampini   mapping->globalend   = 0;
5191bd0b88eSStefano Zampini 
5201bd0b88eSStefano Zampini   /* reset the cached information */
5219566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping->info_procs));
5229566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping->info_numprocs));
5231bd0b88eSStefano Zampini   if (mapping->info_indices) {
5241bd0b88eSStefano Zampini     PetscInt i;
5251bd0b88eSStefano Zampini 
5269566063dSJacob Faibussowitsch     PetscCall(PetscFree((mapping->info_indices)[0]));
52748a46eb9SPierre Jolivet     for (i = 1; i < mapping->info_nproc; i++) PetscCall(PetscFree(mapping->info_indices[i]));
5289566063dSJacob Faibussowitsch     PetscCall(PetscFree(mapping->info_indices));
5291bd0b88eSStefano Zampini   }
5301bd0b88eSStefano Zampini   mapping->info_cached = PETSC_FALSE;
5311bd0b88eSStefano Zampini 
532dbbe0bcdSBarry Smith   PetscTryTypeMethod(mapping, destroy);
53363fa5c83Sstefano_zampini   PetscFunctionReturn(0);
53463fa5c83Sstefano_zampini }
53563fa5c83Sstefano_zampini 
53645b6f7e9SBarry Smith /*@
53745b6f7e9SBarry Smith     ISLocalToGlobalMappingGetBlockSize - Gets the blocksize of the mapping
53845b6f7e9SBarry Smith     ordering and a global parallel ordering.
53945b6f7e9SBarry Smith 
54045b6f7e9SBarry Smith     Not Collective
54145b6f7e9SBarry Smith 
54245b6f7e9SBarry Smith     Input Parameters:
54345b6f7e9SBarry Smith .   mapping - mapping data structure
54445b6f7e9SBarry Smith 
54545b6f7e9SBarry Smith     Output Parameter:
54645b6f7e9SBarry Smith .   bs - the blocksize
54745b6f7e9SBarry Smith 
54845b6f7e9SBarry Smith     Level: advanced
54945b6f7e9SBarry Smith 
550cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`
55145b6f7e9SBarry Smith @*/
552d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping mapping, PetscInt *bs)
553d71ae5a4SJacob Faibussowitsch {
55445b6f7e9SBarry Smith   PetscFunctionBegin;
555cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
55645b6f7e9SBarry Smith   *bs = mapping->bs;
55745b6f7e9SBarry Smith   PetscFunctionReturn(0);
55845b6f7e9SBarry Smith }
55945b6f7e9SBarry Smith 
560ba5bb76aSSatish Balay /*@
56190f02eecSBarry Smith     ISLocalToGlobalMappingCreate - Creates a mapping between a local (0 to n)
56290f02eecSBarry Smith     ordering and a global parallel ordering.
5632362add9SBarry Smith 
56489d82c54SBarry Smith     Not Collective, but communicator may have more than one process
565b9cd556bSLois Curfman McInnes 
5662362add9SBarry Smith     Input Parameters:
56789d82c54SBarry Smith +   comm - MPI communicator
568f0413b6fSBarry Smith .   bs - the block size
56928bc9809SBarry Smith .   n - the number of local elements divided by the block size, or equivalently the number of block indices
57028bc9809SBarry 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
571d5ad8652SBarry Smith -   mode - see PetscCopyMode
5722362add9SBarry Smith 
573a997ad1aSLois Curfman McInnes     Output Parameter:
57490f02eecSBarry Smith .   mapping - new mapping data structure
5752362add9SBarry Smith 
576cab54364SBarry Smith     Level: advanced
577cab54364SBarry Smith 
57895452b02SPatrick Sanan     Notes:
57995452b02SPatrick Sanan     There is one integer value in indices per block and it represents the actual indices bs*idx + j, where j=0,..,bs-1
580413f72f0SBarry Smith 
581cab54364SBarry Smith     For "small" problems when using `ISGlobalToLocalMappingApply()` and `ISGlobalToLocalMappingApplyBlock()`, the `ISLocalToGlobalMappingType`
582cab54364SBarry Smith     of `ISLOCALTOGLOBALMAPPINGBASIC` will be used; this uses more memory but is faster; this approach is not scalable for extremely large mappings.
583413f72f0SBarry Smith 
584cab54364SBarry Smith     For large problems `ISLOCALTOGLOBALMAPPINGHASH` is used, this is scalable.
585cab54364SBarry Smith     Use `ISLocalToGlobalMappingSetType()` or call `ISLocalToGlobalMappingSetFromOptions()` with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
586a997ad1aSLois Curfman McInnes 
587cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingSetFromOptions()`, `ISLOCALTOGLOBALMAPPINGBASIC`, `ISLOCALTOGLOBALMAPPINGHASH`
588db781477SPatrick Sanan           `ISLocalToGlobalMappingSetType()`, `ISLocalToGlobalMappingType`
5892362add9SBarry Smith @*/
590d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingCreate(MPI_Comm comm, PetscInt bs, PetscInt n, const PetscInt indices[], PetscCopyMode mode, ISLocalToGlobalMapping *mapping)
591d71ae5a4SJacob Faibussowitsch {
59232dcc486SBarry Smith   PetscInt *in;
593b46b645bSBarry Smith 
594b46b645bSBarry Smith   PetscFunctionBegin;
595064a246eSJacob Faibussowitsch   if (n) PetscValidIntPointer(indices, 4);
596064a246eSJacob Faibussowitsch   PetscValidPointer(mapping, 6);
597b46b645bSBarry Smith 
5980298fd71SBarry Smith   *mapping = NULL;
5999566063dSJacob Faibussowitsch   PetscCall(ISInitializePackage());
6002362add9SBarry Smith 
6019566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(*mapping, IS_LTOGM_CLASSID, "ISLocalToGlobalMapping", "Local to global mapping", "IS", comm, ISLocalToGlobalMappingDestroy, ISLocalToGlobalMappingView));
602d4bb536fSBarry Smith   (*mapping)->n  = n;
603f0413b6fSBarry Smith   (*mapping)->bs = bs;
604d5ad8652SBarry Smith   if (mode == PETSC_COPY_VALUES) {
6059566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &in));
6069566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(in, indices, n));
607d5ad8652SBarry Smith     (*mapping)->indices         = in;
60871910c26SVaclav Hapla     (*mapping)->dealloc_indices = PETSC_TRUE;
6096389a1a1SBarry Smith   } else if (mode == PETSC_OWN_POINTER) {
6106389a1a1SBarry Smith     (*mapping)->indices         = (PetscInt *)indices;
61171910c26SVaclav Hapla     (*mapping)->dealloc_indices = PETSC_TRUE;
61271910c26SVaclav Hapla   } else if (mode == PETSC_USE_POINTER) {
61371910c26SVaclav Hapla     (*mapping)->indices = (PetscInt *)indices;
6149371c9d4SSatish Balay   } else SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid mode %d", mode);
6153a40ed3dSBarry Smith   PetscFunctionReturn(0);
6162362add9SBarry Smith }
6172362add9SBarry Smith 
618413f72f0SBarry Smith PetscFunctionList ISLocalToGlobalMappingList = NULL;
619413f72f0SBarry Smith 
62090f02eecSBarry Smith /*@
6217e99dc12SLawrence Mitchell    ISLocalToGlobalMappingSetFromOptions - Set mapping options from the options database.
6227e99dc12SLawrence Mitchell 
6237e99dc12SLawrence Mitchell    Not collective
6247e99dc12SLawrence Mitchell 
6257e99dc12SLawrence Mitchell    Input Parameters:
6267e99dc12SLawrence Mitchell .  mapping - mapping data structure
6277e99dc12SLawrence Mitchell 
6287e99dc12SLawrence Mitchell    Level: advanced
6297e99dc12SLawrence Mitchell 
630cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingSetFromOptions()`, `ISLOCALTOGLOBALMAPPINGBASIC`, `ISLOCALTOGLOBALMAPPINGHASH`
631cab54364SBarry Smith           `ISLocalToGlobalMappingSetType()`, `ISLocalToGlobalMappingType`
6327e99dc12SLawrence Mitchell @*/
633d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingSetFromOptions(ISLocalToGlobalMapping mapping)
634d71ae5a4SJacob Faibussowitsch {
635413f72f0SBarry Smith   char                       type[256];
636413f72f0SBarry Smith   ISLocalToGlobalMappingType defaulttype = "Not set";
6377e99dc12SLawrence Mitchell   PetscBool                  flg;
6387e99dc12SLawrence Mitchell 
6397e99dc12SLawrence Mitchell   PetscFunctionBegin;
6407e99dc12SLawrence Mitchell   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
6419566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRegisterAll());
642d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)mapping);
6439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-islocaltoglobalmapping_type", "ISLocalToGlobalMapping method", "ISLocalToGlobalMappingSetType", ISLocalToGlobalMappingList, (char *)(((PetscObject)mapping)->type_name) ? ((PetscObject)mapping)->type_name : defaulttype, type, 256, &flg));
6441baa6e33SBarry Smith   if (flg) PetscCall(ISLocalToGlobalMappingSetType(mapping, type));
645d0609cedSBarry Smith   PetscOptionsEnd();
6467e99dc12SLawrence Mitchell   PetscFunctionReturn(0);
6477e99dc12SLawrence Mitchell }
6487e99dc12SLawrence Mitchell 
6497e99dc12SLawrence Mitchell /*@
65090f02eecSBarry Smith    ISLocalToGlobalMappingDestroy - Destroys a mapping between a local (0 to n)
65190f02eecSBarry Smith    ordering and a global parallel ordering.
65290f02eecSBarry Smith 
6530f5bd95cSBarry Smith    Note Collective
654b9cd556bSLois Curfman McInnes 
65590f02eecSBarry Smith    Input Parameters:
65690f02eecSBarry Smith .  mapping - mapping data structure
65790f02eecSBarry Smith 
658a997ad1aSLois Curfman McInnes    Level: advanced
659a997ad1aSLois Curfman McInnes 
660cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingCreate()`
66190f02eecSBarry Smith @*/
662d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping *mapping)
663d71ae5a4SJacob Faibussowitsch {
6643a40ed3dSBarry Smith   PetscFunctionBegin;
6656bf464f9SBarry Smith   if (!*mapping) PetscFunctionReturn(0);
6666bf464f9SBarry Smith   PetscValidHeaderSpecific((*mapping), IS_LTOGM_CLASSID, 1);
6679371c9d4SSatish Balay   if (--((PetscObject)(*mapping))->refct > 0) {
6689371c9d4SSatish Balay     *mapping = NULL;
6699371c9d4SSatish Balay     PetscFunctionReturn(0);
67071910c26SVaclav Hapla   }
67148a46eb9SPierre Jolivet   if ((*mapping)->dealloc_indices) PetscCall(PetscFree((*mapping)->indices));
6729566063dSJacob Faibussowitsch   PetscCall(PetscFree((*mapping)->info_procs));
6739566063dSJacob Faibussowitsch   PetscCall(PetscFree((*mapping)->info_numprocs));
674268a049cSStefano Zampini   if ((*mapping)->info_indices) {
675268a049cSStefano Zampini     PetscInt i;
676268a049cSStefano Zampini 
6779566063dSJacob Faibussowitsch     PetscCall(PetscFree(((*mapping)->info_indices)[0]));
67848a46eb9SPierre Jolivet     for (i = 1; i < (*mapping)->info_nproc; i++) PetscCall(PetscFree(((*mapping)->info_indices)[i]));
6799566063dSJacob Faibussowitsch     PetscCall(PetscFree((*mapping)->info_indices));
680268a049cSStefano Zampini   }
68148a46eb9SPierre Jolivet   if ((*mapping)->info_nodei) PetscCall(PetscFree(((*mapping)->info_nodei)[0]));
6829566063dSJacob Faibussowitsch   PetscCall(PetscFree2((*mapping)->info_nodec, (*mapping)->info_nodei));
68348a46eb9SPierre Jolivet   if ((*mapping)->ops->destroy) PetscCall((*(*mapping)->ops->destroy)(*mapping));
6849566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(mapping));
6854c8fdceaSLisandro Dalcin   *mapping = NULL;
6863a40ed3dSBarry Smith   PetscFunctionReturn(0);
68790f02eecSBarry Smith }
68890f02eecSBarry Smith 
68990f02eecSBarry Smith /*@
690cab54364SBarry Smith     ISLocalToGlobalMappingApplyIS - Creates from an `IS` in the local numbering
691cab54364SBarry Smith     a new index set using the global numbering defined in an `ISLocalToGlobalMapping`
6923acfe500SLois Curfman McInnes     context.
69390f02eecSBarry Smith 
694c3339decSBarry Smith     Collective
695b9cd556bSLois Curfman McInnes 
69690f02eecSBarry Smith     Input Parameters:
697b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
698b9cd556bSLois Curfman McInnes -   is - index set in local numbering
69990f02eecSBarry Smith 
700cab54364SBarry Smith     Output Parameter:
70190f02eecSBarry Smith .   newis - index set in global numbering
70290f02eecSBarry Smith 
703a997ad1aSLois Curfman McInnes     Level: advanced
704a997ad1aSLois Curfman McInnes 
705cab54364SBarry Smith     Note:
706cab54364SBarry Smith     The output `IS` will have the same communicator of the input `IS`.
707cab54364SBarry Smith 
708cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingCreate()`,
709db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`, `ISGlobalToLocalMappingApply()`
71090f02eecSBarry Smith @*/
711d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping, IS is, IS *newis)
712d71ae5a4SJacob Faibussowitsch {
713e24637baSBarry Smith   PetscInt        n, *idxout;
7145d0c19d7SBarry Smith   const PetscInt *idxin;
7153a40ed3dSBarry Smith 
7163a40ed3dSBarry Smith   PetscFunctionBegin;
7170700a824SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
7180700a824SBarry Smith   PetscValidHeaderSpecific(is, IS_CLASSID, 2);
7194482741eSBarry Smith   PetscValidPointer(newis, 3);
72090f02eecSBarry Smith 
7219566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is, &n));
7229566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is, &idxin));
7239566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &idxout));
7249566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(mapping, n, idxin, idxout));
7259566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is, &idxin));
7269566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)is), n, idxout, PETSC_OWN_POINTER, newis));
7273a40ed3dSBarry Smith   PetscFunctionReturn(0);
72890f02eecSBarry Smith }
72990f02eecSBarry Smith 
730b89cb25eSSatish Balay /*@
7313acfe500SLois Curfman McInnes    ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering
7323acfe500SLois Curfman McInnes    and converts them to the global numbering.
73390f02eecSBarry Smith 
734b9cd556bSLois Curfman McInnes    Not collective
735b9cd556bSLois Curfman McInnes 
736bb25748dSBarry Smith    Input Parameters:
737b9cd556bSLois Curfman McInnes +  mapping - the local to global mapping context
738bb25748dSBarry Smith .  N - number of integers
739b9cd556bSLois Curfman McInnes -  in - input indices in local numbering
740bb25748dSBarry Smith 
741bb25748dSBarry Smith    Output Parameter:
742bb25748dSBarry Smith .  out - indices in global numbering
743bb25748dSBarry Smith 
744a997ad1aSLois Curfman McInnes    Level: advanced
745a997ad1aSLois Curfman McInnes 
746cab54364SBarry Smith    Note:
747cab54364SBarry Smith    The in and out array parameters may be identical.
748cab54364SBarry Smith 
749cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApplyBlock()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingDestroy()`,
750c2e3fba1SPatrick Sanan           `ISLocalToGlobalMappingApplyIS()`, `AOCreateBasic()`, `AOApplicationToPetsc()`,
751db781477SPatrick Sanan           `AOPetscToApplication()`, `ISGlobalToLocalMappingApply()`
752afcb2eb5SJed Brown @*/
753d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping, PetscInt N, const PetscInt in[], PetscInt out[])
754d71ae5a4SJacob Faibussowitsch {
755cbc1caf0SMatthew G. Knepley   PetscInt i, bs, Nmax;
75645b6f7e9SBarry Smith 
75745b6f7e9SBarry Smith   PetscFunctionBegin;
758cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
759cbc1caf0SMatthew G. Knepley   bs   = mapping->bs;
760cbc1caf0SMatthew G. Knepley   Nmax = bs * mapping->n;
76145b6f7e9SBarry Smith   if (bs == 1) {
762cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
76345b6f7e9SBarry Smith     for (i = 0; i < N; i++) {
76445b6f7e9SBarry Smith       if (in[i] < 0) {
76545b6f7e9SBarry Smith         out[i] = in[i];
76645b6f7e9SBarry Smith         continue;
76745b6f7e9SBarry Smith       }
76808401ef6SPierre 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);
76945b6f7e9SBarry Smith       out[i] = idx[in[i]];
77045b6f7e9SBarry Smith     }
77145b6f7e9SBarry Smith   } else {
772cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
77345b6f7e9SBarry Smith     for (i = 0; i < N; i++) {
77445b6f7e9SBarry Smith       if (in[i] < 0) {
77545b6f7e9SBarry Smith         out[i] = in[i];
77645b6f7e9SBarry Smith         continue;
77745b6f7e9SBarry Smith       }
77808401ef6SPierre 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);
77945b6f7e9SBarry Smith       out[i] = idx[in[i] / bs] * bs + (in[i] % bs);
78045b6f7e9SBarry Smith     }
78145b6f7e9SBarry Smith   }
78245b6f7e9SBarry Smith   PetscFunctionReturn(0);
78345b6f7e9SBarry Smith }
78445b6f7e9SBarry Smith 
78545b6f7e9SBarry Smith /*@
7866006e8d2SBarry Smith    ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local block numbering and converts them to the global block numbering
78745b6f7e9SBarry Smith 
78845b6f7e9SBarry Smith    Not collective
78945b6f7e9SBarry Smith 
79045b6f7e9SBarry Smith    Input Parameters:
79145b6f7e9SBarry Smith +  mapping - the local to global mapping context
79245b6f7e9SBarry Smith .  N - number of integers
7936006e8d2SBarry Smith -  in - input indices in local block numbering
79445b6f7e9SBarry Smith 
79545b6f7e9SBarry Smith    Output Parameter:
7966006e8d2SBarry Smith .  out - indices in global block numbering
79745b6f7e9SBarry Smith 
798cab54364SBarry Smith    Level: advanced
799cab54364SBarry Smith 
800cab54364SBarry Smith    Note:
80145b6f7e9SBarry Smith    The in and out array parameters may be identical.
80245b6f7e9SBarry Smith 
8036006e8d2SBarry Smith    Example:
804cab54364SBarry 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
8056006e8d2SBarry Smith      (the first block) would produce 0 and the mapping applied to 1 (the second block) would produce 3.
8066006e8d2SBarry Smith 
807cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingDestroy()`,
808c2e3fba1SPatrick Sanan           `ISLocalToGlobalMappingApplyIS()`, `AOCreateBasic()`, `AOApplicationToPetsc()`,
809db781477SPatrick Sanan           `AOPetscToApplication()`, `ISGlobalToLocalMappingApply()`
81045b6f7e9SBarry Smith @*/
811d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping, PetscInt N, const PetscInt in[], PetscInt out[])
812d71ae5a4SJacob Faibussowitsch {
8138a1f772fSStefano Zampini   PetscInt        i, Nmax;
8148a1f772fSStefano Zampini   const PetscInt *idx;
815d4bb536fSBarry Smith 
816a0d79125SStefano Zampini   PetscFunctionBegin;
817a0d79125SStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
8188a1f772fSStefano Zampini   Nmax = mapping->n;
8198a1f772fSStefano Zampini   idx  = mapping->indices;
820afcb2eb5SJed Brown   for (i = 0; i < N; i++) {
821afcb2eb5SJed Brown     if (in[i] < 0) {
822afcb2eb5SJed Brown       out[i] = in[i];
823afcb2eb5SJed Brown       continue;
824afcb2eb5SJed Brown     }
82508401ef6SPierre 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);
826afcb2eb5SJed Brown     out[i] = idx[in[i]];
827afcb2eb5SJed Brown   }
828afcb2eb5SJed Brown   PetscFunctionReturn(0);
829afcb2eb5SJed Brown }
830d4bb536fSBarry Smith 
8317e99dc12SLawrence Mitchell /*@
832a997ad1aSLois Curfman McInnes     ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers
833a997ad1aSLois Curfman McInnes     specified with a global numbering.
834d4bb536fSBarry Smith 
835b9cd556bSLois Curfman McInnes     Not collective
836b9cd556bSLois Curfman McInnes 
837d4bb536fSBarry Smith     Input Parameters:
838b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
839cab54364SBarry Smith .   type - `IS_GTOLM_MASK` - maps global indices with no local value to -1 in the output list (i.e., mask them)
840cab54364SBarry Smith            `IS_GTOLM_DROP` - drops the indices with no local value from the output list
841d4bb536fSBarry Smith .   n - number of global indices to map
842b9cd556bSLois Curfman McInnes -   idx - global indices to map
843d4bb536fSBarry Smith 
844d4bb536fSBarry Smith     Output Parameters:
845cab54364SBarry Smith +   nout - number of indices in output array (if type == `IS_GTOLM_MASK` then nout = n)
846b9cd556bSLois Curfman McInnes -   idxout - local index of each global index, one must pass in an array long enough
847cab54364SBarry Smith              to hold all the indices. You can call `ISGlobalToLocalMappingApply()` with
8480298fd71SBarry Smith              idxout == NULL to determine the required length (returned in nout)
849cab54364SBarry Smith              and then allocate the required space and call `ISGlobalToLocalMappingApply()`
850e182c471SBarry Smith              a second time to set the values.
851d4bb536fSBarry Smith 
852cab54364SBarry Smith     Level: advanced
853cab54364SBarry Smith 
854b9cd556bSLois Curfman McInnes     Notes:
8550298fd71SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
856d4bb536fSBarry Smith 
857cab54364SBarry Smith     For "small" problems when using `ISGlobalToLocalMappingApply()` and `ISGlobalToLocalMappingApplyBlock()`, the `ISLocalToGlobalMappingType` of `ISLOCALTOGLOBALMAPPINGBASIC` will be used;
858cab54364SBarry 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.
859cab54364SBarry Smith     Use `ISLocalToGlobalMappingSetType()` or call `ISLocalToGlobalMappingSetFromOptions()` with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
8600f5bd95cSBarry Smith 
861cab54364SBarry Smith     Developer Note:
862cab54364SBarry Smith     The manual page states that idx and idxout may be identical but the calling
86332fd6b96SBarry Smith     sequence declares idx as const so it cannot be the same as idxout.
86432fd6b96SBarry Smith 
865cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApply()`, `ISGlobalToLocalMappingApplyBlock()`, `ISLocalToGlobalMappingCreate()`,
866db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`
867d4bb536fSBarry Smith @*/
868d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, PetscInt n, const PetscInt idx[], PetscInt *nout, PetscInt idxout[])
869d71ae5a4SJacob Faibussowitsch {
8709d90f715SBarry Smith   PetscFunctionBegin;
8719d90f715SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
87248a46eb9SPierre Jolivet   if (!mapping->data) PetscCall(ISGlobalToLocalMappingSetUp(mapping));
873dbbe0bcdSBarry Smith   PetscUseTypeMethod(mapping, globaltolocalmappingapply, type, n, idx, nout, idxout);
8749d90f715SBarry Smith   PetscFunctionReturn(0);
8759d90f715SBarry Smith }
8769d90f715SBarry Smith 
877d4fe737eSStefano Zampini /*@
878cab54364SBarry Smith     ISGlobalToLocalMappingApplyIS - Creates from an `IS` in the global numbering
879cab54364SBarry Smith     a new index set using the local numbering defined in an `ISLocalToGlobalMapping`
880d4fe737eSStefano Zampini     context.
881d4fe737eSStefano Zampini 
882d4fe737eSStefano Zampini     Not collective
883d4fe737eSStefano Zampini 
884d4fe737eSStefano Zampini     Input Parameters:
885d4fe737eSStefano Zampini +   mapping - mapping between local and global numbering
886cab54364SBarry Smith .   type - `IS_GTOLM_MASK` - maps global indices with no local value to -1 in the output list (i.e., mask them)
887cab54364SBarry Smith            `IS_GTOLM_DROP` - drops the indices with no local value from the output list
888d4fe737eSStefano Zampini -   is - index set in global numbering
889d4fe737eSStefano Zampini 
890d4fe737eSStefano Zampini     Output Parameters:
891d4fe737eSStefano Zampini .   newis - index set in local numbering
892d4fe737eSStefano Zampini 
893d4fe737eSStefano Zampini     Level: advanced
894d4fe737eSStefano Zampini 
895cab54364SBarry Smith     Note:
896cab54364SBarry Smith     The output `IS` will be sequential, as it encodes a purely local operation
897cab54364SBarry Smith 
898cab54364SBarry Smith .seealso: [](sec_scatter), `ISGlobalToLocalMapping`, `ISGlobalToLocalMappingApply()`, `ISLocalToGlobalMappingCreate()`,
899db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`
900d4fe737eSStefano Zampini @*/
901d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, IS is, IS *newis)
902d71ae5a4SJacob Faibussowitsch {
903d4fe737eSStefano Zampini   PetscInt        n, nout, *idxout;
904d4fe737eSStefano Zampini   const PetscInt *idxin;
905d4fe737eSStefano Zampini 
906d4fe737eSStefano Zampini   PetscFunctionBegin;
907d4fe737eSStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
908d4fe737eSStefano Zampini   PetscValidHeaderSpecific(is, IS_CLASSID, 3);
909d4fe737eSStefano Zampini   PetscValidPointer(newis, 4);
910d4fe737eSStefano Zampini 
9119566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is, &n));
9129566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is, &idxin));
913d4fe737eSStefano Zampini   if (type == IS_GTOLM_MASK) {
9149566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &idxout));
915d4fe737eSStefano Zampini   } else {
9169566063dSJacob Faibussowitsch     PetscCall(ISGlobalToLocalMappingApply(mapping, type, n, idxin, &nout, NULL));
9179566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nout, &idxout));
918d4fe737eSStefano Zampini   }
9199566063dSJacob Faibussowitsch   PetscCall(ISGlobalToLocalMappingApply(mapping, type, n, idxin, &nout, idxout));
9209566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is, &idxin));
9219566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nout, idxout, PETSC_OWN_POINTER, newis));
922d4fe737eSStefano Zampini   PetscFunctionReturn(0);
923d4fe737eSStefano Zampini }
924d4fe737eSStefano Zampini 
9259d90f715SBarry Smith /*@
9269d90f715SBarry Smith     ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers
9279d90f715SBarry Smith     specified with a block global numbering.
9289d90f715SBarry Smith 
9299d90f715SBarry Smith     Not collective
9309d90f715SBarry Smith 
9319d90f715SBarry Smith     Input Parameters:
9329d90f715SBarry Smith +   mapping - mapping between local and global numbering
933cab54364SBarry Smith .   type - `IS_GTOLM_MASK` - maps global indices with no local value to -1 in the output list (i.e., mask them)
934cab54364SBarry Smith            `IS_GTOLM_DROP` - drops the indices with no local value from the output list
9359d90f715SBarry Smith .   n - number of global indices to map
9369d90f715SBarry Smith -   idx - global indices to map
9379d90f715SBarry Smith 
9389d90f715SBarry Smith     Output Parameters:
939cab54364SBarry Smith +   nout - number of indices in output array (if type == `IS_GTOLM_MASK` then nout = n)
9409d90f715SBarry Smith -   idxout - local index of each global index, one must pass in an array long enough
941cab54364SBarry Smith              to hold all the indices. You can call `ISGlobalToLocalMappingApplyBlock()` with
9429d90f715SBarry Smith              idxout == NULL to determine the required length (returned in nout)
943cab54364SBarry Smith              and then allocate the required space and call `ISGlobalToLocalMappingApplyBlock()`
9449d90f715SBarry Smith              a second time to set the values.
9459d90f715SBarry Smith 
946cab54364SBarry Smith     Level: advanced
947cab54364SBarry Smith 
9489d90f715SBarry Smith     Notes:
9499d90f715SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
9509d90f715SBarry Smith 
951cab54364SBarry Smith     For "small" problems when using `ISGlobalToLocalMappingApply()` and `ISGlobalToLocalMappingApplyBlock()`, the `ISLocalToGlobalMappingType` of `ISLOCALTOGLOBALMAPPINGBASIC` will be used;
952cab54364SBarry 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.
953cab54364SBarry Smith     Use `ISLocalToGlobalMappingSetType()` or call `ISLocalToGlobalMappingSetFromOptions()` with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
9549d90f715SBarry Smith 
955cab54364SBarry Smith     Developer Note:
956cab54364SBarry Smith     The manual page states that idx and idxout may be identical but the calling
9579d90f715SBarry Smith     sequence declares idx as const so it cannot be the same as idxout.
9589d90f715SBarry Smith 
959cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApply()`, `ISGlobalToLocalMappingApply()`, `ISLocalToGlobalMappingCreate()`,
960db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`
9619d90f715SBarry Smith @*/
962d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, PetscInt n, const PetscInt idx[], PetscInt *nout, PetscInt idxout[])
963d71ae5a4SJacob Faibussowitsch {
9643a40ed3dSBarry Smith   PetscFunctionBegin;
9650700a824SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
96648a46eb9SPierre Jolivet   if (!mapping->data) PetscCall(ISGlobalToLocalMappingSetUp(mapping));
967dbbe0bcdSBarry Smith   PetscUseTypeMethod(mapping, globaltolocalmappingapplyblock, type, n, idx, nout, idxout);
9683a40ed3dSBarry Smith   PetscFunctionReturn(0);
969d4bb536fSBarry Smith }
97090f02eecSBarry Smith 
97189d82c54SBarry Smith /*@C
9726a818285SBarry Smith     ISLocalToGlobalMappingGetBlockInfo - Gets the neighbor information for each processor and
97389d82c54SBarry Smith      each index shared by more than one processor
97489d82c54SBarry Smith 
975c3339decSBarry Smith     Collective
97689d82c54SBarry Smith 
977f899ff85SJose E. Roman     Input Parameter:
97889d82c54SBarry Smith .   mapping - the mapping from local to global indexing
97989d82c54SBarry Smith 
980d8d19677SJose E. Roman     Output Parameters:
98189d82c54SBarry Smith +   nproc - number of processors that are connected to this one
98289d82c54SBarry Smith .   proc - neighboring processors
98307b52d57SBarry Smith .   numproc - number of indices for each subdomain (processor)
9843463a7baSJed Brown -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
98589d82c54SBarry Smith 
98689d82c54SBarry Smith     Level: advanced
98789d82c54SBarry Smith 
9882cfcea29SBarry Smith     Fortran Usage:
989cab54364SBarry Smith .vb
9902cfcea29SBarry Smith         PetscInt indices[nproc][numprocmax],ierr)
991cab54364SBarry Smith         ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
992cab54364SBarry Smith         ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
993cab54364SBarry Smith .ve
994cab54364SBarry Smith 
995cab54364SBarry Smith    Fortran Note:
996cab54364SBarry Smith         There is no `ISLocalToGlobalMappingRestoreInfo()` in Fortran. You must make sure that procs[], numprocs[] and
9972cfcea29SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
9982cfcea29SBarry Smith 
999cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1000db781477SPatrick Sanan           `ISLocalToGlobalMappingRestoreInfo()`
100189d82c54SBarry Smith @*/
1002d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[])
1003d71ae5a4SJacob Faibussowitsch {
1004268a049cSStefano Zampini   PetscFunctionBegin;
1005268a049cSStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
1006268a049cSStefano Zampini   if (mapping->info_cached) {
1007268a049cSStefano Zampini     *nproc    = mapping->info_nproc;
1008268a049cSStefano Zampini     *procs    = mapping->info_procs;
1009268a049cSStefano Zampini     *numprocs = mapping->info_numprocs;
1010268a049cSStefano Zampini     *indices  = mapping->info_indices;
1011268a049cSStefano Zampini   } else {
10129566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetBlockInfo_Private(mapping, nproc, procs, numprocs, indices));
1013268a049cSStefano Zampini   }
1014268a049cSStefano Zampini   PetscFunctionReturn(0);
1015268a049cSStefano Zampini }
1016268a049cSStefano Zampini 
1017d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[])
1018d71ae5a4SJacob Faibussowitsch {
101997f1f81fSBarry Smith   PetscMPIInt  size, rank, tag1, tag2, tag3, *len, *source, imdex;
102032dcc486SBarry Smith   PetscInt     i, n = mapping->n, Ng, ng, max = 0, *lindices = mapping->indices;
102132dcc486SBarry Smith   PetscInt    *nprocs, *owner, nsends, *sends, j, *starts, nmax, nrecvs, *recvs, proc;
1022c599c493SJunchao Zhang   PetscInt     cnt, scale, *ownedsenders, *nownedsenders, rstart;
102332dcc486SBarry Smith   PetscInt     node, nownedm, nt, *sends2, nsends2, *starts2, *lens2, *dest, nrecvs2, *starts3, *recvs2, k, *bprocs, *tmp;
102432dcc486SBarry Smith   PetscInt     first_procs, first_numprocs, *first_indices;
102589d82c54SBarry Smith   MPI_Request *recv_waits, *send_waits;
102630dcb7c9SBarry Smith   MPI_Status   recv_status, *send_status, *recv_statuses;
1027ce94432eSBarry Smith   MPI_Comm     comm;
1028ace3abfcSBarry Smith   PetscBool    debug = PETSC_FALSE;
102989d82c54SBarry Smith 
103089d82c54SBarry Smith   PetscFunctionBegin;
10319566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mapping, &comm));
10329566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
10339566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
103424cf384cSBarry Smith   if (size == 1) {
103524cf384cSBarry Smith     *nproc = 0;
10360298fd71SBarry Smith     *procs = NULL;
10379566063dSJacob Faibussowitsch     PetscCall(PetscNew(numprocs));
10381e2105dcSBarry Smith     (*numprocs)[0] = 0;
10399566063dSJacob Faibussowitsch     PetscCall(PetscNew(indices));
10400298fd71SBarry Smith     (*indices)[0] = NULL;
1041268a049cSStefano Zampini     /* save info for reuse */
1042268a049cSStefano Zampini     mapping->info_nproc    = *nproc;
1043268a049cSStefano Zampini     mapping->info_procs    = *procs;
1044268a049cSStefano Zampini     mapping->info_numprocs = *numprocs;
1045268a049cSStefano Zampini     mapping->info_indices  = *indices;
1046268a049cSStefano Zampini     mapping->info_cached   = PETSC_TRUE;
104724cf384cSBarry Smith     PetscFunctionReturn(0);
104824cf384cSBarry Smith   }
104924cf384cSBarry Smith 
10509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)mapping)->options, NULL, "-islocaltoglobalmappinggetinfo_debug", &debug, NULL));
105107b52d57SBarry Smith 
10523677ff5aSBarry Smith   /*
10536a818285SBarry Smith     Notes on ISLocalToGlobalMappingGetBlockInfo
10543677ff5aSBarry Smith 
10553677ff5aSBarry Smith     globally owned node - the nodes that have been assigned to this processor in global
10563677ff5aSBarry Smith            numbering, just for this routine.
10573677ff5aSBarry Smith 
10583677ff5aSBarry Smith     nontrivial globally owned node - node assigned to this processor that is on a subdomain
10593677ff5aSBarry Smith            boundary (i.e. is has more than one local owner)
10603677ff5aSBarry Smith 
10613677ff5aSBarry Smith     locally owned node - node that exists on this processors subdomain
10623677ff5aSBarry Smith 
10633677ff5aSBarry Smith     nontrivial locally owned node - node that is not in the interior (i.e. has more than one
10643677ff5aSBarry Smith            local subdomain
10653677ff5aSBarry Smith   */
10669566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetNewTag((PetscObject)mapping, &tag1));
10679566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetNewTag((PetscObject)mapping, &tag2));
10689566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetNewTag((PetscObject)mapping, &tag3));
106989d82c54SBarry Smith 
107089d82c54SBarry Smith   for (i = 0; i < n; i++) {
107189d82c54SBarry Smith     if (lindices[i] > max) max = lindices[i];
107289d82c54SBarry Smith   }
10731c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&max, &Ng, 1, MPIU_INT, MPI_MAX, comm));
107478058e43SBarry Smith   Ng++;
10759566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
10769566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1077bc8ff85bSBarry Smith   scale = Ng / size + 1;
10789371c9d4SSatish Balay   ng    = scale;
10799371c9d4SSatish Balay   if (rank == size - 1) ng = Ng - scale * (size - 1);
10809371c9d4SSatish Balay   ng     = PetscMax(1, ng);
1081caba0dd0SBarry Smith   rstart = scale * rank;
108289d82c54SBarry Smith 
108389d82c54SBarry Smith   /* determine ownership ranges of global indices */
10849566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * size, &nprocs));
10859566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(nprocs, 2 * size));
108689d82c54SBarry Smith 
108789d82c54SBarry Smith   /* determine owners of each local node  */
10889566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &owner));
108989d82c54SBarry Smith   for (i = 0; i < n; i++) {
10903677ff5aSBarry Smith     proc                 = lindices[i] / scale; /* processor that globally owns this index */
109127c402fcSBarry Smith     nprocs[2 * proc + 1] = 1;                   /* processor globally owns at least one of ours */
10923677ff5aSBarry Smith     owner[i]             = proc;
109327c402fcSBarry Smith     nprocs[2 * proc]++; /* count of how many that processor globally owns of ours */
109489d82c54SBarry Smith   }
10959371c9d4SSatish Balay   nsends = 0;
10969371c9d4SSatish Balay   for (i = 0; i < size; i++) nsends += nprocs[2 * i + 1];
10979566063dSJacob Faibussowitsch   PetscCall(PetscInfo(mapping, "Number of global owners for my local data %" PetscInt_FMT "\n", nsends));
109889d82c54SBarry Smith 
109989d82c54SBarry Smith   /* inform other processors of number of messages and max length*/
11009566063dSJacob Faibussowitsch   PetscCall(PetscMaxSum(comm, nprocs, &nmax, &nrecvs));
11019566063dSJacob Faibussowitsch   PetscCall(PetscInfo(mapping, "Number of local owners for my global data %" PetscInt_FMT "\n", nrecvs));
110289d82c54SBarry Smith 
110389d82c54SBarry Smith   /* post receives for owned rows */
11049566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((2 * nrecvs + 1) * (nmax + 1), &recvs));
11059566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs + 1, &recv_waits));
110648a46eb9SPierre Jolivet   for (i = 0; i < nrecvs; i++) PetscCallMPI(MPI_Irecv(recvs + 2 * nmax * i, 2 * nmax, MPIU_INT, MPI_ANY_SOURCE, tag1, comm, recv_waits + i));
110789d82c54SBarry Smith 
110889d82c54SBarry Smith   /* pack messages containing lists of local nodes to owners */
11099566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * n + 1, &sends));
11109566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size + 1, &starts));
111189d82c54SBarry Smith   starts[0] = 0;
1112f6e5521dSKarl Rupp   for (i = 1; i < size; i++) starts[i] = starts[i - 1] + 2 * nprocs[2 * i - 2];
111389d82c54SBarry Smith   for (i = 0; i < n; i++) {
111489d82c54SBarry Smith     sends[starts[owner[i]]++] = lindices[i];
111530dcb7c9SBarry Smith     sends[starts[owner[i]]++] = i;
111689d82c54SBarry Smith   }
11179566063dSJacob Faibussowitsch   PetscCall(PetscFree(owner));
111889d82c54SBarry Smith   starts[0] = 0;
1119f6e5521dSKarl Rupp   for (i = 1; i < size; i++) starts[i] = starts[i - 1] + 2 * nprocs[2 * i - 2];
112089d82c54SBarry Smith 
112189d82c54SBarry Smith   /* send the messages */
11229566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends + 1, &send_waits));
11239566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends + 1, &dest));
112489d82c54SBarry Smith   cnt = 0;
112589d82c54SBarry Smith   for (i = 0; i < size; i++) {
112627c402fcSBarry Smith     if (nprocs[2 * i]) {
11279566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Isend(sends + starts[i], 2 * nprocs[2 * i], MPIU_INT, i, tag1, comm, send_waits + cnt));
112830dcb7c9SBarry Smith       dest[cnt] = i;
112989d82c54SBarry Smith       cnt++;
113089d82c54SBarry Smith     }
113189d82c54SBarry Smith   }
11329566063dSJacob Faibussowitsch   PetscCall(PetscFree(starts));
113389d82c54SBarry Smith 
113489d82c54SBarry Smith   /* wait on receives */
11359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs + 1, &source));
11369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs + 1, &len));
113789d82c54SBarry Smith   cnt = nrecvs;
11389566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(ng + 1, &nownedsenders));
113989d82c54SBarry Smith   while (cnt) {
11409566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitany(nrecvs, recv_waits, &imdex, &recv_status));
114189d82c54SBarry Smith     /* unpack receives into our local space */
11429566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Get_count(&recv_status, MPIU_INT, &len[imdex]));
114389d82c54SBarry Smith     source[imdex] = recv_status.MPI_SOURCE;
114430dcb7c9SBarry Smith     len[imdex]    = len[imdex] / 2;
1145caba0dd0SBarry Smith     /* count how many local owners for each of my global owned indices */
114630dcb7c9SBarry Smith     for (i = 0; i < len[imdex]; i++) nownedsenders[recvs[2 * imdex * nmax + 2 * i] - rstart]++;
114789d82c54SBarry Smith     cnt--;
114889d82c54SBarry Smith   }
11499566063dSJacob Faibussowitsch   PetscCall(PetscFree(recv_waits));
115089d82c54SBarry Smith 
115130dcb7c9SBarry Smith   /* count how many globally owned indices are on an edge multiplied by how many processors own them. */
1152bc8ff85bSBarry Smith   nownedm = 0;
1153bc8ff85bSBarry Smith   for (i = 0; i < ng; i++) {
1154c599c493SJunchao Zhang     if (nownedsenders[i] > 1) nownedm += nownedsenders[i];
1155bc8ff85bSBarry Smith   }
1156bc8ff85bSBarry Smith 
1157bc8ff85bSBarry Smith   /* create single array to contain rank of all local owners of each globally owned index */
11589566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nownedm + 1, &ownedsenders));
11599566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ng + 1, &starts));
1160bc8ff85bSBarry Smith   starts[0] = 0;
1161bc8ff85bSBarry Smith   for (i = 1; i < ng; i++) {
1162bc8ff85bSBarry Smith     if (nownedsenders[i - 1] > 1) starts[i] = starts[i - 1] + nownedsenders[i - 1];
1163bc8ff85bSBarry Smith     else starts[i] = starts[i - 1];
1164bc8ff85bSBarry Smith   }
1165bc8ff85bSBarry Smith 
11666aad120cSJose E. Roman   /* for each nontrivial globally owned node list all arriving processors */
1167bc8ff85bSBarry Smith   for (i = 0; i < nrecvs; i++) {
1168bc8ff85bSBarry Smith     for (j = 0; j < len[i]; j++) {
116930dcb7c9SBarry Smith       node = recvs[2 * i * nmax + 2 * j] - rstart;
1170f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) ownedsenders[starts[node]++] = source[i];
1171bc8ff85bSBarry Smith     }
1172bc8ff85bSBarry Smith   }
1173bc8ff85bSBarry Smith 
117407b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
117530dcb7c9SBarry Smith     starts[0] = 0;
117630dcb7c9SBarry Smith     for (i = 1; i < ng; i++) {
117730dcb7c9SBarry Smith       if (nownedsenders[i - 1] > 1) starts[i] = starts[i - 1] + nownedsenders[i - 1];
117830dcb7c9SBarry Smith       else starts[i] = starts[i - 1];
117930dcb7c9SBarry Smith     }
118030dcb7c9SBarry Smith     for (i = 0; i < ng; i++) {
118130dcb7c9SBarry Smith       if (nownedsenders[i] > 1) {
11829566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(comm, "[%d] global node %" PetscInt_FMT " local owner processors: ", rank, i + rstart));
118348a46eb9SPierre Jolivet         for (j = 0; j < nownedsenders[i]; j++) PetscCall(PetscSynchronizedPrintf(comm, "%" PetscInt_FMT " ", ownedsenders[starts[i] + j]));
11849566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(comm, "\n"));
118530dcb7c9SBarry Smith       }
118630dcb7c9SBarry Smith     }
11879566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
118807b52d57SBarry Smith   } /* -----------------------------------  */
118930dcb7c9SBarry Smith 
11903677ff5aSBarry Smith   /* wait on original sends */
11913a96401aSBarry Smith   if (nsends) {
11929566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nsends, &send_status));
11939566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nsends, send_waits, send_status));
11949566063dSJacob Faibussowitsch     PetscCall(PetscFree(send_status));
11953a96401aSBarry Smith   }
11969566063dSJacob Faibussowitsch   PetscCall(PetscFree(send_waits));
11979566063dSJacob Faibussowitsch   PetscCall(PetscFree(sends));
11989566063dSJacob Faibussowitsch   PetscCall(PetscFree(nprocs));
11993677ff5aSBarry Smith 
12003677ff5aSBarry Smith   /* pack messages to send back to local owners */
120130dcb7c9SBarry Smith   starts[0] = 0;
120230dcb7c9SBarry Smith   for (i = 1; i < ng; i++) {
120330dcb7c9SBarry Smith     if (nownedsenders[i - 1] > 1) starts[i] = starts[i - 1] + nownedsenders[i - 1];
120430dcb7c9SBarry Smith     else starts[i] = starts[i - 1];
120530dcb7c9SBarry Smith   }
120630dcb7c9SBarry Smith   nsends2 = nrecvs;
12079566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends2 + 1, &nprocs)); /* length of each message */
120830dcb7c9SBarry Smith   for (i = 0; i < nrecvs; i++) {
120930dcb7c9SBarry Smith     nprocs[i] = 1;
121030dcb7c9SBarry Smith     for (j = 0; j < len[i]; j++) {
121130dcb7c9SBarry Smith       node = recvs[2 * i * nmax + 2 * j] - rstart;
1212f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) nprocs[i] += 2 + nownedsenders[node];
121330dcb7c9SBarry Smith     }
121430dcb7c9SBarry Smith   }
1215f6e5521dSKarl Rupp   nt = 0;
1216f6e5521dSKarl Rupp   for (i = 0; i < nsends2; i++) nt += nprocs[i];
1217f6e5521dSKarl Rupp 
12189566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt + 1, &sends2));
12199566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends2 + 1, &starts2));
1220f6e5521dSKarl Rupp 
1221f6e5521dSKarl Rupp   starts2[0] = 0;
1222f6e5521dSKarl Rupp   for (i = 1; i < nsends2; i++) starts2[i] = starts2[i - 1] + nprocs[i - 1];
122330dcb7c9SBarry Smith   /*
122430dcb7c9SBarry Smith      Each message is 1 + nprocs[i] long, and consists of
122530dcb7c9SBarry Smith        (0) the number of nodes being sent back
122630dcb7c9SBarry Smith        (1) the local node number,
122730dcb7c9SBarry Smith        (2) the number of processors sharing it,
122830dcb7c9SBarry Smith        (3) the processors sharing it
122930dcb7c9SBarry Smith   */
123030dcb7c9SBarry Smith   for (i = 0; i < nsends2; i++) {
123130dcb7c9SBarry Smith     cnt                = 1;
123230dcb7c9SBarry Smith     sends2[starts2[i]] = 0;
123330dcb7c9SBarry Smith     for (j = 0; j < len[i]; j++) {
123430dcb7c9SBarry Smith       node = recvs[2 * i * nmax + 2 * j] - rstart;
123530dcb7c9SBarry Smith       if (nownedsenders[node] > 1) {
123630dcb7c9SBarry Smith         sends2[starts2[i]]++;
123730dcb7c9SBarry Smith         sends2[starts2[i] + cnt++] = recvs[2 * i * nmax + 2 * j + 1];
123830dcb7c9SBarry Smith         sends2[starts2[i] + cnt++] = nownedsenders[node];
12399566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(&sends2[starts2[i] + cnt], &ownedsenders[starts[node]], nownedsenders[node]));
124030dcb7c9SBarry Smith         cnt += nownedsenders[node];
124130dcb7c9SBarry Smith       }
124230dcb7c9SBarry Smith     }
124330dcb7c9SBarry Smith   }
124430dcb7c9SBarry Smith 
124530dcb7c9SBarry Smith   /* receive the message lengths */
124630dcb7c9SBarry Smith   nrecvs2 = nsends;
12479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs2 + 1, &lens2));
12489566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs2 + 1, &starts3));
12499566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs2 + 1, &recv_waits));
125048a46eb9SPierre Jolivet   for (i = 0; i < nrecvs2; i++) PetscCallMPI(MPI_Irecv(&lens2[i], 1, MPIU_INT, dest[i], tag2, comm, recv_waits + i));
1251d44834fbSBarry Smith 
12528a8e0b3aSBarry Smith   /* send the message lengths */
125348a46eb9SPierre Jolivet   for (i = 0; i < nsends2; i++) PetscCallMPI(MPI_Send(&nprocs[i], 1, MPIU_INT, source[i], tag2, comm));
12548a8e0b3aSBarry Smith 
1255d44834fbSBarry Smith   /* wait on receives of lens */
12560c468ba9SBarry Smith   if (nrecvs2) {
12579566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrecvs2, &recv_statuses));
12589566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrecvs2, recv_waits, recv_statuses));
12599566063dSJacob Faibussowitsch     PetscCall(PetscFree(recv_statuses));
12600c468ba9SBarry Smith   }
12619566063dSJacob Faibussowitsch   PetscCall(PetscFree(recv_waits));
1262d44834fbSBarry Smith 
126330dcb7c9SBarry Smith   starts3[0] = 0;
1264d44834fbSBarry Smith   nt         = 0;
126530dcb7c9SBarry Smith   for (i = 0; i < nrecvs2 - 1; i++) {
126630dcb7c9SBarry Smith     starts3[i + 1] = starts3[i] + lens2[i];
1267d44834fbSBarry Smith     nt += lens2[i];
126830dcb7c9SBarry Smith   }
126976466f69SStefano Zampini   if (nrecvs2) nt += lens2[nrecvs2 - 1];
1270d44834fbSBarry Smith 
12719566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt + 1, &recvs2));
12729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs2 + 1, &recv_waits));
127348a46eb9SPierre Jolivet   for (i = 0; i < nrecvs2; i++) PetscCallMPI(MPI_Irecv(recvs2 + starts3[i], lens2[i], MPIU_INT, dest[i], tag3, comm, recv_waits + i));
127430dcb7c9SBarry Smith 
127530dcb7c9SBarry Smith   /* send the messages */
12769566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends2 + 1, &send_waits));
127748a46eb9SPierre Jolivet   for (i = 0; i < nsends2; i++) PetscCallMPI(MPI_Isend(sends2 + starts2[i], nprocs[i], MPIU_INT, source[i], tag3, comm, send_waits + i));
127830dcb7c9SBarry Smith 
127930dcb7c9SBarry Smith   /* wait on receives */
12800c468ba9SBarry Smith   if (nrecvs2) {
12819566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrecvs2, &recv_statuses));
12829566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrecvs2, recv_waits, recv_statuses));
12839566063dSJacob Faibussowitsch     PetscCall(PetscFree(recv_statuses));
12840c468ba9SBarry Smith   }
12859566063dSJacob Faibussowitsch   PetscCall(PetscFree(recv_waits));
12869566063dSJacob Faibussowitsch   PetscCall(PetscFree(nprocs));
128730dcb7c9SBarry Smith 
128807b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
128930dcb7c9SBarry Smith     cnt = 0;
129030dcb7c9SBarry Smith     for (i = 0; i < nrecvs2; i++) {
129130dcb7c9SBarry Smith       nt = recvs2[cnt++];
129230dcb7c9SBarry Smith       for (j = 0; j < nt; j++) {
12939566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(comm, "[%d] local node %" PetscInt_FMT " number of subdomains %" PetscInt_FMT ": ", rank, recvs2[cnt], recvs2[cnt + 1]));
129448a46eb9SPierre Jolivet         for (k = 0; k < recvs2[cnt + 1]; k++) PetscCall(PetscSynchronizedPrintf(comm, "%" PetscInt_FMT " ", recvs2[cnt + 2 + k]));
129530dcb7c9SBarry Smith         cnt += 2 + recvs2[cnt + 1];
12969566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(comm, "\n"));
129730dcb7c9SBarry Smith       }
129830dcb7c9SBarry Smith     }
12999566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
130007b52d57SBarry Smith   } /* -----------------------------------  */
130130dcb7c9SBarry Smith 
130230dcb7c9SBarry Smith   /* count number subdomains for each local node */
13039566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(size, &nprocs));
130430dcb7c9SBarry Smith   cnt = 0;
130530dcb7c9SBarry Smith   for (i = 0; i < nrecvs2; i++) {
130630dcb7c9SBarry Smith     nt = recvs2[cnt++];
130730dcb7c9SBarry Smith     for (j = 0; j < nt; j++) {
1308f6e5521dSKarl Rupp       for (k = 0; k < recvs2[cnt + 1]; k++) nprocs[recvs2[cnt + 2 + k]]++;
130930dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt + 1];
131030dcb7c9SBarry Smith     }
131130dcb7c9SBarry Smith   }
13129371c9d4SSatish Balay   nt = 0;
13139371c9d4SSatish Balay   for (i = 0; i < size; i++) nt += (nprocs[i] > 0);
131430dcb7c9SBarry Smith   *nproc = nt;
13159566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt + 1, procs));
13169566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt + 1, numprocs));
13179566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt + 1, indices));
13180298fd71SBarry Smith   for (i = 0; i < nt + 1; i++) (*indices)[i] = NULL;
13199566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size, &bprocs));
132030dcb7c9SBarry Smith   cnt = 0;
132130dcb7c9SBarry Smith   for (i = 0; i < size; i++) {
132230dcb7c9SBarry Smith     if (nprocs[i] > 0) {
132330dcb7c9SBarry Smith       bprocs[i]        = cnt;
132430dcb7c9SBarry Smith       (*procs)[cnt]    = i;
132530dcb7c9SBarry Smith       (*numprocs)[cnt] = nprocs[i];
13269566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nprocs[i], &(*indices)[cnt]));
132730dcb7c9SBarry Smith       cnt++;
132830dcb7c9SBarry Smith     }
132930dcb7c9SBarry Smith   }
133030dcb7c9SBarry Smith 
133130dcb7c9SBarry Smith   /* make the list of subdomains for each nontrivial local node */
13329566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(*numprocs, nt));
133330dcb7c9SBarry Smith   cnt = 0;
133430dcb7c9SBarry Smith   for (i = 0; i < nrecvs2; i++) {
133530dcb7c9SBarry Smith     nt = recvs2[cnt++];
133630dcb7c9SBarry Smith     for (j = 0; j < nt; j++) {
1337f6e5521dSKarl Rupp       for (k = 0; k < recvs2[cnt + 1]; k++) (*indices)[bprocs[recvs2[cnt + 2 + k]]][(*numprocs)[bprocs[recvs2[cnt + 2 + k]]]++] = recvs2[cnt];
133830dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt + 1];
133930dcb7c9SBarry Smith     }
134030dcb7c9SBarry Smith   }
13419566063dSJacob Faibussowitsch   PetscCall(PetscFree(bprocs));
13429566063dSJacob Faibussowitsch   PetscCall(PetscFree(recvs2));
134330dcb7c9SBarry Smith 
134407b52d57SBarry Smith   /* sort the node indexing by their global numbers */
134507b52d57SBarry Smith   nt = *nproc;
134607b52d57SBarry Smith   for (i = 0; i < nt; i++) {
13479566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1((*numprocs)[i], &tmp));
1348f6e5521dSKarl Rupp     for (j = 0; j < (*numprocs)[i]; j++) tmp[j] = lindices[(*indices)[i][j]];
13499566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithArray((*numprocs)[i], tmp, (*indices)[i]));
13509566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmp));
135107b52d57SBarry Smith   }
135207b52d57SBarry Smith 
135307b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
135430dcb7c9SBarry Smith     nt = *nproc;
135530dcb7c9SBarry Smith     for (i = 0; i < nt; i++) {
13569566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedPrintf(comm, "[%d] subdomain %" PetscInt_FMT " number of indices %" PetscInt_FMT ": ", rank, (*procs)[i], (*numprocs)[i]));
135748a46eb9SPierre Jolivet       for (j = 0; j < (*numprocs)[i]; j++) PetscCall(PetscSynchronizedPrintf(comm, "%" PetscInt_FMT " ", (*indices)[i][j]));
13589566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedPrintf(comm, "\n"));
135930dcb7c9SBarry Smith     }
13609566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
136107b52d57SBarry Smith   } /* -----------------------------------  */
136230dcb7c9SBarry Smith 
136330dcb7c9SBarry Smith   /* wait on sends */
136430dcb7c9SBarry Smith   if (nsends2) {
13659566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nsends2, &send_status));
13669566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nsends2, send_waits, send_status));
13679566063dSJacob Faibussowitsch     PetscCall(PetscFree(send_status));
136830dcb7c9SBarry Smith   }
136930dcb7c9SBarry Smith 
13709566063dSJacob Faibussowitsch   PetscCall(PetscFree(starts3));
13719566063dSJacob Faibussowitsch   PetscCall(PetscFree(dest));
13729566063dSJacob Faibussowitsch   PetscCall(PetscFree(send_waits));
13733677ff5aSBarry Smith 
13749566063dSJacob Faibussowitsch   PetscCall(PetscFree(nownedsenders));
13759566063dSJacob Faibussowitsch   PetscCall(PetscFree(ownedsenders));
13769566063dSJacob Faibussowitsch   PetscCall(PetscFree(starts));
13779566063dSJacob Faibussowitsch   PetscCall(PetscFree(starts2));
13789566063dSJacob Faibussowitsch   PetscCall(PetscFree(lens2));
137989d82c54SBarry Smith 
13809566063dSJacob Faibussowitsch   PetscCall(PetscFree(source));
13819566063dSJacob Faibussowitsch   PetscCall(PetscFree(len));
13829566063dSJacob Faibussowitsch   PetscCall(PetscFree(recvs));
13839566063dSJacob Faibussowitsch   PetscCall(PetscFree(nprocs));
13849566063dSJacob Faibussowitsch   PetscCall(PetscFree(sends2));
138524cf384cSBarry Smith 
138624cf384cSBarry Smith   /* put the information about myself as the first entry in the list */
138724cf384cSBarry Smith   first_procs    = (*procs)[0];
138824cf384cSBarry Smith   first_numprocs = (*numprocs)[0];
138924cf384cSBarry Smith   first_indices  = (*indices)[0];
139024cf384cSBarry Smith   for (i = 0; i < *nproc; i++) {
139124cf384cSBarry Smith     if ((*procs)[i] == rank) {
139224cf384cSBarry Smith       (*procs)[0]    = (*procs)[i];
139324cf384cSBarry Smith       (*numprocs)[0] = (*numprocs)[i];
139424cf384cSBarry Smith       (*indices)[0]  = (*indices)[i];
139524cf384cSBarry Smith       (*procs)[i]    = first_procs;
139624cf384cSBarry Smith       (*numprocs)[i] = first_numprocs;
139724cf384cSBarry Smith       (*indices)[i]  = first_indices;
139824cf384cSBarry Smith       break;
139924cf384cSBarry Smith     }
140024cf384cSBarry Smith   }
1401268a049cSStefano Zampini 
1402268a049cSStefano Zampini   /* save info for reuse */
1403268a049cSStefano Zampini   mapping->info_nproc    = *nproc;
1404268a049cSStefano Zampini   mapping->info_procs    = *procs;
1405268a049cSStefano Zampini   mapping->info_numprocs = *numprocs;
1406268a049cSStefano Zampini   mapping->info_indices  = *indices;
1407268a049cSStefano Zampini   mapping->info_cached   = PETSC_TRUE;
140889d82c54SBarry Smith   PetscFunctionReturn(0);
140989d82c54SBarry Smith }
141089d82c54SBarry Smith 
14116a818285SBarry Smith /*@C
1412cab54364SBarry Smith     ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by `ISLocalToGlobalMappingGetBlockInfo()`
14136a818285SBarry Smith 
1414c3339decSBarry Smith     Collective
14156a818285SBarry Smith 
1416f899ff85SJose E. Roman     Input Parameter:
14176a818285SBarry Smith .   mapping - the mapping from local to global indexing
14186a818285SBarry Smith 
1419d8d19677SJose E. Roman     Output Parameters:
14206a818285SBarry Smith +   nproc - number of processors that are connected to this one
14216a818285SBarry Smith .   proc - neighboring processors
14226a818285SBarry Smith .   numproc - number of indices for each processor
14236a818285SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
14246a818285SBarry Smith 
14256a818285SBarry Smith     Level: advanced
14266a818285SBarry Smith 
1427cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1428db781477SPatrick Sanan           `ISLocalToGlobalMappingGetInfo()`
14296a818285SBarry Smith @*/
1430d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[])
1431d71ae5a4SJacob Faibussowitsch {
14326a818285SBarry Smith   PetscFunctionBegin;
1433cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
1434268a049cSStefano Zampini   if (mapping->info_free) {
14359566063dSJacob Faibussowitsch     PetscCall(PetscFree(*numprocs));
14366a818285SBarry Smith     if (*indices) {
1437268a049cSStefano Zampini       PetscInt i;
1438268a049cSStefano Zampini 
14399566063dSJacob Faibussowitsch       PetscCall(PetscFree((*indices)[0]));
144048a46eb9SPierre Jolivet       for (i = 1; i < *nproc; i++) PetscCall(PetscFree((*indices)[i]));
14419566063dSJacob Faibussowitsch       PetscCall(PetscFree(*indices));
14426a818285SBarry Smith     }
1443268a049cSStefano Zampini   }
1444268a049cSStefano Zampini   *nproc    = 0;
1445268a049cSStefano Zampini   *procs    = NULL;
1446268a049cSStefano Zampini   *numprocs = NULL;
1447268a049cSStefano Zampini   *indices  = NULL;
14486a818285SBarry Smith   PetscFunctionReturn(0);
14496a818285SBarry Smith }
14506a818285SBarry Smith 
14516a818285SBarry Smith /*@C
14526a818285SBarry Smith     ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each processor and
14536a818285SBarry Smith      each index shared by more than one processor
14546a818285SBarry Smith 
1455c3339decSBarry Smith     Collective
14566a818285SBarry Smith 
1457f899ff85SJose E. Roman     Input Parameter:
14586a818285SBarry Smith .   mapping - the mapping from local to global indexing
14596a818285SBarry Smith 
1460d8d19677SJose E. Roman     Output Parameters:
14616a818285SBarry Smith +   nproc - number of processors that are connected to this one
14626a818285SBarry Smith .   proc - neighboring processors
14636a818285SBarry Smith .   numproc - number of indices for each subdomain (processor)
14646a818285SBarry Smith -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
14656a818285SBarry Smith 
14666a818285SBarry Smith     Level: advanced
14676a818285SBarry Smith 
1468cab54364SBarry Smith     Note:
1469cab54364SBarry Smith     The user needs to call `ISLocalToGlobalMappingRestoreInfo()` when the data is no longer needed.
14701bd0b88eSStefano Zampini 
14716a818285SBarry Smith     Fortran Usage:
1472cab54364SBarry Smith .vb
14736a818285SBarry Smith         PetscInt indices[nproc][numprocmax],ierr)
1474cab54364SBarry Smith         ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
1475cab54364SBarry Smith         ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
1476cab54364SBarry Smith .ve
1477cab54364SBarry Smith 
1478cab54364SBarry Smith     Fortran Note:
1479cab54364SBarry Smith         There is no `ISLocalToGlobalMappingRestoreInfo()` in Fortran. You must make sure that procs[], numprocs[] and
14806a818285SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
14816a818285SBarry Smith 
1482cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1483db781477SPatrick Sanan           `ISLocalToGlobalMappingRestoreInfo()`
14846a818285SBarry Smith @*/
1485d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[])
1486d71ae5a4SJacob Faibussowitsch {
14878a1f772fSStefano Zampini   PetscInt **bindices = NULL, *bnumprocs = NULL, bs, i, j, k;
14886a818285SBarry Smith 
14896a818285SBarry Smith   PetscFunctionBegin;
14906a818285SBarry Smith   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
14918a1f772fSStefano Zampini   bs = mapping->bs;
14929566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockInfo(mapping, nproc, procs, &bnumprocs, &bindices));
1493268a049cSStefano Zampini   if (bs > 1) { /* we need to expand the cached info */
14949566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(*nproc, &*indices));
14959566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(*nproc, &*numprocs));
14966a818285SBarry Smith     for (i = 0; i < *nproc; i++) {
14979566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(bs * bnumprocs[i], &(*indices)[i]));
1498268a049cSStefano Zampini       for (j = 0; j < bnumprocs[i]; j++) {
1499ad540459SPierre Jolivet         for (k = 0; k < bs; k++) (*indices)[i][j * bs + k] = bs * bindices[i][j] + k;
15006a818285SBarry Smith       }
1501268a049cSStefano Zampini       (*numprocs)[i] = bnumprocs[i] * bs;
15026a818285SBarry Smith     }
1503268a049cSStefano Zampini     mapping->info_free = PETSC_TRUE;
1504268a049cSStefano Zampini   } else {
1505268a049cSStefano Zampini     *numprocs = bnumprocs;
1506268a049cSStefano Zampini     *indices  = bindices;
15076a818285SBarry Smith   }
15086a818285SBarry Smith   PetscFunctionReturn(0);
15096a818285SBarry Smith }
15106a818285SBarry Smith 
151107b52d57SBarry Smith /*@C
1512cab54364SBarry Smith     ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by `ISLocalToGlobalMappingGetInfo()`
151389d82c54SBarry Smith 
1514c3339decSBarry Smith     Collective
151507b52d57SBarry Smith 
1516f899ff85SJose E. Roman     Input Parameter:
151707b52d57SBarry Smith .   mapping - the mapping from local to global indexing
151807b52d57SBarry Smith 
1519d8d19677SJose E. Roman     Output Parameters:
152007b52d57SBarry Smith +   nproc - number of processors that are connected to this one
152107b52d57SBarry Smith .   proc - neighboring processors
152207b52d57SBarry Smith .   numproc - number of indices for each processor
152307b52d57SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
152407b52d57SBarry Smith 
152507b52d57SBarry Smith     Level: advanced
152607b52d57SBarry Smith 
1527cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1528db781477SPatrick Sanan           `ISLocalToGlobalMappingGetInfo()`
152907b52d57SBarry Smith @*/
1530d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[])
1531d71ae5a4SJacob Faibussowitsch {
153207b52d57SBarry Smith   PetscFunctionBegin;
15339566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockInfo(mapping, nproc, procs, numprocs, indices));
153407b52d57SBarry Smith   PetscFunctionReturn(0);
153507b52d57SBarry Smith }
153686994e45SJed Brown 
153786994e45SJed Brown /*@C
1538cab54364SBarry Smith     ISLocalToGlobalMappingGetNodeInfo - Gets the neighbor information for each MPI rank
15391bd0b88eSStefano Zampini 
1540c3339decSBarry Smith     Collective
15411bd0b88eSStefano Zampini 
1542f899ff85SJose E. Roman     Input Parameter:
15431bd0b88eSStefano Zampini .   mapping - the mapping from local to global indexing
15441bd0b88eSStefano Zampini 
1545d8d19677SJose E. Roman     Output Parameters:
1546cab54364SBarry Smith +   nnodes - number of local nodes (same `ISLocalToGlobalMappingGetSize()`)
15471bd0b88eSStefano Zampini .   count - number of neighboring processors per node
15481bd0b88eSStefano Zampini -   indices - indices of processes sharing the node (sorted)
15491bd0b88eSStefano Zampini 
15501bd0b88eSStefano Zampini     Level: advanced
15511bd0b88eSStefano Zampini 
1552cab54364SBarry Smith     Note:
1553cab54364SBarry Smith     The user needs to call `ISLocalToGlobalMappingRestoreInfo()` when the data is no longer needed.
15541bd0b88eSStefano Zampini 
1555cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1556db781477SPatrick Sanan           `ISLocalToGlobalMappingGetInfo()`, `ISLocalToGlobalMappingRestoreNodeInfo()`
15571bd0b88eSStefano Zampini @*/
1558d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetNodeInfo(ISLocalToGlobalMapping mapping, PetscInt *nnodes, PetscInt *count[], PetscInt **indices[])
1559d71ae5a4SJacob Faibussowitsch {
15601bd0b88eSStefano Zampini   PetscInt n;
15611bd0b88eSStefano Zampini 
15621bd0b88eSStefano Zampini   PetscFunctionBegin;
15631bd0b88eSStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
15649566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetSize(mapping, &n));
15651bd0b88eSStefano Zampini   if (!mapping->info_nodec) {
15661bd0b88eSStefano Zampini     PetscInt i, m, n_neigh, *neigh, *n_shared, **shared;
15671bd0b88eSStefano Zampini 
15689566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(n + 1, &mapping->info_nodec, n, &mapping->info_nodei));
15699566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetInfo(mapping, &n_neigh, &neigh, &n_shared, &shared));
1570ad540459SPierre Jolivet     for (i = 0; i < n; i++) mapping->info_nodec[i] = 1;
1571071fcb05SBarry Smith     m                      = n;
1572071fcb05SBarry Smith     mapping->info_nodec[n] = 0;
15731bd0b88eSStefano Zampini     for (i = 1; i < n_neigh; i++) {
15741bd0b88eSStefano Zampini       PetscInt j;
15751bd0b88eSStefano Zampini 
15761bd0b88eSStefano Zampini       m += n_shared[i];
15771bd0b88eSStefano Zampini       for (j = 0; j < n_shared[i]; j++) mapping->info_nodec[shared[i][j]] += 1;
15781bd0b88eSStefano Zampini     }
15799566063dSJacob Faibussowitsch     if (n) PetscCall(PetscMalloc1(m, &mapping->info_nodei[0]));
15801bd0b88eSStefano Zampini     for (i = 1; i < n; i++) mapping->info_nodei[i] = mapping->info_nodei[i - 1] + mapping->info_nodec[i - 1];
15819566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(mapping->info_nodec, n));
15829371c9d4SSatish Balay     for (i = 0; i < n; i++) {
15839371c9d4SSatish Balay       mapping->info_nodec[i]    = 1;
15849371c9d4SSatish Balay       mapping->info_nodei[i][0] = neigh[0];
15859371c9d4SSatish Balay     }
15861bd0b88eSStefano Zampini     for (i = 1; i < n_neigh; i++) {
15871bd0b88eSStefano Zampini       PetscInt j;
15881bd0b88eSStefano Zampini 
15891bd0b88eSStefano Zampini       for (j = 0; j < n_shared[i]; j++) {
15901bd0b88eSStefano Zampini         PetscInt k = shared[i][j];
15911bd0b88eSStefano Zampini 
15921bd0b88eSStefano Zampini         mapping->info_nodei[k][mapping->info_nodec[k]] = neigh[i];
15931bd0b88eSStefano Zampini         mapping->info_nodec[k] += 1;
15941bd0b88eSStefano Zampini       }
15951bd0b88eSStefano Zampini     }
15969566063dSJacob Faibussowitsch     for (i = 0; i < n; i++) PetscCall(PetscSortRemoveDupsInt(&mapping->info_nodec[i], mapping->info_nodei[i]));
15979566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreInfo(mapping, &n_neigh, &neigh, &n_shared, &shared));
15981bd0b88eSStefano Zampini   }
15991bd0b88eSStefano Zampini   if (nnodes) *nnodes = n;
16001bd0b88eSStefano Zampini   if (count) *count = mapping->info_nodec;
16011bd0b88eSStefano Zampini   if (indices) *indices = mapping->info_nodei;
16021bd0b88eSStefano Zampini   PetscFunctionReturn(0);
16031bd0b88eSStefano Zampini }
16041bd0b88eSStefano Zampini 
16051bd0b88eSStefano Zampini /*@C
1606cab54364SBarry Smith     ISLocalToGlobalMappingRestoreNodeInfo - Frees the memory allocated by `ISLocalToGlobalMappingGetNodeInfo()`
16071bd0b88eSStefano Zampini 
1608c3339decSBarry Smith     Collective
16091bd0b88eSStefano Zampini 
1610f899ff85SJose E. Roman     Input Parameter:
16111bd0b88eSStefano Zampini .   mapping - the mapping from local to global indexing
16121bd0b88eSStefano Zampini 
1613d8d19677SJose E. Roman     Output Parameters:
16141bd0b88eSStefano Zampini +   nnodes - number of local nodes
16151bd0b88eSStefano Zampini .   count - number of neighboring processors per node
16161bd0b88eSStefano Zampini -   indices - indices of processes sharing the node (sorted)
16171bd0b88eSStefano Zampini 
16181bd0b88eSStefano Zampini     Level: advanced
16191bd0b88eSStefano Zampini 
1620cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1621db781477SPatrick Sanan           `ISLocalToGlobalMappingGetInfo()`
16221bd0b88eSStefano Zampini @*/
1623d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreNodeInfo(ISLocalToGlobalMapping mapping, PetscInt *nnodes, PetscInt *count[], PetscInt **indices[])
1624d71ae5a4SJacob Faibussowitsch {
16251bd0b88eSStefano Zampini   PetscFunctionBegin;
16261bd0b88eSStefano Zampini   PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1);
16271bd0b88eSStefano Zampini   if (nnodes) *nnodes = 0;
16281bd0b88eSStefano Zampini   if (count) *count = NULL;
16291bd0b88eSStefano Zampini   if (indices) *indices = NULL;
16301bd0b88eSStefano Zampini   PetscFunctionReturn(0);
16311bd0b88eSStefano Zampini }
16321bd0b88eSStefano Zampini 
16331bd0b88eSStefano Zampini /*@C
1634107e9a97SBarry Smith    ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped
163586994e45SJed Brown 
163686994e45SJed Brown    Not Collective
163786994e45SJed Brown 
16384165533cSJose E. Roman    Input Parameter:
163986994e45SJed Brown . ltog - local to global mapping
164086994e45SJed Brown 
16414165533cSJose E. Roman    Output Parameter:
1642cab54364SBarry Smith . array - array of indices, the length of this array may be obtained with `ISLocalToGlobalMappingGetSize()`
164386994e45SJed Brown 
164486994e45SJed Brown    Level: advanced
164586994e45SJed Brown 
1646cab54364SBarry Smith    Note:
1647cab54364SBarry Smith     `ISLocalToGlobalMappingGetSize()` returns the length the this array
1648107e9a97SBarry Smith 
1649cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingRestoreIndices()`, `ISLocalToGlobalMappingGetBlockIndices()`, `ISLocalToGlobalMappingRestoreBlockIndices()`
165086994e45SJed Brown @*/
1651d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog, const PetscInt **array)
1652d71ae5a4SJacob Faibussowitsch {
165386994e45SJed Brown   PetscFunctionBegin;
165486994e45SJed Brown   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
165586994e45SJed Brown   PetscValidPointer(array, 2);
165645b6f7e9SBarry Smith   if (ltog->bs == 1) {
165786994e45SJed Brown     *array = ltog->indices;
165845b6f7e9SBarry Smith   } else {
165945b6f7e9SBarry Smith     PetscInt       *jj, k, i, j, n = ltog->n, bs = ltog->bs;
166045b6f7e9SBarry Smith     const PetscInt *ii;
166145b6f7e9SBarry Smith 
16629566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(bs * n, &jj));
166345b6f7e9SBarry Smith     *array = jj;
166445b6f7e9SBarry Smith     k      = 0;
166545b6f7e9SBarry Smith     ii     = ltog->indices;
166645b6f7e9SBarry Smith     for (i = 0; i < n; i++)
16679371c9d4SSatish Balay       for (j = 0; j < bs; j++) jj[k++] = bs * ii[i] + j;
166845b6f7e9SBarry Smith   }
166986994e45SJed Brown   PetscFunctionReturn(0);
167086994e45SJed Brown }
167186994e45SJed Brown 
167286994e45SJed Brown /*@C
1673cab54364SBarry Smith    ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with `ISLocalToGlobalMappingGetIndices()`
167486994e45SJed Brown 
167586994e45SJed Brown    Not Collective
167686994e45SJed Brown 
16774165533cSJose E. Roman    Input Parameters:
167886994e45SJed Brown + ltog - local to global mapping
167986994e45SJed Brown - array - array of indices
168086994e45SJed Brown 
168186994e45SJed Brown    Level: advanced
168286994e45SJed Brown 
1683cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingGetIndices()`
168486994e45SJed Brown @*/
1685d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog, const PetscInt **array)
1686d71ae5a4SJacob Faibussowitsch {
168786994e45SJed Brown   PetscFunctionBegin;
168886994e45SJed Brown   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
168986994e45SJed Brown   PetscValidPointer(array, 2);
1690c9cc58a2SBarry Smith   PetscCheck(ltog->bs != 1 || *array == ltog->indices, PETSC_COMM_SELF, PETSC_ERR_ARG_BADPTR, "Trying to return mismatched pointer");
169145b6f7e9SBarry Smith 
169248a46eb9SPierre Jolivet   if (ltog->bs > 1) PetscCall(PetscFree(*(void **)array));
169345b6f7e9SBarry Smith   PetscFunctionReturn(0);
169445b6f7e9SBarry Smith }
169545b6f7e9SBarry Smith 
169645b6f7e9SBarry Smith /*@C
169745b6f7e9SBarry Smith    ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block
169845b6f7e9SBarry Smith 
169945b6f7e9SBarry Smith    Not Collective
170045b6f7e9SBarry Smith 
17014165533cSJose E. Roman    Input Parameter:
170245b6f7e9SBarry Smith . ltog - local to global mapping
170345b6f7e9SBarry Smith 
17044165533cSJose E. Roman    Output Parameter:
170545b6f7e9SBarry Smith . array - array of indices
170645b6f7e9SBarry Smith 
170745b6f7e9SBarry Smith    Level: advanced
170845b6f7e9SBarry Smith 
1709cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingRestoreBlockIndices()`
171045b6f7e9SBarry Smith @*/
1711d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog, const PetscInt **array)
1712d71ae5a4SJacob Faibussowitsch {
171345b6f7e9SBarry Smith   PetscFunctionBegin;
171445b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
171545b6f7e9SBarry Smith   PetscValidPointer(array, 2);
171645b6f7e9SBarry Smith   *array = ltog->indices;
171745b6f7e9SBarry Smith   PetscFunctionReturn(0);
171845b6f7e9SBarry Smith }
171945b6f7e9SBarry Smith 
172045b6f7e9SBarry Smith /*@C
1721cab54364SBarry Smith    ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with `ISLocalToGlobalMappingGetBlockIndices()`
172245b6f7e9SBarry Smith 
172345b6f7e9SBarry Smith    Not Collective
172445b6f7e9SBarry Smith 
17254165533cSJose E. Roman    Input Parameters:
172645b6f7e9SBarry Smith + ltog - local to global mapping
172745b6f7e9SBarry Smith - array - array of indices
172845b6f7e9SBarry Smith 
172945b6f7e9SBarry Smith    Level: advanced
173045b6f7e9SBarry Smith 
1731cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingGetIndices()`
173245b6f7e9SBarry Smith @*/
1733d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog, const PetscInt **array)
1734d71ae5a4SJacob Faibussowitsch {
173545b6f7e9SBarry Smith   PetscFunctionBegin;
173645b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
173745b6f7e9SBarry Smith   PetscValidPointer(array, 2);
173808401ef6SPierre Jolivet   PetscCheck(*array == ltog->indices, PETSC_COMM_SELF, PETSC_ERR_ARG_BADPTR, "Trying to return mismatched pointer");
17390298fd71SBarry Smith   *array = NULL;
174086994e45SJed Brown   PetscFunctionReturn(0);
174186994e45SJed Brown }
1742f7efa3c7SJed Brown 
1743f7efa3c7SJed Brown /*@C
1744f7efa3c7SJed Brown    ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings
1745f7efa3c7SJed Brown 
1746f7efa3c7SJed Brown    Not Collective
1747f7efa3c7SJed Brown 
17484165533cSJose E. Roman    Input Parameters:
1749f7efa3c7SJed Brown + comm - communicator for the new mapping, must contain the communicator of every mapping to concatenate
1750f7efa3c7SJed Brown . n - number of mappings to concatenate
1751f7efa3c7SJed Brown - ltogs - local to global mappings
1752f7efa3c7SJed Brown 
17534165533cSJose E. Roman    Output Parameter:
1754f7efa3c7SJed Brown . ltogcat - new mapping
1755f7efa3c7SJed Brown 
1756f7efa3c7SJed Brown    Level: advanced
1757f7efa3c7SJed Brown 
1758cab54364SBarry Smith    Note:
1759cab54364SBarry Smith    This currently always returns a mapping with block size of 1
1760cab54364SBarry Smith 
1761cab54364SBarry Smith    Developer Note:
1762cab54364SBarry Smith    If all the input mapping have the same block size we could easily handle that as a special case
1763cab54364SBarry Smith 
1764cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingCreate()`
1765f7efa3c7SJed Brown @*/
1766d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm, PetscInt n, const ISLocalToGlobalMapping ltogs[], ISLocalToGlobalMapping *ltogcat)
1767d71ae5a4SJacob Faibussowitsch {
1768f7efa3c7SJed Brown   PetscInt i, cnt, m, *idx;
1769f7efa3c7SJed Brown 
1770f7efa3c7SJed Brown   PetscFunctionBegin;
177108401ef6SPierre Jolivet   PetscCheck(n >= 0, comm, PETSC_ERR_ARG_OUTOFRANGE, "Must have a non-negative number of mappings, given %" PetscInt_FMT, n);
1772f7efa3c7SJed Brown   if (n > 0) PetscValidPointer(ltogs, 3);
1773f7efa3c7SJed Brown   for (i = 0; i < n; i++) PetscValidHeaderSpecific(ltogs[i], IS_LTOGM_CLASSID, 3);
1774f7efa3c7SJed Brown   PetscValidPointer(ltogcat, 4);
1775f7efa3c7SJed Brown   for (cnt = 0, i = 0; i < n; i++) {
17769566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(ltogs[i], &m));
1777f7efa3c7SJed Brown     cnt += m;
1778f7efa3c7SJed Brown   }
17799566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(cnt, &idx));
1780f7efa3c7SJed Brown   for (cnt = 0, i = 0; i < n; i++) {
1781f7efa3c7SJed Brown     const PetscInt *subidx;
17829566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(ltogs[i], &m));
17839566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetIndices(ltogs[i], &subidx));
17849566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(&idx[cnt], subidx, m));
17859566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreIndices(ltogs[i], &subidx));
1786f7efa3c7SJed Brown     cnt += m;
1787f7efa3c7SJed Brown   }
17889566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(comm, 1, cnt, idx, PETSC_OWN_POINTER, ltogcat));
1789f7efa3c7SJed Brown   PetscFunctionReturn(0);
1790f7efa3c7SJed Brown }
179104a59952SBarry Smith 
1792413f72f0SBarry Smith /*MC
1793cab54364SBarry Smith       ISLOCALTOGLOBALMAPPINGBASIC - basic implementation of the `ISLocalToGlobalMapping` object. When `ISGlobalToLocalMappingApply()` is
1794413f72f0SBarry Smith                                     used this is good for only small and moderate size problems.
1795413f72f0SBarry Smith 
1796413f72f0SBarry Smith    Options Database Keys:
1797a2b725a8SWilliam Gropp .   -islocaltoglobalmapping_type basic - select this method
1798413f72f0SBarry Smith 
1799413f72f0SBarry Smith    Level: beginner
1800413f72f0SBarry Smith 
1801cab54364SBarry Smith    Developer Note:
1802cab54364SBarry Smith    This stores all the mapping information on each MPI rank.
1803cab54364SBarry Smith 
1804cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`, `ISLOCALTOGLOBALMAPPINGHASH`
1805413f72f0SBarry Smith M*/
1806d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Basic(ISLocalToGlobalMapping ltog)
1807d71ae5a4SJacob Faibussowitsch {
1808413f72f0SBarry Smith   PetscFunctionBegin;
1809413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Basic;
1810413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Basic;
1811413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Basic;
1812413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Basic;
1813413f72f0SBarry Smith   PetscFunctionReturn(0);
1814413f72f0SBarry Smith }
1815413f72f0SBarry Smith 
1816413f72f0SBarry Smith /*MC
1817cab54364SBarry Smith       ISLOCALTOGLOBALMAPPINGHASH - hash implementation of the `ISLocalToGlobalMapping` object. When `ISGlobalToLocalMappingApply()` is
1818ed56e8eaSBarry Smith                                     used this is good for large memory problems.
1819413f72f0SBarry Smith 
1820413f72f0SBarry Smith    Options Database Keys:
1821a2b725a8SWilliam Gropp .   -islocaltoglobalmapping_type hash - select this method
1822413f72f0SBarry Smith 
1823413f72f0SBarry Smith    Level: beginner
1824413f72f0SBarry Smith 
1825cab54364SBarry Smith    Note:
1826cab54364SBarry Smith     This is selected automatically for large problems if the user does not set the type.
1827cab54364SBarry Smith 
1828cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`, `ISLOCALTOGLOBALMAPPINGBASIC`
1829413f72f0SBarry Smith M*/
1830d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Hash(ISLocalToGlobalMapping ltog)
1831d71ae5a4SJacob Faibussowitsch {
1832413f72f0SBarry Smith   PetscFunctionBegin;
1833413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Hash;
1834413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Hash;
1835413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Hash;
1836413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Hash;
1837413f72f0SBarry Smith   PetscFunctionReturn(0);
1838413f72f0SBarry Smith }
1839413f72f0SBarry Smith 
1840413f72f0SBarry Smith /*@C
1841cab54364SBarry Smith     ISLocalToGlobalMappingRegister -  Registers a method for applying a global to local mapping with an `ISLocalToGlobalMapping`
1842413f72f0SBarry Smith 
1843413f72f0SBarry Smith    Not Collective
1844413f72f0SBarry Smith 
1845413f72f0SBarry Smith    Input Parameters:
1846413f72f0SBarry Smith +  sname - name of a new method
1847413f72f0SBarry Smith -  routine_create - routine to create method context
1848413f72f0SBarry Smith 
1849413f72f0SBarry Smith    Sample usage:
1850413f72f0SBarry Smith .vb
1851413f72f0SBarry Smith    ISLocalToGlobalMappingRegister("my_mapper",MyCreate);
1852413f72f0SBarry Smith .ve
1853413f72f0SBarry Smith 
1854ed56e8eaSBarry Smith    Then, your mapping can be chosen with the procedural interface via
1855413f72f0SBarry Smith $     ISLocalToGlobalMappingSetType(ltog,"my_mapper")
1856413f72f0SBarry Smith    or at runtime via the option
1857ed56e8eaSBarry Smith $     -islocaltoglobalmapping_type my_mapper
1858413f72f0SBarry Smith 
1859413f72f0SBarry Smith    Level: advanced
1860413f72f0SBarry Smith 
1861cab54364SBarry Smith    Note:
1862cab54364SBarry Smith    `ISLocalToGlobalMappingRegister()` may be called multiple times to add several user-defined mappings.
1863413f72f0SBarry Smith 
1864cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingRegisterAll()`, `ISLocalToGlobalMappingRegisterDestroy()`, `ISLOCALTOGLOBALMAPPINGBASIC`,
1865cab54364SBarry Smith           `ISLOCALTOGLOBALMAPPINGHASH`, `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApply()`
1866413f72f0SBarry Smith @*/
1867d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRegister(const char sname[], PetscErrorCode (*function)(ISLocalToGlobalMapping))
1868d71ae5a4SJacob Faibussowitsch {
1869413f72f0SBarry Smith   PetscFunctionBegin;
18709566063dSJacob Faibussowitsch   PetscCall(ISInitializePackage());
18719566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&ISLocalToGlobalMappingList, sname, function));
1872413f72f0SBarry Smith   PetscFunctionReturn(0);
1873413f72f0SBarry Smith }
1874413f72f0SBarry Smith 
1875413f72f0SBarry Smith /*@C
1876cab54364SBarry Smith    ISLocalToGlobalMappingSetType - Sets the implementation type `ISLocalToGlobalMapping` will use
1877413f72f0SBarry Smith 
1878c3339decSBarry Smith    Logically Collective
1879413f72f0SBarry Smith 
1880413f72f0SBarry Smith    Input Parameters:
1881cab54364SBarry Smith +  ltog - the `ISLocalToGlobalMapping` object
1882413f72f0SBarry Smith -  type - a known method
1883413f72f0SBarry Smith 
1884413f72f0SBarry Smith    Options Database Key:
1885cab54364SBarry Smith .  -islocaltoglobalmapping_type  <method> - Sets the method; use -help for a list of available methods (for instance, basic or hash)
1886cab54364SBarry Smith 
1887cab54364SBarry Smith   Level: intermediate
1888413f72f0SBarry Smith 
1889413f72f0SBarry Smith    Notes:
1890413f72f0SBarry Smith    See "petsc/include/petscis.h" for available methods
1891413f72f0SBarry Smith 
1892cab54364SBarry Smith   Normally, it is best to use the `ISLocalToGlobalMappingSetFromOptions()` command and
1893cab54364SBarry Smith   then set the `ISLocalToGlobalMappingType` from the options database rather than by using
1894413f72f0SBarry Smith   this routine.
1895413f72f0SBarry Smith 
1896cab54364SBarry Smith   Developer Note:
1897cab54364SBarry Smith   `ISLocalToGlobalMappingRegister()` is used to add new types to `ISLocalToGlobalMappingList` from which they
1898cab54364SBarry Smith   are accessed by `ISLocalToGlobalMappingSetType()`.
1899413f72f0SBarry Smith 
1900cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingRegister()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingGetType()`
1901413f72f0SBarry Smith @*/
1902d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType type)
1903d71ae5a4SJacob Faibussowitsch {
1904413f72f0SBarry Smith   PetscBool match;
19055f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(ISLocalToGlobalMapping) = NULL;
1906413f72f0SBarry Smith 
1907413f72f0SBarry Smith   PetscFunctionBegin;
1908413f72f0SBarry Smith   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
1909a0d79125SStefano Zampini   if (type) PetscValidCharPointer(type, 2);
1910413f72f0SBarry Smith 
19119566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ltog, type, &match));
1912413f72f0SBarry Smith   if (match) PetscFunctionReturn(0);
1913413f72f0SBarry Smith 
1914a0d79125SStefano Zampini   /* L2G maps defer type setup at globaltolocal calls, allow passing NULL here */
1915a0d79125SStefano Zampini   if (type) {
19169566063dSJacob Faibussowitsch     PetscCall(PetscFunctionListFind(ISLocalToGlobalMappingList, type, &r));
1917a0d79125SStefano Zampini     PetscCheck(r, PetscObjectComm((PetscObject)ltog), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested ISLocalToGlobalMapping type %s", type);
1918a0d79125SStefano Zampini   }
1919413f72f0SBarry Smith   /* Destroy the previous private LTOG context */
1920dbbe0bcdSBarry Smith   PetscTryTypeMethod(ltog, destroy);
1921413f72f0SBarry Smith   ltog->ops->destroy = NULL;
1922dbbe0bcdSBarry Smith 
19239566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)ltog, type));
19249566063dSJacob Faibussowitsch   if (r) PetscCall((*r)(ltog));
1925a0d79125SStefano Zampini   PetscFunctionReturn(0);
1926a0d79125SStefano Zampini }
1927a0d79125SStefano Zampini 
1928a0d79125SStefano Zampini /*@C
1929cab54364SBarry Smith    ISLocalToGlobalMappingGetType - Get the type of the `ISLocalToGlobalMapping`
1930a0d79125SStefano Zampini 
1931a0d79125SStefano Zampini    Not Collective
1932a0d79125SStefano Zampini 
1933a0d79125SStefano Zampini    Input Parameter:
1934cab54364SBarry Smith .  ltog - the `ISLocalToGlobalMapping` object
1935a0d79125SStefano Zampini 
1936a0d79125SStefano Zampini    Output Parameter:
1937a0d79125SStefano Zampini .  type - the type
1938a0d79125SStefano Zampini 
193949762cbcSSatish Balay    Level: intermediate
194049762cbcSSatish Balay 
1941cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingRegister()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`
1942a0d79125SStefano Zampini @*/
1943d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType *type)
1944d71ae5a4SJacob Faibussowitsch {
1945a0d79125SStefano Zampini   PetscFunctionBegin;
1946a0d79125SStefano Zampini   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1);
1947a0d79125SStefano Zampini   PetscValidPointer(type, 2);
1948a0d79125SStefano Zampini   *type = ((PetscObject)ltog)->type_name;
1949413f72f0SBarry Smith   PetscFunctionReturn(0);
1950413f72f0SBarry Smith }
1951413f72f0SBarry Smith 
1952413f72f0SBarry Smith PetscBool ISLocalToGlobalMappingRegisterAllCalled = PETSC_FALSE;
1953413f72f0SBarry Smith 
1954413f72f0SBarry Smith /*@C
1955cab54364SBarry Smith   ISLocalToGlobalMappingRegisterAll - Registers all of the local to global mapping components in the `IS` package.
1956413f72f0SBarry Smith 
1957413f72f0SBarry Smith   Not Collective
1958413f72f0SBarry Smith 
1959413f72f0SBarry Smith   Level: advanced
1960413f72f0SBarry Smith 
1961cab54364SBarry Smith .seealso: [](sec_scatter), `ISRegister()`, `ISLocalToGlobalRegister()`
1962413f72f0SBarry Smith @*/
1963d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRegisterAll(void)
1964d71ae5a4SJacob Faibussowitsch {
1965413f72f0SBarry Smith   PetscFunctionBegin;
1966413f72f0SBarry Smith   if (ISLocalToGlobalMappingRegisterAllCalled) PetscFunctionReturn(0);
1967413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAllCalled = PETSC_TRUE;
19689566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGBASIC, ISLocalToGlobalMappingCreate_Basic));
19699566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGHASH, ISLocalToGlobalMappingCreate_Hash));
1970413f72f0SBarry Smith   PetscFunctionReturn(0);
1971413f72f0SBarry Smith }
1972