xref: /petsc/src/vec/is/utils/isltog.c (revision dbbe0bcd3f3a8fbab5a45420dc06f8387e5764c6)
12362add9SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/isimpl.h>    /*I "petscis.h"  I*/
3e8f14785SLisandro Dalcin #include <petsc/private/hashmapi.h>
40c312b8eSJed Brown #include <petscsf.h>
5665c2dedSJed Brown #include <petscviewer.h>
62362add9SBarry Smith 
77087cfbeSBarry Smith PetscClassId IS_LTOGM_CLASSID;
8268a049cSStefano Zampini static PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping,PetscInt*,PetscInt**,PetscInt**,PetscInt***);
98e58c17dSMatthew Knepley 
10413f72f0SBarry Smith typedef struct {
11413f72f0SBarry Smith   PetscInt *globals;
12413f72f0SBarry Smith } ISLocalToGlobalMapping_Basic;
13413f72f0SBarry Smith 
14413f72f0SBarry Smith typedef struct {
15e8f14785SLisandro Dalcin   PetscHMapI globalht;
16413f72f0SBarry Smith } ISLocalToGlobalMapping_Hash;
17413f72f0SBarry Smith 
186528b96dSMatthew G. Knepley /*@C
196528b96dSMatthew G. Knepley   ISGetPointRange - Returns a description of the points in an IS suitable for traversal
20413f72f0SBarry Smith 
216528b96dSMatthew G. Knepley   Not collective
226528b96dSMatthew G. Knepley 
236528b96dSMatthew G. Knepley   Input Parameter:
246528b96dSMatthew G. Knepley . pointIS - The IS object
256528b96dSMatthew G. Knepley 
266528b96dSMatthew G. Knepley   Output Parameters:
276528b96dSMatthew G. Knepley + pStart - The first index, see notes
286528b96dSMatthew G. Knepley . pEnd   - One past the last index, see notes
296528b96dSMatthew G. Knepley - points - The indices, see notes
306528b96dSMatthew G. Knepley 
316528b96dSMatthew G. Knepley   Notes:
326528b96dSMatthew G. Knepley   If the IS contains contiguous indices in an ISSTRIDE, then the indices are contained in [pStart, pEnd) and points = NULL. Otherwise, pStart = 0, pEnd = numIndices, and points is an array of the indices. This supports the following pattern
336528b96dSMatthew G. Knepley $ ISGetPointRange(is, &pStart, &pEnd, &points);
346528b96dSMatthew G. Knepley $ for (p = pStart; p < pEnd; ++p) {
356528b96dSMatthew G. Knepley $   const PetscInt point = points ? points[p] : p;
366528b96dSMatthew G. Knepley $ }
376528b96dSMatthew G. Knepley $ ISRestorePointRange(is, &pstart, &pEnd, &points);
386528b96dSMatthew G. Knepley 
396528b96dSMatthew G. Knepley   Level: intermediate
406528b96dSMatthew G. Knepley 
41db781477SPatrick Sanan .seealso: `ISRestorePointRange()`, `ISGetPointSubrange()`, `ISGetIndices()`, `ISCreateStride()`
426528b96dSMatthew G. Knepley @*/
439305a4c7SMatthew G. Knepley PetscErrorCode ISGetPointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points)
449305a4c7SMatthew G. Knepley {
459305a4c7SMatthew G. Knepley   PetscInt       numCells, step = 1;
469305a4c7SMatthew G. Knepley   PetscBool      isStride;
479305a4c7SMatthew G. Knepley 
489305a4c7SMatthew G. Knepley   PetscFunctionBeginHot;
499305a4c7SMatthew G. Knepley   *pStart = 0;
509305a4c7SMatthew G. Knepley   *points = NULL;
519566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(pointIS, &numCells));
529566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject) pointIS, ISSTRIDE, &isStride));
539566063dSJacob Faibussowitsch   if (isStride) PetscCall(ISStrideGetInfo(pointIS, pStart, &step));
549305a4c7SMatthew G. Knepley   *pEnd   = *pStart + numCells;
559566063dSJacob Faibussowitsch   if (!isStride || step != 1) PetscCall(ISGetIndices(pointIS, points));
569305a4c7SMatthew G. Knepley   PetscFunctionReturn(0);
579305a4c7SMatthew G. Knepley }
589305a4c7SMatthew G. Knepley 
596528b96dSMatthew G. Knepley /*@C
606528b96dSMatthew G. Knepley   ISRestorePointRange - Destroys the traversal description
616528b96dSMatthew G. Knepley 
626528b96dSMatthew G. Knepley   Not collective
636528b96dSMatthew G. Knepley 
646528b96dSMatthew G. Knepley   Input Parameters:
656528b96dSMatthew G. Knepley + pointIS - The IS object
666528b96dSMatthew G. Knepley . pStart  - The first index, from ISGetPointRange()
676528b96dSMatthew G. Knepley . pEnd    - One past the last index, from ISGetPointRange()
686528b96dSMatthew G. Knepley - points  - The indices, from ISGetPointRange()
696528b96dSMatthew G. Knepley 
706528b96dSMatthew G. Knepley   Notes:
716528b96dSMatthew G. Knepley   If the IS contains contiguous indices in an ISSTRIDE, then the indices are contained in [pStart, pEnd) and points = NULL. Otherwise, pStart = 0, pEnd = numIndices, and points is an array of the indices. This supports the following pattern
726528b96dSMatthew G. Knepley $ ISGetPointRange(is, &pStart, &pEnd, &points);
736528b96dSMatthew G. Knepley $ for (p = pStart; p < pEnd; ++p) {
746528b96dSMatthew G. Knepley $   const PetscInt point = points ? points[p] : p;
756528b96dSMatthew G. Knepley $ }
766528b96dSMatthew G. Knepley $ ISRestorePointRange(is, &pstart, &pEnd, &points);
776528b96dSMatthew G. Knepley 
786528b96dSMatthew G. Knepley   Level: intermediate
796528b96dSMatthew G. Knepley 
80db781477SPatrick Sanan .seealso: `ISGetPointRange()`, `ISGetPointSubrange()`, `ISGetIndices()`, `ISCreateStride()`
816528b96dSMatthew G. Knepley @*/
829305a4c7SMatthew G. Knepley PetscErrorCode ISRestorePointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points)
839305a4c7SMatthew G. Knepley {
849305a4c7SMatthew G. Knepley   PetscInt       step = 1;
859305a4c7SMatthew G. Knepley   PetscBool      isStride;
869305a4c7SMatthew G. Knepley 
879305a4c7SMatthew G. Knepley   PetscFunctionBeginHot;
889566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject) pointIS, ISSTRIDE, &isStride));
899566063dSJacob Faibussowitsch   if (isStride) PetscCall(ISStrideGetInfo(pointIS, pStart, &step));
909566063dSJacob Faibussowitsch   if (!isStride || step != 1) PetscCall(ISGetIndices(pointIS, points));
919305a4c7SMatthew G. Knepley   PetscFunctionReturn(0);
929305a4c7SMatthew G. Knepley }
939305a4c7SMatthew G. Knepley 
946528b96dSMatthew G. Knepley /*@C
956528b96dSMatthew G. Knepley   ISGetPointSubrange - Configures the input IS to be a subrange for the traversal information given
966528b96dSMatthew G. Knepley 
976528b96dSMatthew G. Knepley   Not collective
986528b96dSMatthew G. Knepley 
996528b96dSMatthew G. Knepley   Input Parameters:
1006528b96dSMatthew G. Knepley + subpointIS - The IS object to be configured
1016528b96dSMatthew G. Knepley . pStar   t  - The first index of the subrange
1026528b96dSMatthew G. Knepley . pEnd       - One past the last index for the subrange
1036528b96dSMatthew G. Knepley - points     - The indices for the entire range, from ISGetPointRange()
1046528b96dSMatthew G. Knepley 
1056528b96dSMatthew G. Knepley   Output Parameters:
1066528b96dSMatthew G. Knepley . subpointIS - The IS object now configured to be a subrange
1076528b96dSMatthew G. Knepley 
1086528b96dSMatthew G. Knepley   Notes:
1096528b96dSMatthew G. Knepley   The input IS will now respond properly to calls to ISGetPointRange() and return the subrange.
1106528b96dSMatthew G. Knepley 
1116528b96dSMatthew G. Knepley   Level: intermediate
1126528b96dSMatthew G. Knepley 
113db781477SPatrick Sanan .seealso: `ISGetPointRange()`, `ISRestorePointRange()`, `ISGetIndices()`, `ISCreateStride()`
1146528b96dSMatthew G. Knepley @*/
1159305a4c7SMatthew G. Knepley PetscErrorCode ISGetPointSubrange(IS subpointIS, PetscInt pStart, PetscInt pEnd, const PetscInt *points)
1169305a4c7SMatthew G. Knepley {
1179305a4c7SMatthew G. Knepley   PetscFunctionBeginHot;
1189305a4c7SMatthew G. Knepley   if (points) {
1199566063dSJacob Faibussowitsch     PetscCall(ISSetType(subpointIS, ISGENERAL));
1209566063dSJacob Faibussowitsch     PetscCall(ISGeneralSetIndices(subpointIS, pEnd-pStart, &points[pStart], PETSC_USE_POINTER));
1219305a4c7SMatthew G. Knepley   } else {
1229566063dSJacob Faibussowitsch     PetscCall(ISSetType(subpointIS, ISSTRIDE));
1239566063dSJacob Faibussowitsch     PetscCall(ISStrideSetStride(subpointIS, pEnd-pStart, pStart, 1));
1249305a4c7SMatthew G. Knepley   }
1259305a4c7SMatthew G. Knepley   PetscFunctionReturn(0);
1269305a4c7SMatthew G. Knepley }
1279305a4c7SMatthew G. Knepley 
128413f72f0SBarry Smith /* -----------------------------------------------------------------------------------------*/
129413f72f0SBarry Smith 
130413f72f0SBarry Smith /*
131413f72f0SBarry Smith     Creates the global mapping information in the ISLocalToGlobalMapping structure
132413f72f0SBarry Smith 
133413f72f0SBarry Smith     If the user has not selected how to handle the global to local mapping then use HASH for "large" problems
134413f72f0SBarry Smith */
135413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp(ISLocalToGlobalMapping mapping)
136413f72f0SBarry Smith {
137413f72f0SBarry Smith   PetscInt       i,*idx = mapping->indices,n = mapping->n,end,start;
138413f72f0SBarry Smith 
139413f72f0SBarry Smith   PetscFunctionBegin;
140413f72f0SBarry Smith   if (mapping->data) PetscFunctionReturn(0);
141413f72f0SBarry Smith   end   = 0;
142413f72f0SBarry Smith   start = PETSC_MAX_INT;
143413f72f0SBarry Smith 
144413f72f0SBarry Smith   for (i=0; i<n; i++) {
145413f72f0SBarry Smith     if (idx[i] < 0) continue;
146413f72f0SBarry Smith     if (idx[i] < start) start = idx[i];
147413f72f0SBarry Smith     if (idx[i] > end)   end   = idx[i];
148413f72f0SBarry Smith   }
149413f72f0SBarry Smith   if (start > end) {start = 0; end = -1;}
150413f72f0SBarry Smith   mapping->globalstart = start;
151413f72f0SBarry Smith   mapping->globalend   = end;
152413f72f0SBarry Smith   if (!((PetscObject)mapping)->type_name) {
153413f72f0SBarry Smith     if ((end - start) > PetscMax(4*n,1000000)) {
1549566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingSetType(mapping,ISLOCALTOGLOBALMAPPINGHASH));
155413f72f0SBarry Smith     } else {
1569566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingSetType(mapping,ISLOCALTOGLOBALMAPPINGBASIC));
157413f72f0SBarry Smith     }
158413f72f0SBarry Smith   }
159*dbbe0bcdSBarry Smith   PetscTryTypeMethod(mapping,globaltolocalmappingsetup);
160413f72f0SBarry Smith   PetscFunctionReturn(0);
161413f72f0SBarry Smith }
162413f72f0SBarry Smith 
163413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp_Basic(ISLocalToGlobalMapping mapping)
164413f72f0SBarry Smith {
165413f72f0SBarry Smith   PetscInt                    i,*idx = mapping->indices,n = mapping->n,end,start,*globals;
166413f72f0SBarry Smith   ISLocalToGlobalMapping_Basic *map;
167413f72f0SBarry Smith 
168413f72f0SBarry Smith   PetscFunctionBegin;
169413f72f0SBarry Smith   start            = mapping->globalstart;
170413f72f0SBarry Smith   end              = mapping->globalend;
1719566063dSJacob Faibussowitsch   PetscCall(PetscNew(&map));
1729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(end-start+2,&globals));
173413f72f0SBarry Smith   map->globals     = globals;
174413f72f0SBarry Smith   for (i=0; i<end-start+1; i++) globals[i] = -1;
175413f72f0SBarry Smith   for (i=0; i<n; i++) {
176413f72f0SBarry Smith     if (idx[i] < 0) continue;
177413f72f0SBarry Smith     globals[idx[i] - start] = i;
178413f72f0SBarry Smith   }
179413f72f0SBarry Smith   mapping->data = (void*)map;
1809566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)mapping,(end-start+1)*sizeof(PetscInt)));
181413f72f0SBarry Smith   PetscFunctionReturn(0);
182413f72f0SBarry Smith }
183413f72f0SBarry Smith 
184413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp_Hash(ISLocalToGlobalMapping mapping)
185413f72f0SBarry Smith {
186413f72f0SBarry Smith   PetscInt                    i,*idx = mapping->indices,n = mapping->n;
187413f72f0SBarry Smith   ISLocalToGlobalMapping_Hash *map;
188413f72f0SBarry Smith 
189413f72f0SBarry Smith   PetscFunctionBegin;
1909566063dSJacob Faibussowitsch   PetscCall(PetscNew(&map));
1919566063dSJacob Faibussowitsch   PetscCall(PetscHMapICreate(&map->globalht));
192413f72f0SBarry Smith   for (i=0; i<n; i++) {
193413f72f0SBarry Smith     if (idx[i] < 0) continue;
1949566063dSJacob Faibussowitsch     PetscCall(PetscHMapISet(map->globalht,idx[i],i));
195413f72f0SBarry Smith   }
196413f72f0SBarry Smith   mapping->data = (void*)map;
1979566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)mapping,2*n*sizeof(PetscInt)));
198413f72f0SBarry Smith   PetscFunctionReturn(0);
199413f72f0SBarry Smith }
200413f72f0SBarry Smith 
201413f72f0SBarry Smith static PetscErrorCode ISLocalToGlobalMappingDestroy_Basic(ISLocalToGlobalMapping mapping)
202413f72f0SBarry Smith {
203413f72f0SBarry Smith   ISLocalToGlobalMapping_Basic *map  = (ISLocalToGlobalMapping_Basic *)mapping->data;
204413f72f0SBarry Smith 
205413f72f0SBarry Smith   PetscFunctionBegin;
206413f72f0SBarry Smith   if (!map) PetscFunctionReturn(0);
2079566063dSJacob Faibussowitsch   PetscCall(PetscFree(map->globals));
2089566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping->data));
209413f72f0SBarry Smith   PetscFunctionReturn(0);
210413f72f0SBarry Smith }
211413f72f0SBarry Smith 
212413f72f0SBarry Smith static PetscErrorCode ISLocalToGlobalMappingDestroy_Hash(ISLocalToGlobalMapping mapping)
213413f72f0SBarry Smith {
214413f72f0SBarry Smith   ISLocalToGlobalMapping_Hash *map  = (ISLocalToGlobalMapping_Hash*)mapping->data;
215413f72f0SBarry Smith 
216413f72f0SBarry Smith   PetscFunctionBegin;
217413f72f0SBarry Smith   if (!map) PetscFunctionReturn(0);
2189566063dSJacob Faibussowitsch   PetscCall(PetscHMapIDestroy(&map->globalht));
2199566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping->data));
220413f72f0SBarry Smith   PetscFunctionReturn(0);
221413f72f0SBarry Smith }
222413f72f0SBarry Smith 
223413f72f0SBarry Smith #define GTOLTYPE _Basic
224413f72f0SBarry Smith #define GTOLNAME _Basic
225541bf97eSAdrian Croucher #define GTOLBS mapping->bs
226413f72f0SBarry Smith #define GTOL(g, local) do {                  \
227413f72f0SBarry Smith     local = map->globals[g/bs - start];      \
2280040bde1SJunchao Zhang     if (local >= 0) local = bs*local + (g % bs); \
229413f72f0SBarry Smith   } while (0)
230541bf97eSAdrian Croucher 
231413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
232413f72f0SBarry Smith 
233413f72f0SBarry Smith #define GTOLTYPE _Basic
234413f72f0SBarry Smith #define GTOLNAME Block_Basic
235541bf97eSAdrian Croucher #define GTOLBS 1
236413f72f0SBarry Smith #define GTOL(g, local) do {                  \
237413f72f0SBarry Smith     local = map->globals[g - start];         \
238413f72f0SBarry Smith   } while (0)
239413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
240413f72f0SBarry Smith 
241413f72f0SBarry Smith #define GTOLTYPE _Hash
242413f72f0SBarry Smith #define GTOLNAME _Hash
243541bf97eSAdrian Croucher #define GTOLBS mapping->bs
244413f72f0SBarry Smith #define GTOL(g, local) do {                         \
245e8f14785SLisandro Dalcin     (void)PetscHMapIGet(map->globalht,g/bs,&local); \
2460040bde1SJunchao Zhang     if (local >= 0) local = bs*local + (g % bs);   \
247413f72f0SBarry Smith    } while (0)
248413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
249413f72f0SBarry Smith 
250413f72f0SBarry Smith #define GTOLTYPE _Hash
251413f72f0SBarry Smith #define GTOLNAME Block_Hash
252541bf97eSAdrian Croucher #define GTOLBS 1
253413f72f0SBarry Smith #define GTOL(g, local) do {                         \
254e8f14785SLisandro Dalcin     (void)PetscHMapIGet(map->globalht,g,&local);    \
255413f72f0SBarry Smith   } while (0)
256413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h>
257413f72f0SBarry Smith 
2586658fb44Sstefano_zampini /*@
2596658fb44Sstefano_zampini     ISLocalToGlobalMappingDuplicate - Duplicates the local to global mapping object
2606658fb44Sstefano_zampini 
2616658fb44Sstefano_zampini     Not Collective
2626658fb44Sstefano_zampini 
2636658fb44Sstefano_zampini     Input Parameter:
2646658fb44Sstefano_zampini .   ltog - local to global mapping
2656658fb44Sstefano_zampini 
2666658fb44Sstefano_zampini     Output Parameter:
2676658fb44Sstefano_zampini .   nltog - the duplicated local to global mapping
2686658fb44Sstefano_zampini 
2696658fb44Sstefano_zampini     Level: advanced
2706658fb44Sstefano_zampini 
271db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`
2726658fb44Sstefano_zampini @*/
2736658fb44Sstefano_zampini PetscErrorCode  ISLocalToGlobalMappingDuplicate(ISLocalToGlobalMapping ltog,ISLocalToGlobalMapping* nltog)
2746658fb44Sstefano_zampini {
275a0d79125SStefano Zampini   ISLocalToGlobalMappingType l2gtype;
2766658fb44Sstefano_zampini 
2776658fb44Sstefano_zampini   PetscFunctionBegin;
2786658fb44Sstefano_zampini   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
2799566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)ltog),ltog->bs,ltog->n,ltog->indices,PETSC_COPY_VALUES,nltog));
2809566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetType(ltog,&l2gtype));
2819566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingSetType(*nltog,l2gtype));
2826658fb44Sstefano_zampini   PetscFunctionReturn(0);
2836658fb44Sstefano_zampini }
2846658fb44Sstefano_zampini 
285565245c5SBarry Smith /*@
286107e9a97SBarry Smith     ISLocalToGlobalMappingGetSize - Gets the local size of a local to global mapping
2873b9aefa3SBarry Smith 
2883b9aefa3SBarry Smith     Not Collective
2893b9aefa3SBarry Smith 
2903b9aefa3SBarry Smith     Input Parameter:
2913b9aefa3SBarry Smith .   ltog - local to global mapping
2923b9aefa3SBarry Smith 
2933b9aefa3SBarry Smith     Output Parameter:
294107e9a97SBarry Smith .   n - the number of entries in the local mapping, ISLocalToGlobalMappingGetIndices() returns an array of this length
2953b9aefa3SBarry Smith 
2963b9aefa3SBarry Smith     Level: advanced
2973b9aefa3SBarry Smith 
298db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`
2993b9aefa3SBarry Smith @*/
3007087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping mapping,PetscInt *n)
3013b9aefa3SBarry Smith {
3023b9aefa3SBarry Smith   PetscFunctionBegin;
3030700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
3044482741eSBarry Smith   PetscValidIntPointer(n,2);
305107e9a97SBarry Smith   *n = mapping->bs*mapping->n;
3063b9aefa3SBarry Smith   PetscFunctionReturn(0);
3073b9aefa3SBarry Smith }
3083b9aefa3SBarry Smith 
3095a5d4f66SBarry Smith /*@C
310fe2efc57SMark    ISLocalToGlobalMappingViewFromOptions - View from Options
311fe2efc57SMark 
312fe2efc57SMark    Collective on ISLocalToGlobalMapping
313fe2efc57SMark 
314fe2efc57SMark    Input Parameters:
315fe2efc57SMark +  A - the local to global mapping object
316736c3998SJose E. Roman .  obj - Optional object
317736c3998SJose E. Roman -  name - command line option
318fe2efc57SMark 
319fe2efc57SMark    Level: intermediate
320db781477SPatrick Sanan .seealso: `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingView`, `PetscObjectViewFromOptions()`, `ISLocalToGlobalMappingCreate()`
321fe2efc57SMark @*/
322fe2efc57SMark PetscErrorCode  ISLocalToGlobalMappingViewFromOptions(ISLocalToGlobalMapping A,PetscObject obj,const char name[])
323fe2efc57SMark {
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 
341db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`
3425a5d4f66SBarry Smith @*/
3437087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingView(ISLocalToGlobalMapping mapping,PetscViewer viewer)
3445a5d4f66SBarry Smith {
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);
3513050cee2SBarry Smith   if (!viewer) {
3529566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mapping),&viewer));
3533050cee2SBarry Smith   }
3540700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
3555a5d4f66SBarry Smith 
3569566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mapping),&rank));
3579566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
35832077d6dSBarry Smith   if (iascii) {
3599566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)mapping,viewer));
3609566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
3615a5d4f66SBarry Smith     for (i=0; i<mapping->n; i++) {
3629566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %" PetscInt_FMT " %" PetscInt_FMT "\n",rank,i,mapping->indices[i]));
3636831982aSBarry Smith     }
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 
38295452b02SPatrick Sanan     Notes:
38395452b02SPatrick Sanan     the block size of the IS determines the block size of the mapping
384a997ad1aSLois Curfman McInnes     Level: advanced
385a997ad1aSLois Curfman McInnes 
386db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetFromOptions()`
3872bdab257SBarry Smith @*/
3887087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingCreateIS(IS is,ISLocalToGlobalMapping *mapping)
3892bdab257SBarry Smith {
3903bbf0e92SBarry Smith   PetscInt       n,bs;
3915d0c19d7SBarry Smith   const PetscInt *indices;
3922bdab257SBarry Smith   MPI_Comm       comm;
3933bbf0e92SBarry Smith   PetscBool      isblock;
3943a40ed3dSBarry Smith 
3953a40ed3dSBarry Smith   PetscFunctionBegin;
3960700a824SBarry Smith   PetscValidHeaderSpecific(is,IS_CLASSID,1);
3974482741eSBarry Smith   PetscValidPointer(mapping,2);
3982bdab257SBarry Smith 
3999566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)is,&comm));
4009566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is,&n));
4019566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)is,ISBLOCK,&isblock));
4026006e8d2SBarry Smith   if (!isblock) {
4039566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(is,&indices));
4049566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(comm,1,n,indices,PETSC_COPY_VALUES,mapping));
4059566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(is,&indices));
4066006e8d2SBarry Smith   } else {
4079566063dSJacob Faibussowitsch     PetscCall(ISGetBlockSize(is,&bs));
4089566063dSJacob Faibussowitsch     PetscCall(ISBlockGetIndices(is,&indices));
4099566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(comm,bs,n/bs,indices,PETSC_COPY_VALUES,mapping));
4109566063dSJacob Faibussowitsch     PetscCall(ISBlockRestoreIndices(is,&indices));
4116006e8d2SBarry Smith   }
4123a40ed3dSBarry Smith   PetscFunctionReturn(0);
4132bdab257SBarry Smith }
4145a5d4f66SBarry Smith 
415a4d96a55SJed Brown /*@C
416a4d96a55SJed Brown     ISLocalToGlobalMappingCreateSF - Creates a mapping between a local (0 to n)
417a4d96a55SJed Brown     ordering and a global parallel ordering.
418a4d96a55SJed Brown 
419a4d96a55SJed Brown     Collective
420a4d96a55SJed Brown 
421d8d19677SJose E. Roman     Input Parameters:
422a4d96a55SJed Brown +   sf - star forest mapping contiguous local indices to (rank, offset)
4239a535bafSVaclav Hapla -   start - first global index on this process, or PETSC_DECIDE to compute contiguous global numbering automatically
424a4d96a55SJed Brown 
425a4d96a55SJed Brown     Output Parameter:
426a4d96a55SJed Brown .   mapping - new mapping data structure
427a4d96a55SJed Brown 
428a4d96a55SJed Brown     Level: advanced
429a4d96a55SJed Brown 
4309a535bafSVaclav Hapla     Notes:
4319a535bafSVaclav Hapla     If any processor calls this with start = PETSC_DECIDE then all processors must, otherwise the program will hang.
4329a535bafSVaclav Hapla 
433db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingSetFromOptions()`
434a4d96a55SJed Brown @*/
435a4d96a55SJed Brown PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF sf,PetscInt start,ISLocalToGlobalMapping *mapping)
436a4d96a55SJed Brown {
437a4d96a55SJed Brown   PetscInt       i,maxlocal,nroots,nleaves,*globals,*ltog;
438a4d96a55SJed Brown   MPI_Comm       comm;
439a4d96a55SJed Brown 
440a4d96a55SJed Brown   PetscFunctionBegin;
441a4d96a55SJed Brown   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
442a4d96a55SJed Brown   PetscValidPointer(mapping,3);
4439566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sf,&comm));
44441f4c31fSVaclav Hapla   PetscCall(PetscSFGetGraph(sf,&nroots,&nleaves,NULL,NULL));
4459a535bafSVaclav Hapla   if (start == PETSC_DECIDE) {
4469a535bafSVaclav Hapla     start = 0;
4479566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Exscan(&nroots,&start,1,MPIU_INT,MPI_SUM,comm));
44841f4c31fSVaclav Hapla   } else PetscCheck(start >= 0,comm, PETSC_ERR_ARG_OUTOFRANGE, "start must be nonnegative or PETSC_DECIDE");
44941f4c31fSVaclav Hapla   PetscCall(PetscSFGetLeafRange(sf, NULL, &maxlocal));
45041f4c31fSVaclav Hapla   ++maxlocal;
4519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nroots,&globals));
4529566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxlocal,&ltog));
453a4d96a55SJed Brown   for (i=0; i<nroots; i++) globals[i] = start + i;
454a4d96a55SJed Brown   for (i=0; i<maxlocal; i++) ltog[i] = -1;
4559566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sf,MPIU_INT,globals,ltog,MPI_REPLACE));
4569566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sf,MPIU_INT,globals,ltog,MPI_REPLACE));
4579566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(comm,1,maxlocal,ltog,PETSC_OWN_POINTER,mapping));
4589566063dSJacob Faibussowitsch   PetscCall(PetscFree(globals));
459a4d96a55SJed Brown   PetscFunctionReturn(0);
460a4d96a55SJed Brown }
461b46b645bSBarry Smith 
46263fa5c83Sstefano_zampini /*@
46363fa5c83Sstefano_zampini     ISLocalToGlobalMappingSetBlockSize - Sets the blocksize of the mapping
46463fa5c83Sstefano_zampini 
46563fa5c83Sstefano_zampini     Not collective
46663fa5c83Sstefano_zampini 
46763fa5c83Sstefano_zampini     Input Parameters:
468a2b725a8SWilliam Gropp +   mapping - mapping data structure
469a2b725a8SWilliam Gropp -   bs - the blocksize
47063fa5c83Sstefano_zampini 
47163fa5c83Sstefano_zampini     Level: advanced
47263fa5c83Sstefano_zampini 
473db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`
47463fa5c83Sstefano_zampini @*/
47563fa5c83Sstefano_zampini PetscErrorCode  ISLocalToGlobalMappingSetBlockSize(ISLocalToGlobalMapping mapping,PetscInt bs)
47663fa5c83Sstefano_zampini {
477a59f3c4dSstefano_zampini   PetscInt       *nid;
478a59f3c4dSstefano_zampini   const PetscInt *oid;
479a59f3c4dSstefano_zampini   PetscInt       i,cn,on,obs,nn;
48063fa5c83Sstefano_zampini 
48163fa5c83Sstefano_zampini   PetscFunctionBegin;
48263fa5c83Sstefano_zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
48308401ef6SPierre Jolivet   PetscCheck(bs >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid block size %" PetscInt_FMT,bs);
48463fa5c83Sstefano_zampini   if (bs == mapping->bs) PetscFunctionReturn(0);
48563fa5c83Sstefano_zampini   on  = mapping->n;
48663fa5c83Sstefano_zampini   obs = mapping->bs;
48763fa5c83Sstefano_zampini   oid = mapping->indices;
48863fa5c83Sstefano_zampini   nn  = (on*obs)/bs;
48908401ef6SPierre 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);
490a59f3c4dSstefano_zampini 
4919566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nn,&nid));
4929566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetIndices(mapping,&oid));
493a59f3c4dSstefano_zampini   for (i=0;i<nn;i++) {
494a59f3c4dSstefano_zampini     PetscInt j;
495a59f3c4dSstefano_zampini     for (j=0,cn=0;j<bs-1;j++) {
496a59f3c4dSstefano_zampini       if (oid[i*bs+j] < 0) { cn++; continue; }
49708401ef6SPierre 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]);
498a59f3c4dSstefano_zampini     }
499a59f3c4dSstefano_zampini     if (oid[i*bs+j] < 0) cn++;
5008b7cb0e6Sstefano_zampini     if (cn) {
50108401ef6SPierre 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);
502a59f3c4dSstefano_zampini       nid[i] = -1;
5038b7cb0e6Sstefano_zampini     } else {
504a59f3c4dSstefano_zampini       nid[i] = oid[i*bs]/bs;
50563fa5c83Sstefano_zampini     }
50663fa5c83Sstefano_zampini   }
5079566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreIndices(mapping,&oid));
508a59f3c4dSstefano_zampini 
50963fa5c83Sstefano_zampini   mapping->n           = nn;
51063fa5c83Sstefano_zampini   mapping->bs          = bs;
5119566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping->indices));
51263fa5c83Sstefano_zampini   mapping->indices     = nid;
513c9345713Sstefano_zampini   mapping->globalstart = 0;
514c9345713Sstefano_zampini   mapping->globalend   = 0;
5151bd0b88eSStefano Zampini 
5161bd0b88eSStefano Zampini   /* reset the cached information */
5179566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping->info_procs));
5189566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping->info_numprocs));
5191bd0b88eSStefano Zampini   if (mapping->info_indices) {
5201bd0b88eSStefano Zampini     PetscInt i;
5211bd0b88eSStefano Zampini 
5229566063dSJacob Faibussowitsch     PetscCall(PetscFree((mapping->info_indices)[0]));
5231bd0b88eSStefano Zampini     for (i=1; i<mapping->info_nproc; i++) {
5249566063dSJacob Faibussowitsch       PetscCall(PetscFree(mapping->info_indices[i]));
5251bd0b88eSStefano Zampini     }
5269566063dSJacob Faibussowitsch     PetscCall(PetscFree(mapping->info_indices));
5271bd0b88eSStefano Zampini   }
5281bd0b88eSStefano Zampini   mapping->info_cached = PETSC_FALSE;
5291bd0b88eSStefano Zampini 
530*dbbe0bcdSBarry Smith   PetscTryTypeMethod(mapping,destroy);
53163fa5c83Sstefano_zampini   PetscFunctionReturn(0);
53263fa5c83Sstefano_zampini }
53363fa5c83Sstefano_zampini 
53445b6f7e9SBarry Smith /*@
53545b6f7e9SBarry Smith     ISLocalToGlobalMappingGetBlockSize - Gets the blocksize of the mapping
53645b6f7e9SBarry Smith     ordering and a global parallel ordering.
53745b6f7e9SBarry Smith 
53845b6f7e9SBarry Smith     Not Collective
53945b6f7e9SBarry Smith 
54045b6f7e9SBarry Smith     Input Parameters:
54145b6f7e9SBarry Smith .   mapping - mapping data structure
54245b6f7e9SBarry Smith 
54345b6f7e9SBarry Smith     Output Parameter:
54445b6f7e9SBarry Smith .   bs - the blocksize
54545b6f7e9SBarry Smith 
54645b6f7e9SBarry Smith     Level: advanced
54745b6f7e9SBarry Smith 
548db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`
54945b6f7e9SBarry Smith @*/
55045b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping mapping,PetscInt *bs)
55145b6f7e9SBarry Smith {
55245b6f7e9SBarry Smith   PetscFunctionBegin;
553cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
55445b6f7e9SBarry Smith   *bs = mapping->bs;
55545b6f7e9SBarry Smith   PetscFunctionReturn(0);
55645b6f7e9SBarry Smith }
55745b6f7e9SBarry Smith 
558ba5bb76aSSatish Balay /*@
55990f02eecSBarry Smith     ISLocalToGlobalMappingCreate - Creates a mapping between a local (0 to n)
56090f02eecSBarry Smith     ordering and a global parallel ordering.
5612362add9SBarry Smith 
56289d82c54SBarry Smith     Not Collective, but communicator may have more than one process
563b9cd556bSLois Curfman McInnes 
5642362add9SBarry Smith     Input Parameters:
56589d82c54SBarry Smith +   comm - MPI communicator
566f0413b6fSBarry Smith .   bs - the block size
56728bc9809SBarry Smith .   n - the number of local elements divided by the block size, or equivalently the number of block indices
56828bc9809SBarry 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
569d5ad8652SBarry Smith -   mode - see PetscCopyMode
5702362add9SBarry Smith 
571a997ad1aSLois Curfman McInnes     Output Parameter:
57290f02eecSBarry Smith .   mapping - new mapping data structure
5732362add9SBarry Smith 
57495452b02SPatrick Sanan     Notes:
57595452b02SPatrick Sanan     There is one integer value in indices per block and it represents the actual indices bs*idx + j, where j=0,..,bs-1
576413f72f0SBarry Smith 
5779a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
578413f72f0SBarry 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.
579413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
580413f72f0SBarry Smith 
581a997ad1aSLois Curfman McInnes     Level: advanced
582a997ad1aSLois Curfman McInnes 
583db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingSetFromOptions()`, `ISLOCALTOGLOBALMAPPINGBASIC`, `ISLOCALTOGLOBALMAPPINGHASH`
584db781477SPatrick Sanan           `ISLocalToGlobalMappingSetType()`, `ISLocalToGlobalMappingType`
5852362add9SBarry Smith @*/
58660c7cefcSBarry Smith PetscErrorCode  ISLocalToGlobalMappingCreate(MPI_Comm comm,PetscInt bs,PetscInt n,const PetscInt indices[],PetscCopyMode mode,ISLocalToGlobalMapping *mapping)
5872362add9SBarry Smith {
58832dcc486SBarry Smith   PetscInt       *in;
589b46b645bSBarry Smith 
590b46b645bSBarry Smith   PetscFunctionBegin;
591064a246eSJacob Faibussowitsch   if (n) PetscValidIntPointer(indices,4);
592064a246eSJacob Faibussowitsch   PetscValidPointer(mapping,6);
593b46b645bSBarry Smith 
5940298fd71SBarry Smith   *mapping = NULL;
5959566063dSJacob Faibussowitsch   PetscCall(ISInitializePackage());
5962362add9SBarry Smith 
5979566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(*mapping,IS_LTOGM_CLASSID,"ISLocalToGlobalMapping","Local to global mapping","IS",comm,ISLocalToGlobalMappingDestroy,ISLocalToGlobalMappingView));
598d4bb536fSBarry Smith   (*mapping)->n  = n;
599f0413b6fSBarry Smith   (*mapping)->bs = bs;
600d5ad8652SBarry Smith   if (mode == PETSC_COPY_VALUES) {
6019566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n,&in));
6029566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(in,indices,n));
603d5ad8652SBarry Smith     (*mapping)->indices = in;
60471910c26SVaclav Hapla     (*mapping)->dealloc_indices = PETSC_TRUE;
6059566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt)));
6066389a1a1SBarry Smith   } else if (mode == PETSC_OWN_POINTER) {
6076389a1a1SBarry Smith     (*mapping)->indices = (PetscInt*)indices;
60871910c26SVaclav Hapla     (*mapping)->dealloc_indices = PETSC_TRUE;
6099566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt)));
61071910c26SVaclav Hapla   } else if (mode == PETSC_USE_POINTER) {
61171910c26SVaclav Hapla     (*mapping)->indices = (PetscInt*)indices;
6126389a1a1SBarry Smith   }
61398921bdaSJacob Faibussowitsch   else SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid mode %d", mode);
6143a40ed3dSBarry Smith   PetscFunctionReturn(0);
6152362add9SBarry Smith }
6162362add9SBarry Smith 
617413f72f0SBarry Smith PetscFunctionList ISLocalToGlobalMappingList = NULL;
618413f72f0SBarry Smith 
61990f02eecSBarry Smith /*@
6207e99dc12SLawrence Mitchell    ISLocalToGlobalMappingSetFromOptions - Set mapping options from the options database.
6217e99dc12SLawrence Mitchell 
6227e99dc12SLawrence Mitchell    Not collective
6237e99dc12SLawrence Mitchell 
6247e99dc12SLawrence Mitchell    Input Parameters:
6257e99dc12SLawrence Mitchell .  mapping - mapping data structure
6267e99dc12SLawrence Mitchell 
6277e99dc12SLawrence Mitchell    Level: advanced
6287e99dc12SLawrence Mitchell 
6297e99dc12SLawrence Mitchell @*/
6307e99dc12SLawrence Mitchell PetscErrorCode ISLocalToGlobalMappingSetFromOptions(ISLocalToGlobalMapping mapping)
6317e99dc12SLawrence Mitchell {
632413f72f0SBarry Smith   char                       type[256];
633413f72f0SBarry Smith   ISLocalToGlobalMappingType defaulttype = "Not set";
6347e99dc12SLawrence Mitchell   PetscBool                  flg;
6357e99dc12SLawrence Mitchell 
6367e99dc12SLawrence Mitchell   PetscFunctionBegin;
6377e99dc12SLawrence Mitchell   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
6389566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRegisterAll());
639d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)mapping);
6409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-islocaltoglobalmapping_type","ISLocalToGlobalMapping method","ISLocalToGlobalMappingSetType",ISLocalToGlobalMappingList,(char*)(((PetscObject)mapping)->type_name) ? ((PetscObject)mapping)->type_name : defaulttype,type,256,&flg));
6411baa6e33SBarry Smith   if (flg) PetscCall(ISLocalToGlobalMappingSetType(mapping,type));
642d0609cedSBarry Smith   PetscOptionsEnd();
6437e99dc12SLawrence Mitchell   PetscFunctionReturn(0);
6447e99dc12SLawrence Mitchell }
6457e99dc12SLawrence Mitchell 
6467e99dc12SLawrence Mitchell /*@
64790f02eecSBarry Smith    ISLocalToGlobalMappingDestroy - Destroys a mapping between a local (0 to n)
64890f02eecSBarry Smith    ordering and a global parallel ordering.
64990f02eecSBarry Smith 
6500f5bd95cSBarry Smith    Note Collective
651b9cd556bSLois Curfman McInnes 
65290f02eecSBarry Smith    Input Parameters:
65390f02eecSBarry Smith .  mapping - mapping data structure
65490f02eecSBarry Smith 
655a997ad1aSLois Curfman McInnes    Level: advanced
656a997ad1aSLois Curfman McInnes 
657db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`
65890f02eecSBarry Smith @*/
6596bf464f9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping *mapping)
66090f02eecSBarry Smith {
6613a40ed3dSBarry Smith   PetscFunctionBegin;
6626bf464f9SBarry Smith   if (!*mapping) PetscFunctionReturn(0);
6636bf464f9SBarry Smith   PetscValidHeaderSpecific((*mapping),IS_LTOGM_CLASSID,1);
6644c8fdceaSLisandro Dalcin   if (--((PetscObject)(*mapping))->refct > 0) {*mapping = NULL;PetscFunctionReturn(0);}
66571910c26SVaclav Hapla   if ((*mapping)->dealloc_indices) {
6669566063dSJacob Faibussowitsch     PetscCall(PetscFree((*mapping)->indices));
66771910c26SVaclav Hapla   }
6689566063dSJacob Faibussowitsch   PetscCall(PetscFree((*mapping)->info_procs));
6699566063dSJacob Faibussowitsch   PetscCall(PetscFree((*mapping)->info_numprocs));
670268a049cSStefano Zampini   if ((*mapping)->info_indices) {
671268a049cSStefano Zampini     PetscInt i;
672268a049cSStefano Zampini 
6739566063dSJacob Faibussowitsch     PetscCall(PetscFree(((*mapping)->info_indices)[0]));
674268a049cSStefano Zampini     for (i=1; i<(*mapping)->info_nproc; i++) {
6759566063dSJacob Faibussowitsch       PetscCall(PetscFree(((*mapping)->info_indices)[i]));
676268a049cSStefano Zampini     }
6779566063dSJacob Faibussowitsch     PetscCall(PetscFree((*mapping)->info_indices));
678268a049cSStefano Zampini   }
6791bd0b88eSStefano Zampini   if ((*mapping)->info_nodei) {
6809566063dSJacob Faibussowitsch     PetscCall(PetscFree(((*mapping)->info_nodei)[0]));
6811bd0b88eSStefano Zampini   }
6829566063dSJacob Faibussowitsch   PetscCall(PetscFree2((*mapping)->info_nodec,(*mapping)->info_nodei));
683413f72f0SBarry Smith   if ((*mapping)->ops->destroy) {
6849566063dSJacob Faibussowitsch     PetscCall((*(*mapping)->ops->destroy)(*mapping));
685413f72f0SBarry Smith   }
6869566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(mapping));
6874c8fdceaSLisandro Dalcin   *mapping = NULL;
6883a40ed3dSBarry Smith   PetscFunctionReturn(0);
68990f02eecSBarry Smith }
69090f02eecSBarry Smith 
69190f02eecSBarry Smith /*@
6923acfe500SLois Curfman McInnes     ISLocalToGlobalMappingApplyIS - Creates from an IS in the local numbering
6933acfe500SLois Curfman McInnes     a new index set using the global numbering defined in an ISLocalToGlobalMapping
6943acfe500SLois Curfman McInnes     context.
69590f02eecSBarry Smith 
6964cb36875SStefano Zampini     Collective on is
697b9cd556bSLois Curfman McInnes 
69890f02eecSBarry Smith     Input Parameters:
699b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
700b9cd556bSLois Curfman McInnes -   is - index set in local numbering
70190f02eecSBarry Smith 
70290f02eecSBarry Smith     Output Parameters:
70390f02eecSBarry Smith .   newis - index set in global numbering
70490f02eecSBarry Smith 
7054cb36875SStefano Zampini     Notes:
7064cb36875SStefano Zampini     The output IS will have the same communicator of the input IS.
7074cb36875SStefano Zampini 
708a997ad1aSLois Curfman McInnes     Level: advanced
709a997ad1aSLois Curfman McInnes 
710db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingCreate()`,
711db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`, `ISGlobalToLocalMappingApply()`
71290f02eecSBarry Smith @*/
7137087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping,IS is,IS *newis)
71490f02eecSBarry Smith {
715e24637baSBarry Smith   PetscInt       n,*idxout;
7165d0c19d7SBarry Smith   const PetscInt *idxin;
7173a40ed3dSBarry Smith 
7183a40ed3dSBarry Smith   PetscFunctionBegin;
7190700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
7200700a824SBarry Smith   PetscValidHeaderSpecific(is,IS_CLASSID,2);
7214482741eSBarry Smith   PetscValidPointer(newis,3);
72290f02eecSBarry Smith 
7239566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is,&n));
7249566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is,&idxin));
7259566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n,&idxout));
7269566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(mapping,n,idxin,idxout));
7279566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is,&idxin));
7289566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)is),n,idxout,PETSC_OWN_POINTER,newis));
7293a40ed3dSBarry Smith   PetscFunctionReturn(0);
73090f02eecSBarry Smith }
73190f02eecSBarry Smith 
732b89cb25eSSatish Balay /*@
7333acfe500SLois Curfman McInnes    ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering
7343acfe500SLois Curfman McInnes    and converts them to the global numbering.
73590f02eecSBarry Smith 
736b9cd556bSLois Curfman McInnes    Not collective
737b9cd556bSLois Curfman McInnes 
738bb25748dSBarry Smith    Input Parameters:
739b9cd556bSLois Curfman McInnes +  mapping - the local to global mapping context
740bb25748dSBarry Smith .  N - number of integers
741b9cd556bSLois Curfman McInnes -  in - input indices in local numbering
742bb25748dSBarry Smith 
743bb25748dSBarry Smith    Output Parameter:
744bb25748dSBarry Smith .  out - indices in global numbering
745bb25748dSBarry Smith 
746b9cd556bSLois Curfman McInnes    Notes:
747b9cd556bSLois Curfman McInnes    The in and out array parameters may be identical.
748d4bb536fSBarry Smith 
749a997ad1aSLois Curfman McInnes    Level: advanced
750a997ad1aSLois Curfman McInnes 
751c2e3fba1SPatrick Sanan .seealso: `ISLocalToGlobalMappingApplyBlock()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingDestroy()`,
752c2e3fba1SPatrick Sanan           `ISLocalToGlobalMappingApplyIS()`, `AOCreateBasic()`, `AOApplicationToPetsc()`,
753db781477SPatrick Sanan           `AOPetscToApplication()`, `ISGlobalToLocalMappingApply()`
754bb25748dSBarry Smith 
755afcb2eb5SJed Brown @*/
756afcb2eb5SJed Brown PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
757afcb2eb5SJed Brown {
758cbc1caf0SMatthew G. Knepley   PetscInt i,bs,Nmax;
75945b6f7e9SBarry Smith 
76045b6f7e9SBarry Smith   PetscFunctionBegin;
761cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
762cbc1caf0SMatthew G. Knepley   bs   = mapping->bs;
763cbc1caf0SMatthew G. Knepley   Nmax = bs*mapping->n;
76445b6f7e9SBarry Smith   if (bs == 1) {
765cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
76645b6f7e9SBarry Smith     for (i=0; i<N; i++) {
76745b6f7e9SBarry Smith       if (in[i] < 0) {
76845b6f7e9SBarry Smith         out[i] = in[i];
76945b6f7e9SBarry Smith         continue;
77045b6f7e9SBarry Smith       }
77108401ef6SPierre 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);
77245b6f7e9SBarry Smith       out[i] = idx[in[i]];
77345b6f7e9SBarry Smith     }
77445b6f7e9SBarry Smith   } else {
775cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
77645b6f7e9SBarry Smith     for (i=0; i<N; i++) {
77745b6f7e9SBarry Smith       if (in[i] < 0) {
77845b6f7e9SBarry Smith         out[i] = in[i];
77945b6f7e9SBarry Smith         continue;
78045b6f7e9SBarry Smith       }
78108401ef6SPierre 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);
78245b6f7e9SBarry Smith       out[i] = idx[in[i]/bs]*bs + (in[i] % bs);
78345b6f7e9SBarry Smith     }
78445b6f7e9SBarry Smith   }
78545b6f7e9SBarry Smith   PetscFunctionReturn(0);
78645b6f7e9SBarry Smith }
78745b6f7e9SBarry Smith 
78845b6f7e9SBarry Smith /*@
7896006e8d2SBarry Smith    ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local block numbering and converts them to the global block numbering
79045b6f7e9SBarry Smith 
79145b6f7e9SBarry Smith    Not collective
79245b6f7e9SBarry Smith 
79345b6f7e9SBarry Smith    Input Parameters:
79445b6f7e9SBarry Smith +  mapping - the local to global mapping context
79545b6f7e9SBarry Smith .  N - number of integers
7966006e8d2SBarry Smith -  in - input indices in local block numbering
79745b6f7e9SBarry Smith 
79845b6f7e9SBarry Smith    Output Parameter:
7996006e8d2SBarry Smith .  out - indices in global block numbering
80045b6f7e9SBarry Smith 
80145b6f7e9SBarry Smith    Notes:
80245b6f7e9SBarry Smith    The in and out array parameters may be identical.
80345b6f7e9SBarry Smith 
8046006e8d2SBarry Smith    Example:
8056006e8d2SBarry 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
8066006e8d2SBarry Smith      (the first block) would produce 0 and the mapping applied to 1 (the second block) would produce 3.
8076006e8d2SBarry Smith 
80845b6f7e9SBarry Smith    Level: advanced
80945b6f7e9SBarry Smith 
810c2e3fba1SPatrick Sanan .seealso: `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingDestroy()`,
811c2e3fba1SPatrick Sanan           `ISLocalToGlobalMappingApplyIS()`, `AOCreateBasic()`, `AOApplicationToPetsc()`,
812db781477SPatrick Sanan           `AOPetscToApplication()`, `ISGlobalToLocalMappingApply()`
81345b6f7e9SBarry Smith 
81445b6f7e9SBarry Smith @*/
81545b6f7e9SBarry Smith PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
81645b6f7e9SBarry Smith {
8178a1f772fSStefano Zampini   PetscInt       i,Nmax;
8188a1f772fSStefano Zampini   const PetscInt *idx;
819d4bb536fSBarry Smith 
820a0d79125SStefano Zampini   PetscFunctionBegin;
821a0d79125SStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
8228a1f772fSStefano Zampini   Nmax = mapping->n;
8238a1f772fSStefano Zampini   idx = mapping->indices;
824afcb2eb5SJed Brown   for (i=0; i<N; i++) {
825afcb2eb5SJed Brown     if (in[i] < 0) {
826afcb2eb5SJed Brown       out[i] = in[i];
827afcb2eb5SJed Brown       continue;
828afcb2eb5SJed Brown     }
82908401ef6SPierre 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);
830afcb2eb5SJed Brown     out[i] = idx[in[i]];
831afcb2eb5SJed Brown   }
832afcb2eb5SJed Brown   PetscFunctionReturn(0);
833afcb2eb5SJed Brown }
834d4bb536fSBarry Smith 
8357e99dc12SLawrence Mitchell /*@
836a997ad1aSLois Curfman McInnes     ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers
837a997ad1aSLois Curfman McInnes     specified with a global numbering.
838d4bb536fSBarry Smith 
839b9cd556bSLois Curfman McInnes     Not collective
840b9cd556bSLois Curfman McInnes 
841d4bb536fSBarry Smith     Input Parameters:
842b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
8430040bde1SJunchao Zhang .   type - IS_GTOLM_MASK - maps global indices with no local value to -1 in the output list (i.e., mask them)
844d4bb536fSBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
845d4bb536fSBarry Smith .   n - number of global indices to map
846b9cd556bSLois Curfman McInnes -   idx - global indices to map
847d4bb536fSBarry Smith 
848d4bb536fSBarry Smith     Output Parameters:
849b9cd556bSLois Curfman McInnes +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
850b9cd556bSLois Curfman McInnes -   idxout - local index of each global index, one must pass in an array long enough
851e182c471SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApply() with
8520298fd71SBarry Smith              idxout == NULL to determine the required length (returned in nout)
853e182c471SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApply()
854e182c471SBarry Smith              a second time to set the values.
855d4bb536fSBarry Smith 
856b9cd556bSLois Curfman McInnes     Notes:
8570298fd71SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
858d4bb536fSBarry Smith 
8599a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
860413f72f0SBarry 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.
861413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
8620f5bd95cSBarry Smith 
863a997ad1aSLois Curfman McInnes     Level: advanced
864a997ad1aSLois Curfman McInnes 
86532fd6b96SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
86632fd6b96SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
86732fd6b96SBarry Smith 
868db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingApply()`, `ISGlobalToLocalMappingApplyBlock()`, `ISLocalToGlobalMappingCreate()`,
869db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`
870d4bb536fSBarry Smith @*/
871413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
872d4bb536fSBarry Smith {
8739d90f715SBarry Smith   PetscFunctionBegin;
8749d90f715SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
875413f72f0SBarry Smith   if (!mapping->data) {
8769566063dSJacob Faibussowitsch     PetscCall(ISGlobalToLocalMappingSetUp(mapping));
8779d90f715SBarry Smith   }
878*dbbe0bcdSBarry Smith   PetscUseTypeMethod(mapping,globaltolocalmappingapply ,type,n,idx,nout,idxout);
8799d90f715SBarry Smith   PetscFunctionReturn(0);
8809d90f715SBarry Smith }
8819d90f715SBarry Smith 
882d4fe737eSStefano Zampini /*@
883d4fe737eSStefano Zampini     ISGlobalToLocalMappingApplyIS - Creates from an IS in the global numbering
884d4fe737eSStefano Zampini     a new index set using the local numbering defined in an ISLocalToGlobalMapping
885d4fe737eSStefano Zampini     context.
886d4fe737eSStefano Zampini 
887d4fe737eSStefano Zampini     Not collective
888d4fe737eSStefano Zampini 
889d4fe737eSStefano Zampini     Input Parameters:
890d4fe737eSStefano Zampini +   mapping - mapping between local and global numbering
8910040bde1SJunchao Zhang .   type - IS_GTOLM_MASK - maps global indices with no local value to -1 in the output list (i.e., mask them)
8922785b321SStefano Zampini            IS_GTOLM_DROP - drops the indices with no local value from the output list
893d4fe737eSStefano Zampini -   is - index set in global numbering
894d4fe737eSStefano Zampini 
895d4fe737eSStefano Zampini     Output Parameters:
896d4fe737eSStefano Zampini .   newis - index set in local numbering
897d4fe737eSStefano Zampini 
8984cb36875SStefano Zampini     Notes:
8994cb36875SStefano Zampini     The output IS will be sequential, as it encodes a purely local operation
9004cb36875SStefano Zampini 
901d4fe737eSStefano Zampini     Level: advanced
902d4fe737eSStefano Zampini 
903db781477SPatrick Sanan .seealso: `ISGlobalToLocalMappingApply()`, `ISLocalToGlobalMappingCreate()`,
904db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`
905d4fe737eSStefano Zampini @*/
906413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,IS is,IS *newis)
907d4fe737eSStefano Zampini {
908d4fe737eSStefano Zampini   PetscInt       n,nout,*idxout;
909d4fe737eSStefano Zampini   const PetscInt *idxin;
910d4fe737eSStefano Zampini 
911d4fe737eSStefano Zampini   PetscFunctionBegin;
912d4fe737eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
913d4fe737eSStefano Zampini   PetscValidHeaderSpecific(is,IS_CLASSID,3);
914d4fe737eSStefano Zampini   PetscValidPointer(newis,4);
915d4fe737eSStefano Zampini 
9169566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is,&n));
9179566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is,&idxin));
918d4fe737eSStefano Zampini   if (type == IS_GTOLM_MASK) {
9199566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n,&idxout));
920d4fe737eSStefano Zampini   } else {
9219566063dSJacob Faibussowitsch     PetscCall(ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,NULL));
9229566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nout,&idxout));
923d4fe737eSStefano Zampini   }
9249566063dSJacob Faibussowitsch   PetscCall(ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,idxout));
9259566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is,&idxin));
9269566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,nout,idxout,PETSC_OWN_POINTER,newis));
927d4fe737eSStefano Zampini   PetscFunctionReturn(0);
928d4fe737eSStefano Zampini }
929d4fe737eSStefano Zampini 
9309d90f715SBarry Smith /*@
9319d90f715SBarry Smith     ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers
9329d90f715SBarry Smith     specified with a block global numbering.
9339d90f715SBarry Smith 
9349d90f715SBarry Smith     Not collective
9359d90f715SBarry Smith 
9369d90f715SBarry Smith     Input Parameters:
9379d90f715SBarry Smith +   mapping - mapping between local and global numbering
9380040bde1SJunchao Zhang .   type - IS_GTOLM_MASK - maps global indices with no local value to -1 in the output list (i.e., mask them)
9399d90f715SBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
9409d90f715SBarry Smith .   n - number of global indices to map
9419d90f715SBarry Smith -   idx - global indices to map
9429d90f715SBarry Smith 
9439d90f715SBarry Smith     Output Parameters:
9449d90f715SBarry Smith +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
9459d90f715SBarry Smith -   idxout - local index of each global index, one must pass in an array long enough
9469d90f715SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApplyBlock() with
9479d90f715SBarry Smith              idxout == NULL to determine the required length (returned in nout)
9489d90f715SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApplyBlock()
9499d90f715SBarry Smith              a second time to set the values.
9509d90f715SBarry Smith 
9519d90f715SBarry Smith     Notes:
9529d90f715SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
9539d90f715SBarry Smith 
9549a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
955413f72f0SBarry 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.
956413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
9579d90f715SBarry Smith 
9589d90f715SBarry Smith     Level: advanced
9599d90f715SBarry Smith 
9609d90f715SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
9619d90f715SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
9629d90f715SBarry Smith 
963db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingApply()`, `ISGlobalToLocalMappingApply()`, `ISLocalToGlobalMappingCreate()`,
964db781477SPatrick Sanan           `ISLocalToGlobalMappingDestroy()`
9659d90f715SBarry Smith @*/
966413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,
9679d90f715SBarry Smith                                                  PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
9689d90f715SBarry Smith {
9693a40ed3dSBarry Smith   PetscFunctionBegin;
9700700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
971413f72f0SBarry Smith   if (!mapping->data) {
9729566063dSJacob Faibussowitsch     PetscCall(ISGlobalToLocalMappingSetUp(mapping));
973d4bb536fSBarry Smith   }
974*dbbe0bcdSBarry Smith   PetscUseTypeMethod(mapping,globaltolocalmappingapplyblock ,type,n,idx,nout,idxout);
9753a40ed3dSBarry Smith   PetscFunctionReturn(0);
976d4bb536fSBarry Smith }
97790f02eecSBarry Smith 
97889d82c54SBarry Smith /*@C
9796a818285SBarry Smith     ISLocalToGlobalMappingGetBlockInfo - Gets the neighbor information for each processor and
98089d82c54SBarry Smith      each index shared by more than one processor
98189d82c54SBarry Smith 
98289d82c54SBarry Smith     Collective on ISLocalToGlobalMapping
98389d82c54SBarry Smith 
984f899ff85SJose E. Roman     Input Parameter:
98589d82c54SBarry Smith .   mapping - the mapping from local to global indexing
98689d82c54SBarry Smith 
987d8d19677SJose E. Roman     Output Parameters:
98889d82c54SBarry Smith +   nproc - number of processors that are connected to this one
98989d82c54SBarry Smith .   proc - neighboring processors
99007b52d57SBarry Smith .   numproc - number of indices for each subdomain (processor)
9913463a7baSJed Brown -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
99289d82c54SBarry Smith 
99389d82c54SBarry Smith     Level: advanced
99489d82c54SBarry Smith 
9952cfcea29SBarry Smith     Fortran Usage:
9962cfcea29SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
9972cfcea29SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
9982cfcea29SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
9992cfcea29SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
10002cfcea29SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
10012cfcea29SBarry Smith 
1002db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1003db781477SPatrick Sanan           `ISLocalToGlobalMappingRestoreInfo()`
100489d82c54SBarry Smith @*/
10056a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
100689d82c54SBarry Smith {
1007268a049cSStefano Zampini   PetscFunctionBegin;
1008268a049cSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1009268a049cSStefano Zampini   if (mapping->info_cached) {
1010268a049cSStefano Zampini     *nproc    = mapping->info_nproc;
1011268a049cSStefano Zampini     *procs    = mapping->info_procs;
1012268a049cSStefano Zampini     *numprocs = mapping->info_numprocs;
1013268a049cSStefano Zampini     *indices  = mapping->info_indices;
1014268a049cSStefano Zampini   } else {
10159566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetBlockInfo_Private(mapping,nproc,procs,numprocs,indices));
1016268a049cSStefano Zampini   }
1017268a049cSStefano Zampini   PetscFunctionReturn(0);
1018268a049cSStefano Zampini }
1019268a049cSStefano Zampini 
1020268a049cSStefano Zampini static PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
1021268a049cSStefano Zampini {
102297f1f81fSBarry Smith   PetscMPIInt    size,rank,tag1,tag2,tag3,*len,*source,imdex;
102332dcc486SBarry Smith   PetscInt       i,n = mapping->n,Ng,ng,max = 0,*lindices = mapping->indices;
102432dcc486SBarry Smith   PetscInt       *nprocs,*owner,nsends,*sends,j,*starts,nmax,nrecvs,*recvs,proc;
1025c599c493SJunchao Zhang   PetscInt       cnt,scale,*ownedsenders,*nownedsenders,rstart;
102632dcc486SBarry Smith   PetscInt       node,nownedm,nt,*sends2,nsends2,*starts2,*lens2,*dest,nrecvs2,*starts3,*recvs2,k,*bprocs,*tmp;
102732dcc486SBarry Smith   PetscInt       first_procs,first_numprocs,*first_indices;
102889d82c54SBarry Smith   MPI_Request    *recv_waits,*send_waits;
102930dcb7c9SBarry Smith   MPI_Status     recv_status,*send_status,*recv_statuses;
1030ce94432eSBarry Smith   MPI_Comm       comm;
1031ace3abfcSBarry Smith   PetscBool      debug = PETSC_FALSE;
103289d82c54SBarry Smith 
103389d82c54SBarry Smith   PetscFunctionBegin;
10349566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mapping,&comm));
10359566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
10369566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
103724cf384cSBarry Smith   if (size == 1) {
103824cf384cSBarry Smith     *nproc         = 0;
10390298fd71SBarry Smith     *procs         = NULL;
10409566063dSJacob Faibussowitsch     PetscCall(PetscNew(numprocs));
10411e2105dcSBarry Smith     (*numprocs)[0] = 0;
10429566063dSJacob Faibussowitsch     PetscCall(PetscNew(indices));
10430298fd71SBarry Smith     (*indices)[0]  = NULL;
1044268a049cSStefano Zampini     /* save info for reuse */
1045268a049cSStefano Zampini     mapping->info_nproc = *nproc;
1046268a049cSStefano Zampini     mapping->info_procs = *procs;
1047268a049cSStefano Zampini     mapping->info_numprocs = *numprocs;
1048268a049cSStefano Zampini     mapping->info_indices = *indices;
1049268a049cSStefano Zampini     mapping->info_cached = PETSC_TRUE;
105024cf384cSBarry Smith     PetscFunctionReturn(0);
105124cf384cSBarry Smith   }
105224cf384cSBarry Smith 
10539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)mapping)->options,NULL,"-islocaltoglobalmappinggetinfo_debug",&debug,NULL));
105407b52d57SBarry Smith 
10553677ff5aSBarry Smith   /*
10566a818285SBarry Smith     Notes on ISLocalToGlobalMappingGetBlockInfo
10573677ff5aSBarry Smith 
10583677ff5aSBarry Smith     globally owned node - the nodes that have been assigned to this processor in global
10593677ff5aSBarry Smith            numbering, just for this routine.
10603677ff5aSBarry Smith 
10613677ff5aSBarry Smith     nontrivial globally owned node - node assigned to this processor that is on a subdomain
10623677ff5aSBarry Smith            boundary (i.e. is has more than one local owner)
10633677ff5aSBarry Smith 
10643677ff5aSBarry Smith     locally owned node - node that exists on this processors subdomain
10653677ff5aSBarry Smith 
10663677ff5aSBarry Smith     nontrivial locally owned node - node that is not in the interior (i.e. has more than one
10673677ff5aSBarry Smith            local subdomain
10683677ff5aSBarry Smith   */
10699566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetNewTag((PetscObject)mapping,&tag1));
10709566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetNewTag((PetscObject)mapping,&tag2));
10719566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetNewTag((PetscObject)mapping,&tag3));
107289d82c54SBarry Smith 
107389d82c54SBarry Smith   for (i=0; i<n; i++) {
107489d82c54SBarry Smith     if (lindices[i] > max) max = lindices[i];
107589d82c54SBarry Smith   }
10761c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&max,&Ng,1,MPIU_INT,MPI_MAX,comm));
107778058e43SBarry Smith   Ng++;
10789566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
10799566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
1080bc8ff85bSBarry Smith   scale  = Ng/size + 1;
1081a2e34c3dSBarry Smith   ng     = scale; if (rank == size-1) ng = Ng - scale*(size-1); ng = PetscMax(1,ng);
1082caba0dd0SBarry Smith   rstart = scale*rank;
108389d82c54SBarry Smith 
108489d82c54SBarry Smith   /* determine ownership ranges of global indices */
10859566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2*size,&nprocs));
10869566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(nprocs,2*size));
108789d82c54SBarry Smith 
108889d82c54SBarry Smith   /* determine owners of each local node  */
10899566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n,&owner));
109089d82c54SBarry Smith   for (i=0; i<n; i++) {
10913677ff5aSBarry Smith     proc             = lindices[i]/scale; /* processor that globally owns this index */
109227c402fcSBarry Smith     nprocs[2*proc+1] = 1;                 /* processor globally owns at least one of ours */
10933677ff5aSBarry Smith     owner[i]         = proc;
109427c402fcSBarry Smith     nprocs[2*proc]++;                     /* count of how many that processor globally owns of ours */
109589d82c54SBarry Smith   }
109627c402fcSBarry Smith   nsends = 0; 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));
110689d82c54SBarry Smith   for (i=0; i<nrecvs; i++) {
11079566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Irecv(recvs+2*nmax*i,2*nmax,MPIU_INT,MPI_ANY_SOURCE,tag1,comm,recv_waits+i));
110889d82c54SBarry Smith   }
110989d82c54SBarry Smith 
111089d82c54SBarry Smith   /* pack messages containing lists of local nodes to owners */
11119566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2*n+1,&sends));
11129566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size+1,&starts));
111389d82c54SBarry Smith   starts[0] = 0;
1114f6e5521dSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
111589d82c54SBarry Smith   for (i=0; i<n; i++) {
111689d82c54SBarry Smith     sends[starts[owner[i]]++] = lindices[i];
111730dcb7c9SBarry Smith     sends[starts[owner[i]]++] = i;
111889d82c54SBarry Smith   }
11199566063dSJacob Faibussowitsch   PetscCall(PetscFree(owner));
112089d82c54SBarry Smith   starts[0] = 0;
1121f6e5521dSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
112289d82c54SBarry Smith 
112389d82c54SBarry Smith   /* send the messages */
11249566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends+1,&send_waits));
11259566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends+1,&dest));
112689d82c54SBarry Smith   cnt = 0;
112789d82c54SBarry Smith   for (i=0; i<size; i++) {
112827c402fcSBarry Smith     if (nprocs[2*i]) {
11299566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Isend(sends+starts[i],2*nprocs[2*i],MPIU_INT,i,tag1,comm,send_waits+cnt));
113030dcb7c9SBarry Smith       dest[cnt] = i;
113189d82c54SBarry Smith       cnt++;
113289d82c54SBarry Smith     }
113389d82c54SBarry Smith   }
11349566063dSJacob Faibussowitsch   PetscCall(PetscFree(starts));
113589d82c54SBarry Smith 
113689d82c54SBarry Smith   /* wait on receives */
11379566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs+1,&source));
11389566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs+1,&len));
113989d82c54SBarry Smith   cnt  = nrecvs;
11409566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(ng+1,&nownedsenders));
114189d82c54SBarry Smith   while (cnt) {
11429566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status));
114389d82c54SBarry Smith     /* unpack receives into our local space */
11449566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Get_count(&recv_status,MPIU_INT,&len[imdex]));
114589d82c54SBarry Smith     source[imdex] = recv_status.MPI_SOURCE;
114630dcb7c9SBarry Smith     len[imdex]    = len[imdex]/2;
1147caba0dd0SBarry Smith     /* count how many local owners for each of my global owned indices */
114830dcb7c9SBarry Smith     for (i=0; i<len[imdex]; i++) nownedsenders[recvs[2*imdex*nmax+2*i]-rstart]++;
114989d82c54SBarry Smith     cnt--;
115089d82c54SBarry Smith   }
11519566063dSJacob Faibussowitsch   PetscCall(PetscFree(recv_waits));
115289d82c54SBarry Smith 
115330dcb7c9SBarry Smith   /* count how many globally owned indices are on an edge multiplied by how many processors own them. */
1154bc8ff85bSBarry Smith   nownedm = 0;
1155bc8ff85bSBarry Smith   for (i=0; i<ng; i++) {
1156c599c493SJunchao Zhang     if (nownedsenders[i] > 1) nownedm += nownedsenders[i];
1157bc8ff85bSBarry Smith   }
1158bc8ff85bSBarry Smith 
1159bc8ff85bSBarry Smith   /* create single array to contain rank of all local owners of each globally owned index */
11609566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nownedm+1,&ownedsenders));
11619566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ng+1,&starts));
1162bc8ff85bSBarry Smith   starts[0] = 0;
1163bc8ff85bSBarry Smith   for (i=1; i<ng; i++) {
1164bc8ff85bSBarry Smith     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
1165bc8ff85bSBarry Smith     else starts[i] = starts[i-1];
1166bc8ff85bSBarry Smith   }
1167bc8ff85bSBarry Smith 
11686aad120cSJose E. Roman   /* for each nontrivial globally owned node list all arriving processors */
1169bc8ff85bSBarry Smith   for (i=0; i<nrecvs; i++) {
1170bc8ff85bSBarry Smith     for (j=0; j<len[i]; j++) {
117130dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
1172f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) ownedsenders[starts[node]++] = source[i];
1173bc8ff85bSBarry Smith     }
1174bc8ff85bSBarry Smith   }
1175bc8ff85bSBarry Smith 
117607b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
117730dcb7c9SBarry Smith     starts[0] = 0;
117830dcb7c9SBarry Smith     for (i=1; i<ng; i++) {
117930dcb7c9SBarry Smith       if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
118030dcb7c9SBarry Smith       else starts[i] = starts[i-1];
118130dcb7c9SBarry Smith     }
118230dcb7c9SBarry Smith     for (i=0; i<ng; i++) {
118330dcb7c9SBarry Smith       if (nownedsenders[i] > 1) {
11849566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(comm,"[%d] global node %" PetscInt_FMT " local owner processors: ",rank,i+rstart));
118530dcb7c9SBarry Smith         for (j=0; j<nownedsenders[i]; j++) {
11869566063dSJacob Faibussowitsch           PetscCall(PetscSynchronizedPrintf(comm,"%" PetscInt_FMT " ",ownedsenders[starts[i]+j]));
118730dcb7c9SBarry Smith         }
11889566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(comm,"\n"));
118930dcb7c9SBarry Smith       }
119030dcb7c9SBarry Smith     }
11919566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm,PETSC_STDOUT));
119207b52d57SBarry Smith   } /* -----------------------------------  */
119330dcb7c9SBarry Smith 
11943677ff5aSBarry Smith   /* wait on original sends */
11953a96401aSBarry Smith   if (nsends) {
11969566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nsends,&send_status));
11979566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nsends,send_waits,send_status));
11989566063dSJacob Faibussowitsch     PetscCall(PetscFree(send_status));
11993a96401aSBarry Smith   }
12009566063dSJacob Faibussowitsch   PetscCall(PetscFree(send_waits));
12019566063dSJacob Faibussowitsch   PetscCall(PetscFree(sends));
12029566063dSJacob Faibussowitsch   PetscCall(PetscFree(nprocs));
12033677ff5aSBarry Smith 
12043677ff5aSBarry Smith   /* pack messages to send back to local owners */
120530dcb7c9SBarry Smith   starts[0] = 0;
120630dcb7c9SBarry Smith   for (i=1; i<ng; i++) {
120730dcb7c9SBarry Smith     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
120830dcb7c9SBarry Smith     else starts[i] = starts[i-1];
120930dcb7c9SBarry Smith   }
121030dcb7c9SBarry Smith   nsends2 = nrecvs;
12119566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends2+1,&nprocs)); /* length of each message */
121230dcb7c9SBarry Smith   for (i=0; i<nrecvs; i++) {
121330dcb7c9SBarry Smith     nprocs[i] = 1;
121430dcb7c9SBarry Smith     for (j=0; j<len[i]; j++) {
121530dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
1216f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) nprocs[i] += 2 + nownedsenders[node];
121730dcb7c9SBarry Smith     }
121830dcb7c9SBarry Smith   }
1219f6e5521dSKarl Rupp   nt = 0;
1220f6e5521dSKarl Rupp   for (i=0; i<nsends2; i++) nt += nprocs[i];
1221f6e5521dSKarl Rupp 
12229566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt+1,&sends2));
12239566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends2+1,&starts2));
1224f6e5521dSKarl Rupp 
1225f6e5521dSKarl Rupp   starts2[0] = 0;
1226f6e5521dSKarl Rupp   for (i=1; i<nsends2; i++) starts2[i] = starts2[i-1] + nprocs[i-1];
122730dcb7c9SBarry Smith   /*
122830dcb7c9SBarry Smith      Each message is 1 + nprocs[i] long, and consists of
122930dcb7c9SBarry Smith        (0) the number of nodes being sent back
123030dcb7c9SBarry Smith        (1) the local node number,
123130dcb7c9SBarry Smith        (2) the number of processors sharing it,
123230dcb7c9SBarry Smith        (3) the processors sharing it
123330dcb7c9SBarry Smith   */
123430dcb7c9SBarry Smith   for (i=0; i<nsends2; i++) {
123530dcb7c9SBarry Smith     cnt = 1;
123630dcb7c9SBarry Smith     sends2[starts2[i]] = 0;
123730dcb7c9SBarry Smith     for (j=0; j<len[i]; j++) {
123830dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
123930dcb7c9SBarry Smith       if (nownedsenders[node] > 1) {
124030dcb7c9SBarry Smith         sends2[starts2[i]]++;
124130dcb7c9SBarry Smith         sends2[starts2[i]+cnt++] = recvs[2*i*nmax+2*j+1];
124230dcb7c9SBarry Smith         sends2[starts2[i]+cnt++] = nownedsenders[node];
12439566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(&sends2[starts2[i]+cnt],&ownedsenders[starts[node]],nownedsenders[node]));
124430dcb7c9SBarry Smith         cnt += nownedsenders[node];
124530dcb7c9SBarry Smith       }
124630dcb7c9SBarry Smith     }
124730dcb7c9SBarry Smith   }
124830dcb7c9SBarry Smith 
124930dcb7c9SBarry Smith   /* receive the message lengths */
125030dcb7c9SBarry Smith   nrecvs2 = nsends;
12519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs2+1,&lens2));
12529566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs2+1,&starts3));
12539566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs2+1,&recv_waits));
125430dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
12559566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Irecv(&lens2[i],1,MPIU_INT,dest[i],tag2,comm,recv_waits+i));
125630dcb7c9SBarry Smith   }
1257d44834fbSBarry Smith 
12588a8e0b3aSBarry Smith   /* send the message lengths */
12598a8e0b3aSBarry Smith   for (i=0; i<nsends2; i++) {
12609566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(&nprocs[i],1,MPIU_INT,source[i],tag2,comm));
12618a8e0b3aSBarry Smith   }
12628a8e0b3aSBarry Smith 
1263d44834fbSBarry Smith   /* wait on receives of lens */
12640c468ba9SBarry Smith   if (nrecvs2) {
12659566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrecvs2,&recv_statuses));
12669566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrecvs2,recv_waits,recv_statuses));
12679566063dSJacob Faibussowitsch     PetscCall(PetscFree(recv_statuses));
12680c468ba9SBarry Smith   }
12699566063dSJacob Faibussowitsch   PetscCall(PetscFree(recv_waits));
1270d44834fbSBarry Smith 
127130dcb7c9SBarry Smith   starts3[0] = 0;
1272d44834fbSBarry Smith   nt         = 0;
127330dcb7c9SBarry Smith   for (i=0; i<nrecvs2-1; i++) {
127430dcb7c9SBarry Smith     starts3[i+1] = starts3[i] + lens2[i];
1275d44834fbSBarry Smith     nt          += lens2[i];
127630dcb7c9SBarry Smith   }
127776466f69SStefano Zampini   if (nrecvs2) nt += lens2[nrecvs2-1];
1278d44834fbSBarry Smith 
12799566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt+1,&recvs2));
12809566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs2+1,&recv_waits));
128152b72c4aSBarry Smith   for (i=0; i<nrecvs2; i++) {
12829566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Irecv(recvs2+starts3[i],lens2[i],MPIU_INT,dest[i],tag3,comm,recv_waits+i));
128330dcb7c9SBarry Smith   }
128430dcb7c9SBarry Smith 
128530dcb7c9SBarry Smith   /* send the messages */
12869566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends2+1,&send_waits));
128730dcb7c9SBarry Smith   for (i=0; i<nsends2; i++) {
12889566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Isend(sends2+starts2[i],nprocs[i],MPIU_INT,source[i],tag3,comm,send_waits+i));
128930dcb7c9SBarry Smith   }
129030dcb7c9SBarry Smith 
129130dcb7c9SBarry Smith   /* wait on receives */
12920c468ba9SBarry Smith   if (nrecvs2) {
12939566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrecvs2,&recv_statuses));
12949566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrecvs2,recv_waits,recv_statuses));
12959566063dSJacob Faibussowitsch     PetscCall(PetscFree(recv_statuses));
12960c468ba9SBarry Smith   }
12979566063dSJacob Faibussowitsch   PetscCall(PetscFree(recv_waits));
12989566063dSJacob Faibussowitsch   PetscCall(PetscFree(nprocs));
129930dcb7c9SBarry Smith 
130007b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
130130dcb7c9SBarry Smith     cnt = 0;
130230dcb7c9SBarry Smith     for (i=0; i<nrecvs2; i++) {
130330dcb7c9SBarry Smith       nt = recvs2[cnt++];
130430dcb7c9SBarry Smith       for (j=0; j<nt; j++) {
13059566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(comm,"[%d] local node %" PetscInt_FMT " number of subdomains %" PetscInt_FMT ": ",rank,recvs2[cnt],recvs2[cnt+1]));
130630dcb7c9SBarry Smith         for (k=0; k<recvs2[cnt+1]; k++) {
13079566063dSJacob Faibussowitsch           PetscCall(PetscSynchronizedPrintf(comm,"%" PetscInt_FMT " ",recvs2[cnt+2+k]));
130830dcb7c9SBarry Smith         }
130930dcb7c9SBarry Smith         cnt += 2 + recvs2[cnt+1];
13109566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(comm,"\n"));
131130dcb7c9SBarry Smith       }
131230dcb7c9SBarry Smith     }
13139566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm,PETSC_STDOUT));
131407b52d57SBarry Smith   } /* -----------------------------------  */
131530dcb7c9SBarry Smith 
131630dcb7c9SBarry Smith   /* count number subdomains for each local node */
13179566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(size,&nprocs));
131830dcb7c9SBarry Smith   cnt  = 0;
131930dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
132030dcb7c9SBarry Smith     nt = recvs2[cnt++];
132130dcb7c9SBarry Smith     for (j=0; j<nt; j++) {
1322f6e5521dSKarl Rupp       for (k=0; k<recvs2[cnt+1]; k++) nprocs[recvs2[cnt+2+k]]++;
132330dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt+1];
132430dcb7c9SBarry Smith     }
132530dcb7c9SBarry Smith   }
132630dcb7c9SBarry Smith   nt = 0; for (i=0; i<size; i++) nt += (nprocs[i] > 0);
132730dcb7c9SBarry Smith   *nproc    = nt;
13289566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt+1,procs));
13299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt+1,numprocs));
13309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt+1,indices));
13310298fd71SBarry Smith   for (i=0;i<nt+1;i++) (*indices)[i]=NULL;
13329566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size,&bprocs));
133330dcb7c9SBarry Smith   cnt  = 0;
133430dcb7c9SBarry Smith   for (i=0; i<size; i++) {
133530dcb7c9SBarry Smith     if (nprocs[i] > 0) {
133630dcb7c9SBarry Smith       bprocs[i]        = cnt;
133730dcb7c9SBarry Smith       (*procs)[cnt]    = i;
133830dcb7c9SBarry Smith       (*numprocs)[cnt] = nprocs[i];
13399566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nprocs[i],&(*indices)[cnt]));
134030dcb7c9SBarry Smith       cnt++;
134130dcb7c9SBarry Smith     }
134230dcb7c9SBarry Smith   }
134330dcb7c9SBarry Smith 
134430dcb7c9SBarry Smith   /* make the list of subdomains for each nontrivial local node */
13459566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(*numprocs,nt));
134630dcb7c9SBarry Smith   cnt  = 0;
134730dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
134830dcb7c9SBarry Smith     nt = recvs2[cnt++];
134930dcb7c9SBarry Smith     for (j=0; j<nt; j++) {
1350f6e5521dSKarl Rupp       for (k=0; k<recvs2[cnt+1]; k++) (*indices)[bprocs[recvs2[cnt+2+k]]][(*numprocs)[bprocs[recvs2[cnt+2+k]]]++] = recvs2[cnt];
135130dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt+1];
135230dcb7c9SBarry Smith     }
135330dcb7c9SBarry Smith   }
13549566063dSJacob Faibussowitsch   PetscCall(PetscFree(bprocs));
13559566063dSJacob Faibussowitsch   PetscCall(PetscFree(recvs2));
135630dcb7c9SBarry Smith 
135707b52d57SBarry Smith   /* sort the node indexing by their global numbers */
135807b52d57SBarry Smith   nt = *nproc;
135907b52d57SBarry Smith   for (i=0; i<nt; i++) {
13609566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1((*numprocs)[i],&tmp));
1361f6e5521dSKarl Rupp     for (j=0; j<(*numprocs)[i]; j++) tmp[j] = lindices[(*indices)[i][j]];
13629566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithArray((*numprocs)[i],tmp,(*indices)[i]));
13639566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmp));
136407b52d57SBarry Smith   }
136507b52d57SBarry Smith 
136607b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
136730dcb7c9SBarry Smith     nt = *nproc;
136830dcb7c9SBarry Smith     for (i=0; i<nt; i++) {
13699566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedPrintf(comm,"[%d] subdomain %" PetscInt_FMT " number of indices %" PetscInt_FMT ": ",rank,(*procs)[i],(*numprocs)[i]));
137030dcb7c9SBarry Smith       for (j=0; j<(*numprocs)[i]; j++) {
13719566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(comm,"%" PetscInt_FMT " ",(*indices)[i][j]));
137230dcb7c9SBarry Smith       }
13739566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedPrintf(comm,"\n"));
137430dcb7c9SBarry Smith     }
13759566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm,PETSC_STDOUT));
137607b52d57SBarry Smith   } /* -----------------------------------  */
137730dcb7c9SBarry Smith 
137830dcb7c9SBarry Smith   /* wait on sends */
137930dcb7c9SBarry Smith   if (nsends2) {
13809566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nsends2,&send_status));
13819566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nsends2,send_waits,send_status));
13829566063dSJacob Faibussowitsch     PetscCall(PetscFree(send_status));
138330dcb7c9SBarry Smith   }
138430dcb7c9SBarry Smith 
13859566063dSJacob Faibussowitsch   PetscCall(PetscFree(starts3));
13869566063dSJacob Faibussowitsch   PetscCall(PetscFree(dest));
13879566063dSJacob Faibussowitsch   PetscCall(PetscFree(send_waits));
13883677ff5aSBarry Smith 
13899566063dSJacob Faibussowitsch   PetscCall(PetscFree(nownedsenders));
13909566063dSJacob Faibussowitsch   PetscCall(PetscFree(ownedsenders));
13919566063dSJacob Faibussowitsch   PetscCall(PetscFree(starts));
13929566063dSJacob Faibussowitsch   PetscCall(PetscFree(starts2));
13939566063dSJacob Faibussowitsch   PetscCall(PetscFree(lens2));
139489d82c54SBarry Smith 
13959566063dSJacob Faibussowitsch   PetscCall(PetscFree(source));
13969566063dSJacob Faibussowitsch   PetscCall(PetscFree(len));
13979566063dSJacob Faibussowitsch   PetscCall(PetscFree(recvs));
13989566063dSJacob Faibussowitsch   PetscCall(PetscFree(nprocs));
13999566063dSJacob Faibussowitsch   PetscCall(PetscFree(sends2));
140024cf384cSBarry Smith 
140124cf384cSBarry Smith   /* put the information about myself as the first entry in the list */
140224cf384cSBarry Smith   first_procs    = (*procs)[0];
140324cf384cSBarry Smith   first_numprocs = (*numprocs)[0];
140424cf384cSBarry Smith   first_indices  = (*indices)[0];
140524cf384cSBarry Smith   for (i=0; i<*nproc; i++) {
140624cf384cSBarry Smith     if ((*procs)[i] == rank) {
140724cf384cSBarry Smith       (*procs)[0]    = (*procs)[i];
140824cf384cSBarry Smith       (*numprocs)[0] = (*numprocs)[i];
140924cf384cSBarry Smith       (*indices)[0]  = (*indices)[i];
141024cf384cSBarry Smith       (*procs)[i]    = first_procs;
141124cf384cSBarry Smith       (*numprocs)[i] = first_numprocs;
141224cf384cSBarry Smith       (*indices)[i]  = first_indices;
141324cf384cSBarry Smith       break;
141424cf384cSBarry Smith     }
141524cf384cSBarry Smith   }
1416268a049cSStefano Zampini 
1417268a049cSStefano Zampini   /* save info for reuse */
1418268a049cSStefano Zampini   mapping->info_nproc = *nproc;
1419268a049cSStefano Zampini   mapping->info_procs = *procs;
1420268a049cSStefano Zampini   mapping->info_numprocs = *numprocs;
1421268a049cSStefano Zampini   mapping->info_indices = *indices;
1422268a049cSStefano Zampini   mapping->info_cached = PETSC_TRUE;
142389d82c54SBarry Smith   PetscFunctionReturn(0);
142489d82c54SBarry Smith }
142589d82c54SBarry Smith 
14266a818285SBarry Smith /*@C
14276a818285SBarry Smith     ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by ISLocalToGlobalMappingGetBlockInfo()
14286a818285SBarry Smith 
14296a818285SBarry Smith     Collective on ISLocalToGlobalMapping
14306a818285SBarry Smith 
1431f899ff85SJose E. Roman     Input Parameter:
14326a818285SBarry Smith .   mapping - the mapping from local to global indexing
14336a818285SBarry Smith 
1434d8d19677SJose E. Roman     Output Parameters:
14356a818285SBarry Smith +   nproc - number of processors that are connected to this one
14366a818285SBarry Smith .   proc - neighboring processors
14376a818285SBarry Smith .   numproc - number of indices for each processor
14386a818285SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
14396a818285SBarry Smith 
14406a818285SBarry Smith     Level: advanced
14416a818285SBarry Smith 
1442db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1443db781477SPatrick Sanan           `ISLocalToGlobalMappingGetInfo()`
14446a818285SBarry Smith @*/
14456a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
14466a818285SBarry Smith {
14476a818285SBarry Smith   PetscFunctionBegin;
1448cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1449268a049cSStefano Zampini   if (mapping->info_free) {
14509566063dSJacob Faibussowitsch     PetscCall(PetscFree(*numprocs));
14516a818285SBarry Smith     if (*indices) {
1452268a049cSStefano Zampini       PetscInt i;
1453268a049cSStefano Zampini 
14549566063dSJacob Faibussowitsch       PetscCall(PetscFree((*indices)[0]));
14556a818285SBarry Smith       for (i=1; i<*nproc; i++) {
14569566063dSJacob Faibussowitsch         PetscCall(PetscFree((*indices)[i]));
14576a818285SBarry Smith       }
14589566063dSJacob Faibussowitsch       PetscCall(PetscFree(*indices));
14596a818285SBarry Smith     }
1460268a049cSStefano Zampini   }
1461268a049cSStefano Zampini   *nproc    = 0;
1462268a049cSStefano Zampini   *procs    = NULL;
1463268a049cSStefano Zampini   *numprocs = NULL;
1464268a049cSStefano Zampini   *indices  = NULL;
14656a818285SBarry Smith   PetscFunctionReturn(0);
14666a818285SBarry Smith }
14676a818285SBarry Smith 
14686a818285SBarry Smith /*@C
14696a818285SBarry Smith     ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each processor and
14706a818285SBarry Smith      each index shared by more than one processor
14716a818285SBarry Smith 
14726a818285SBarry Smith     Collective on ISLocalToGlobalMapping
14736a818285SBarry Smith 
1474f899ff85SJose E. Roman     Input Parameter:
14756a818285SBarry Smith .   mapping - the mapping from local to global indexing
14766a818285SBarry Smith 
1477d8d19677SJose E. Roman     Output Parameters:
14786a818285SBarry Smith +   nproc - number of processors that are connected to this one
14796a818285SBarry Smith .   proc - neighboring processors
14806a818285SBarry Smith .   numproc - number of indices for each subdomain (processor)
14816a818285SBarry Smith -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
14826a818285SBarry Smith 
14836a818285SBarry Smith     Level: advanced
14846a818285SBarry Smith 
14851bd0b88eSStefano Zampini     Notes: The user needs to call ISLocalToGlobalMappingRestoreInfo when the data is no longer needed.
14861bd0b88eSStefano Zampini 
14876a818285SBarry Smith     Fortran Usage:
14886a818285SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
14896a818285SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
14906a818285SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
14916a818285SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
14926a818285SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
14936a818285SBarry Smith 
1494db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1495db781477SPatrick Sanan           `ISLocalToGlobalMappingRestoreInfo()`
14966a818285SBarry Smith @*/
14976a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
14986a818285SBarry Smith {
14998a1f772fSStefano Zampini   PetscInt       **bindices = NULL,*bnumprocs = NULL,bs,i,j,k;
15006a818285SBarry Smith 
15016a818285SBarry Smith   PetscFunctionBegin;
15026a818285SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
15038a1f772fSStefano Zampini   bs = mapping->bs;
15049566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockInfo(mapping,nproc,procs,&bnumprocs,&bindices));
1505268a049cSStefano Zampini   if (bs > 1) { /* we need to expand the cached info */
15069566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(*nproc,&*indices));
15079566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(*nproc,&*numprocs));
15086a818285SBarry Smith     for (i=0; i<*nproc; i++) {
15099566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(bs*bnumprocs[i],&(*indices)[i]));
1510268a049cSStefano Zampini       for (j=0; j<bnumprocs[i]; j++) {
15116a818285SBarry Smith         for (k=0; k<bs; k++) {
15126a818285SBarry Smith           (*indices)[i][j*bs+k] = bs*bindices[i][j] + k;
15136a818285SBarry Smith         }
15146a818285SBarry Smith       }
1515268a049cSStefano Zampini       (*numprocs)[i] = bnumprocs[i]*bs;
15166a818285SBarry Smith     }
1517268a049cSStefano Zampini     mapping->info_free = PETSC_TRUE;
1518268a049cSStefano Zampini   } else {
1519268a049cSStefano Zampini     *numprocs = bnumprocs;
1520268a049cSStefano Zampini     *indices  = bindices;
15216a818285SBarry Smith   }
15226a818285SBarry Smith   PetscFunctionReturn(0);
15236a818285SBarry Smith }
15246a818285SBarry Smith 
152507b52d57SBarry Smith /*@C
152607b52d57SBarry Smith     ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by ISLocalToGlobalMappingGetInfo()
152789d82c54SBarry Smith 
152807b52d57SBarry Smith     Collective on ISLocalToGlobalMapping
152907b52d57SBarry Smith 
1530f899ff85SJose E. Roman     Input Parameter:
153107b52d57SBarry Smith .   mapping - the mapping from local to global indexing
153207b52d57SBarry Smith 
1533d8d19677SJose E. Roman     Output Parameters:
153407b52d57SBarry Smith +   nproc - number of processors that are connected to this one
153507b52d57SBarry Smith .   proc - neighboring processors
153607b52d57SBarry Smith .   numproc - number of indices for each processor
153707b52d57SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
153807b52d57SBarry Smith 
153907b52d57SBarry Smith     Level: advanced
154007b52d57SBarry Smith 
1541db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1542db781477SPatrick Sanan           `ISLocalToGlobalMappingGetInfo()`
154307b52d57SBarry Smith @*/
15447087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
154507b52d57SBarry Smith {
154607b52d57SBarry Smith   PetscFunctionBegin;
15479566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockInfo(mapping,nproc,procs,numprocs,indices));
154807b52d57SBarry Smith   PetscFunctionReturn(0);
154907b52d57SBarry Smith }
155086994e45SJed Brown 
155186994e45SJed Brown /*@C
15521bd0b88eSStefano Zampini     ISLocalToGlobalMappingGetNodeInfo - Gets the neighbor information for each node
15531bd0b88eSStefano Zampini 
15541bd0b88eSStefano Zampini     Collective on ISLocalToGlobalMapping
15551bd0b88eSStefano Zampini 
1556f899ff85SJose E. Roman     Input Parameter:
15571bd0b88eSStefano Zampini .   mapping - the mapping from local to global indexing
15581bd0b88eSStefano Zampini 
1559d8d19677SJose E. Roman     Output Parameters:
15601bd0b88eSStefano Zampini +   nnodes - number of local nodes (same ISLocalToGlobalMappingGetSize())
15611bd0b88eSStefano Zampini .   count - number of neighboring processors per node
15621bd0b88eSStefano Zampini -   indices - indices of processes sharing the node (sorted)
15631bd0b88eSStefano Zampini 
15641bd0b88eSStefano Zampini     Level: advanced
15651bd0b88eSStefano Zampini 
15661bd0b88eSStefano Zampini     Notes: The user needs to call ISLocalToGlobalMappingRestoreInfo when the data is no longer needed.
15671bd0b88eSStefano Zampini 
1568db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1569db781477SPatrick Sanan           `ISLocalToGlobalMappingGetInfo()`, `ISLocalToGlobalMappingRestoreNodeInfo()`
15701bd0b88eSStefano Zampini @*/
15711bd0b88eSStefano Zampini PetscErrorCode  ISLocalToGlobalMappingGetNodeInfo(ISLocalToGlobalMapping mapping,PetscInt *nnodes,PetscInt *count[],PetscInt **indices[])
15721bd0b88eSStefano Zampini {
15731bd0b88eSStefano Zampini   PetscInt       n;
15741bd0b88eSStefano Zampini 
15751bd0b88eSStefano Zampini   PetscFunctionBegin;
15761bd0b88eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
15779566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetSize(mapping,&n));
15781bd0b88eSStefano Zampini   if (!mapping->info_nodec) {
15791bd0b88eSStefano Zampini     PetscInt i,m,n_neigh,*neigh,*n_shared,**shared;
15801bd0b88eSStefano Zampini 
15819566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(n+1,&mapping->info_nodec,n,&mapping->info_nodei));
15829566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetInfo(mapping,&n_neigh,&neigh,&n_shared,&shared));
1583071fcb05SBarry Smith     for (i=0;i<n;i++) { mapping->info_nodec[i] = 1;}
1584071fcb05SBarry Smith     m = n;
1585071fcb05SBarry Smith     mapping->info_nodec[n] = 0;
15861bd0b88eSStefano Zampini     for (i=1;i<n_neigh;i++) {
15871bd0b88eSStefano Zampini       PetscInt j;
15881bd0b88eSStefano Zampini 
15891bd0b88eSStefano Zampini       m += n_shared[i];
15901bd0b88eSStefano Zampini       for (j=0;j<n_shared[i];j++) mapping->info_nodec[shared[i][j]] += 1;
15911bd0b88eSStefano Zampini     }
15929566063dSJacob Faibussowitsch     if (n) PetscCall(PetscMalloc1(m,&mapping->info_nodei[0]));
15931bd0b88eSStefano Zampini     for (i=1;i<n;i++) mapping->info_nodei[i] = mapping->info_nodei[i-1] + mapping->info_nodec[i-1];
15949566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(mapping->info_nodec,n));
15951bd0b88eSStefano Zampini     for (i=0;i<n;i++) { mapping->info_nodec[i] = 1; mapping->info_nodei[i][0] = neigh[0]; }
15961bd0b88eSStefano Zampini     for (i=1;i<n_neigh;i++) {
15971bd0b88eSStefano Zampini       PetscInt j;
15981bd0b88eSStefano Zampini 
15991bd0b88eSStefano Zampini       for (j=0;j<n_shared[i];j++) {
16001bd0b88eSStefano Zampini         PetscInt k = shared[i][j];
16011bd0b88eSStefano Zampini 
16021bd0b88eSStefano Zampini         mapping->info_nodei[k][mapping->info_nodec[k]] = neigh[i];
16031bd0b88eSStefano Zampini         mapping->info_nodec[k] += 1;
16041bd0b88eSStefano Zampini       }
16051bd0b88eSStefano Zampini     }
16069566063dSJacob Faibussowitsch     for (i=0;i<n;i++) PetscCall(PetscSortRemoveDupsInt(&mapping->info_nodec[i],mapping->info_nodei[i]));
16079566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreInfo(mapping,&n_neigh,&neigh,&n_shared,&shared));
16081bd0b88eSStefano Zampini   }
16091bd0b88eSStefano Zampini   if (nnodes)  *nnodes  = n;
16101bd0b88eSStefano Zampini   if (count)   *count   = mapping->info_nodec;
16111bd0b88eSStefano Zampini   if (indices) *indices = mapping->info_nodei;
16121bd0b88eSStefano Zampini   PetscFunctionReturn(0);
16131bd0b88eSStefano Zampini }
16141bd0b88eSStefano Zampini 
16151bd0b88eSStefano Zampini /*@C
16161bd0b88eSStefano Zampini     ISLocalToGlobalMappingRestoreNodeInfo - Frees the memory allocated by ISLocalToGlobalMappingGetNodeInfo()
16171bd0b88eSStefano Zampini 
16181bd0b88eSStefano Zampini     Collective on ISLocalToGlobalMapping
16191bd0b88eSStefano Zampini 
1620f899ff85SJose E. Roman     Input Parameter:
16211bd0b88eSStefano Zampini .   mapping - the mapping from local to global indexing
16221bd0b88eSStefano Zampini 
1623d8d19677SJose E. Roman     Output Parameters:
16241bd0b88eSStefano Zampini +   nnodes - number of local nodes
16251bd0b88eSStefano Zampini .   count - number of neighboring processors per node
16261bd0b88eSStefano Zampini -   indices - indices of processes sharing the node (sorted)
16271bd0b88eSStefano Zampini 
16281bd0b88eSStefano Zampini     Level: advanced
16291bd0b88eSStefano Zampini 
1630db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`,
1631db781477SPatrick Sanan           `ISLocalToGlobalMappingGetInfo()`
16321bd0b88eSStefano Zampini @*/
16331bd0b88eSStefano Zampini PetscErrorCode  ISLocalToGlobalMappingRestoreNodeInfo(ISLocalToGlobalMapping mapping,PetscInt *nnodes,PetscInt *count[],PetscInt **indices[])
16341bd0b88eSStefano Zampini {
16351bd0b88eSStefano Zampini   PetscFunctionBegin;
16361bd0b88eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
16371bd0b88eSStefano Zampini   if (nnodes)  *nnodes  = 0;
16381bd0b88eSStefano Zampini   if (count)   *count   = NULL;
16391bd0b88eSStefano Zampini   if (indices) *indices = NULL;
16401bd0b88eSStefano Zampini   PetscFunctionReturn(0);
16411bd0b88eSStefano Zampini }
16421bd0b88eSStefano Zampini 
16431bd0b88eSStefano Zampini /*@C
1644107e9a97SBarry Smith    ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped
164586994e45SJed Brown 
164686994e45SJed Brown    Not Collective
164786994e45SJed Brown 
16484165533cSJose E. Roman    Input Parameter:
164986994e45SJed Brown . ltog - local to global mapping
165086994e45SJed Brown 
16514165533cSJose E. Roman    Output Parameter:
1652565245c5SBarry Smith . array - array of indices, the length of this array may be obtained with ISLocalToGlobalMappingGetSize()
165386994e45SJed Brown 
165486994e45SJed Brown    Level: advanced
165586994e45SJed Brown 
165695452b02SPatrick Sanan    Notes:
165795452b02SPatrick Sanan     ISLocalToGlobalMappingGetSize() returns the length the this array
1658107e9a97SBarry Smith 
1659db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingRestoreIndices()`, `ISLocalToGlobalMappingGetBlockIndices()`, `ISLocalToGlobalMappingRestoreBlockIndices()`
166086994e45SJed Brown @*/
16617087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
166286994e45SJed Brown {
166386994e45SJed Brown   PetscFunctionBegin;
166486994e45SJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
166586994e45SJed Brown   PetscValidPointer(array,2);
166645b6f7e9SBarry Smith   if (ltog->bs == 1) {
166786994e45SJed Brown     *array = ltog->indices;
166845b6f7e9SBarry Smith   } else {
166945b6f7e9SBarry Smith     PetscInt       *jj,k,i,j,n = ltog->n, bs = ltog->bs;
167045b6f7e9SBarry Smith     const PetscInt *ii;
167145b6f7e9SBarry Smith 
16729566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(bs*n,&jj));
167345b6f7e9SBarry Smith     *array = jj;
167445b6f7e9SBarry Smith     k    = 0;
167545b6f7e9SBarry Smith     ii   = ltog->indices;
167645b6f7e9SBarry Smith     for (i=0; i<n; i++)
167745b6f7e9SBarry Smith       for (j=0; j<bs; j++)
167845b6f7e9SBarry Smith         jj[k++] = bs*ii[i] + j;
167945b6f7e9SBarry Smith   }
168086994e45SJed Brown   PetscFunctionReturn(0);
168186994e45SJed Brown }
168286994e45SJed Brown 
168386994e45SJed Brown /*@C
1684193a2b41SJulian Andrej    ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with ISLocalToGlobalMappingGetIndices()
168586994e45SJed Brown 
168686994e45SJed Brown    Not Collective
168786994e45SJed Brown 
16884165533cSJose E. Roman    Input Parameters:
168986994e45SJed Brown + ltog - local to global mapping
169086994e45SJed Brown - array - array of indices
169186994e45SJed Brown 
169286994e45SJed Brown    Level: advanced
169386994e45SJed Brown 
1694db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingGetIndices()`
169586994e45SJed Brown @*/
16967087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
169786994e45SJed Brown {
169886994e45SJed Brown   PetscFunctionBegin;
169986994e45SJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
170086994e45SJed Brown   PetscValidPointer(array,2);
1701c9cc58a2SBarry Smith   PetscCheck(ltog->bs != 1 || *array == ltog->indices,PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
170245b6f7e9SBarry Smith 
170345b6f7e9SBarry Smith   if (ltog->bs > 1) {
17049566063dSJacob Faibussowitsch     PetscCall(PetscFree(*(void**)array));
170545b6f7e9SBarry Smith   }
170645b6f7e9SBarry Smith   PetscFunctionReturn(0);
170745b6f7e9SBarry Smith }
170845b6f7e9SBarry Smith 
170945b6f7e9SBarry Smith /*@C
171045b6f7e9SBarry Smith    ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block
171145b6f7e9SBarry Smith 
171245b6f7e9SBarry Smith    Not Collective
171345b6f7e9SBarry Smith 
17144165533cSJose E. Roman    Input Parameter:
171545b6f7e9SBarry Smith . ltog - local to global mapping
171645b6f7e9SBarry Smith 
17174165533cSJose E. Roman    Output Parameter:
171845b6f7e9SBarry Smith . array - array of indices
171945b6f7e9SBarry Smith 
172045b6f7e9SBarry Smith    Level: advanced
172145b6f7e9SBarry Smith 
1722db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingRestoreBlockIndices()`
172345b6f7e9SBarry Smith @*/
172445b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
172545b6f7e9SBarry Smith {
172645b6f7e9SBarry Smith   PetscFunctionBegin;
172745b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
172845b6f7e9SBarry Smith   PetscValidPointer(array,2);
172945b6f7e9SBarry Smith   *array = ltog->indices;
173045b6f7e9SBarry Smith   PetscFunctionReturn(0);
173145b6f7e9SBarry Smith }
173245b6f7e9SBarry Smith 
173345b6f7e9SBarry Smith /*@C
173445b6f7e9SBarry Smith    ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with ISLocalToGlobalMappingGetBlockIndices()
173545b6f7e9SBarry Smith 
173645b6f7e9SBarry Smith    Not Collective
173745b6f7e9SBarry Smith 
17384165533cSJose E. Roman    Input Parameters:
173945b6f7e9SBarry Smith + ltog - local to global mapping
174045b6f7e9SBarry Smith - array - array of indices
174145b6f7e9SBarry Smith 
174245b6f7e9SBarry Smith    Level: advanced
174345b6f7e9SBarry Smith 
1744db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingGetIndices()`
174545b6f7e9SBarry Smith @*/
174645b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
174745b6f7e9SBarry Smith {
174845b6f7e9SBarry Smith   PetscFunctionBegin;
174945b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
175045b6f7e9SBarry Smith   PetscValidPointer(array,2);
175108401ef6SPierre Jolivet   PetscCheck(*array == ltog->indices,PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
17520298fd71SBarry Smith   *array = NULL;
175386994e45SJed Brown   PetscFunctionReturn(0);
175486994e45SJed Brown }
1755f7efa3c7SJed Brown 
1756f7efa3c7SJed Brown /*@C
1757f7efa3c7SJed Brown    ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings
1758f7efa3c7SJed Brown 
1759f7efa3c7SJed Brown    Not Collective
1760f7efa3c7SJed Brown 
17614165533cSJose E. Roman    Input Parameters:
1762f7efa3c7SJed Brown + comm - communicator for the new mapping, must contain the communicator of every mapping to concatenate
1763f7efa3c7SJed Brown . n - number of mappings to concatenate
1764f7efa3c7SJed Brown - ltogs - local to global mappings
1765f7efa3c7SJed Brown 
17664165533cSJose E. Roman    Output Parameter:
1767f7efa3c7SJed Brown . ltogcat - new mapping
1768f7efa3c7SJed Brown 
17699d90f715SBarry Smith    Note: this currently always returns a mapping with block size of 1
17709d90f715SBarry Smith 
17719d90f715SBarry Smith    Developer Note: If all the input mapping have the same block size we could easily handle that as a special case
17729d90f715SBarry Smith 
1773f7efa3c7SJed Brown    Level: advanced
1774f7efa3c7SJed Brown 
1775db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`
1776f7efa3c7SJed Brown @*/
1777f7efa3c7SJed Brown PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm,PetscInt n,const ISLocalToGlobalMapping ltogs[],ISLocalToGlobalMapping *ltogcat)
1778f7efa3c7SJed Brown {
1779f7efa3c7SJed Brown   PetscInt       i,cnt,m,*idx;
1780f7efa3c7SJed Brown 
1781f7efa3c7SJed Brown   PetscFunctionBegin;
178208401ef6SPierre Jolivet   PetscCheck(n >= 0,comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have a non-negative number of mappings, given %" PetscInt_FMT,n);
1783f7efa3c7SJed Brown   if (n > 0) PetscValidPointer(ltogs,3);
1784f7efa3c7SJed Brown   for (i=0; i<n; i++) PetscValidHeaderSpecific(ltogs[i],IS_LTOGM_CLASSID,3);
1785f7efa3c7SJed Brown   PetscValidPointer(ltogcat,4);
1786f7efa3c7SJed Brown   for (cnt=0,i=0; i<n; i++) {
17879566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(ltogs[i],&m));
1788f7efa3c7SJed Brown     cnt += m;
1789f7efa3c7SJed Brown   }
17909566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(cnt,&idx));
1791f7efa3c7SJed Brown   for (cnt=0,i=0; i<n; i++) {
1792f7efa3c7SJed Brown     const PetscInt *subidx;
17939566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(ltogs[i],&m));
17949566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetIndices(ltogs[i],&subidx));
17959566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(&idx[cnt],subidx,m));
17969566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreIndices(ltogs[i],&subidx));
1797f7efa3c7SJed Brown     cnt += m;
1798f7efa3c7SJed Brown   }
17999566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(comm,1,cnt,idx,PETSC_OWN_POINTER,ltogcat));
1800f7efa3c7SJed Brown   PetscFunctionReturn(0);
1801f7efa3c7SJed Brown }
180204a59952SBarry Smith 
1803413f72f0SBarry Smith /*MC
1804413f72f0SBarry Smith       ISLOCALTOGLOBALMAPPINGBASIC - basic implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is
1805413f72f0SBarry Smith                                     used this is good for only small and moderate size problems.
1806413f72f0SBarry Smith 
1807413f72f0SBarry Smith    Options Database Keys:
1808a2b725a8SWilliam Gropp .   -islocaltoglobalmapping_type basic - select this method
1809413f72f0SBarry Smith 
1810413f72f0SBarry Smith    Level: beginner
1811413f72f0SBarry Smith 
1812db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`, `ISLOCALTOGLOBALMAPPINGHASH`
1813413f72f0SBarry Smith M*/
1814413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Basic(ISLocalToGlobalMapping ltog)
1815413f72f0SBarry Smith {
1816413f72f0SBarry Smith   PetscFunctionBegin;
1817413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Basic;
1818413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Basic;
1819413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Basic;
1820413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Basic;
1821413f72f0SBarry Smith   PetscFunctionReturn(0);
1822413f72f0SBarry Smith }
1823413f72f0SBarry Smith 
1824413f72f0SBarry Smith /*MC
1825413f72f0SBarry Smith       ISLOCALTOGLOBALMAPPINGHASH - hash implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is
1826ed56e8eaSBarry Smith                                     used this is good for large memory problems.
1827413f72f0SBarry Smith 
1828413f72f0SBarry Smith    Options Database Keys:
1829a2b725a8SWilliam Gropp .   -islocaltoglobalmapping_type hash - select this method
1830413f72f0SBarry Smith 
183195452b02SPatrick Sanan    Notes:
183295452b02SPatrick Sanan     This is selected automatically for large problems if the user does not set the type.
1833ed56e8eaSBarry Smith 
1834413f72f0SBarry Smith    Level: beginner
1835413f72f0SBarry Smith 
1836db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`, `ISLOCALTOGLOBALMAPPINGHASH`
1837413f72f0SBarry Smith M*/
1838413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Hash(ISLocalToGlobalMapping ltog)
1839413f72f0SBarry Smith {
1840413f72f0SBarry Smith   PetscFunctionBegin;
1841413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Hash;
1842413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Hash;
1843413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Hash;
1844413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Hash;
1845413f72f0SBarry Smith   PetscFunctionReturn(0);
1846413f72f0SBarry Smith }
1847413f72f0SBarry Smith 
1848413f72f0SBarry Smith /*@C
1849413f72f0SBarry Smith     ISLocalToGlobalMappingRegister -  Adds a method for applying a global to local mapping with an ISLocalToGlobalMapping
1850413f72f0SBarry Smith 
1851413f72f0SBarry Smith    Not Collective
1852413f72f0SBarry Smith 
1853413f72f0SBarry Smith    Input Parameters:
1854413f72f0SBarry Smith +  sname - name of a new method
1855413f72f0SBarry Smith -  routine_create - routine to create method context
1856413f72f0SBarry Smith 
1857413f72f0SBarry Smith    Notes:
1858ed56e8eaSBarry Smith    ISLocalToGlobalMappingRegister() may be called multiple times to add several user-defined mappings.
1859413f72f0SBarry Smith 
1860413f72f0SBarry Smith    Sample usage:
1861413f72f0SBarry Smith .vb
1862413f72f0SBarry Smith    ISLocalToGlobalMappingRegister("my_mapper",MyCreate);
1863413f72f0SBarry Smith .ve
1864413f72f0SBarry Smith 
1865ed56e8eaSBarry Smith    Then, your mapping can be chosen with the procedural interface via
1866413f72f0SBarry Smith $     ISLocalToGlobalMappingSetType(ltog,"my_mapper")
1867413f72f0SBarry Smith    or at runtime via the option
1868ed56e8eaSBarry Smith $     -islocaltoglobalmapping_type my_mapper
1869413f72f0SBarry Smith 
1870413f72f0SBarry Smith    Level: advanced
1871413f72f0SBarry Smith 
1872db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingRegisterAll()`, `ISLocalToGlobalMappingRegisterDestroy()`, `ISLOCALTOGLOBALMAPPINGBASIC`, `ISLOCALTOGLOBALMAPPINGHASH`
1873413f72f0SBarry Smith 
1874413f72f0SBarry Smith @*/
1875413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRegister(const char sname[],PetscErrorCode (*function)(ISLocalToGlobalMapping))
1876413f72f0SBarry Smith {
1877413f72f0SBarry Smith   PetscFunctionBegin;
18789566063dSJacob Faibussowitsch   PetscCall(ISInitializePackage());
18799566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&ISLocalToGlobalMappingList,sname,function));
1880413f72f0SBarry Smith   PetscFunctionReturn(0);
1881413f72f0SBarry Smith }
1882413f72f0SBarry Smith 
1883413f72f0SBarry Smith /*@C
1884ed56e8eaSBarry Smith    ISLocalToGlobalMappingSetType - Builds ISLocalToGlobalMapping for a particular global to local mapping approach.
1885413f72f0SBarry Smith 
1886413f72f0SBarry Smith    Logically Collective on ISLocalToGlobalMapping
1887413f72f0SBarry Smith 
1888413f72f0SBarry Smith    Input Parameters:
1889413f72f0SBarry Smith +  ltog - the ISLocalToGlobalMapping object
1890413f72f0SBarry Smith -  type - a known method
1891413f72f0SBarry Smith 
1892413f72f0SBarry Smith    Options Database Key:
1893ed56e8eaSBarry Smith .  -islocaltoglobalmapping_type  <method> - Sets the method; use -help for a list
1894413f72f0SBarry Smith     of available methods (for instance, basic or hash)
1895413f72f0SBarry Smith 
1896413f72f0SBarry Smith    Notes:
1897413f72f0SBarry Smith    See "petsc/include/petscis.h" for available methods
1898413f72f0SBarry Smith 
1899413f72f0SBarry Smith   Normally, it is best to use the ISLocalToGlobalMappingSetFromOptions() command and
1900413f72f0SBarry Smith   then set the ISLocalToGlobalMapping type from the options database rather than by using
1901413f72f0SBarry Smith   this routine.
1902413f72f0SBarry Smith 
1903413f72f0SBarry Smith   Level: intermediate
1904413f72f0SBarry Smith 
1905413f72f0SBarry Smith   Developer Note: ISLocalToGlobalMappingRegister() is used to add new types to ISLocalToGlobalMappingList from which they
1906413f72f0SBarry Smith   are accessed by ISLocalToGlobalMappingSetType().
1907413f72f0SBarry Smith 
1908db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingRegister()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingGetType()`
1909413f72f0SBarry Smith @*/
1910413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType type)
1911413f72f0SBarry Smith {
1912413f72f0SBarry Smith   PetscBool      match;
19135f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(ISLocalToGlobalMapping) = NULL;
1914413f72f0SBarry Smith 
1915413f72f0SBarry Smith   PetscFunctionBegin;
1916413f72f0SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
1917a0d79125SStefano Zampini   if (type) PetscValidCharPointer(type,2);
1918413f72f0SBarry Smith 
19199566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ltog,type,&match));
1920413f72f0SBarry Smith   if (match) PetscFunctionReturn(0);
1921413f72f0SBarry Smith 
1922a0d79125SStefano Zampini   /* L2G maps defer type setup at globaltolocal calls, allow passing NULL here */
1923a0d79125SStefano Zampini   if (type) {
19249566063dSJacob Faibussowitsch     PetscCall(PetscFunctionListFind(ISLocalToGlobalMappingList,type,&r));
1925a0d79125SStefano Zampini     PetscCheck(r,PetscObjectComm((PetscObject)ltog),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested ISLocalToGlobalMapping type %s",type);
1926a0d79125SStefano Zampini   }
1927413f72f0SBarry Smith   /* Destroy the previous private LTOG context */
1928*dbbe0bcdSBarry Smith   PetscTryTypeMethod(ltog,destroy);
1929413f72f0SBarry Smith   ltog->ops->destroy = NULL;
1930*dbbe0bcdSBarry Smith 
19319566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)ltog,type));
19329566063dSJacob Faibussowitsch   if (r) PetscCall((*r)(ltog));
1933a0d79125SStefano Zampini   PetscFunctionReturn(0);
1934a0d79125SStefano Zampini }
1935a0d79125SStefano Zampini 
1936a0d79125SStefano Zampini /*@C
1937a0d79125SStefano Zampini    ISLocalToGlobalMappingGetType - Get the type of the l2g map
1938a0d79125SStefano Zampini 
1939a0d79125SStefano Zampini    Not Collective
1940a0d79125SStefano Zampini 
1941a0d79125SStefano Zampini    Input Parameter:
1942a0d79125SStefano Zampini .  ltog - the ISLocalToGlobalMapping object
1943a0d79125SStefano Zampini 
1944a0d79125SStefano Zampini    Output Parameter:
1945a0d79125SStefano Zampini .  type - the type
1946a0d79125SStefano Zampini 
194749762cbcSSatish Balay    Level: intermediate
194849762cbcSSatish Balay 
1949db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingRegister()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`
1950a0d79125SStefano Zampini @*/
1951a0d79125SStefano Zampini PetscErrorCode  ISLocalToGlobalMappingGetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType *type)
1952a0d79125SStefano Zampini {
1953a0d79125SStefano Zampini   PetscFunctionBegin;
1954a0d79125SStefano Zampini   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
1955a0d79125SStefano Zampini   PetscValidPointer(type,2);
1956a0d79125SStefano Zampini   *type = ((PetscObject)ltog)->type_name;
1957413f72f0SBarry Smith   PetscFunctionReturn(0);
1958413f72f0SBarry Smith }
1959413f72f0SBarry Smith 
1960413f72f0SBarry Smith PetscBool ISLocalToGlobalMappingRegisterAllCalled = PETSC_FALSE;
1961413f72f0SBarry Smith 
1962413f72f0SBarry Smith /*@C
1963413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAll - Registers all of the local to global mapping components in the IS package.
1964413f72f0SBarry Smith 
1965413f72f0SBarry Smith   Not Collective
1966413f72f0SBarry Smith 
1967413f72f0SBarry Smith   Level: advanced
1968413f72f0SBarry Smith 
1969db781477SPatrick Sanan .seealso: `ISRegister()`, `ISLocalToGlobalRegister()`
1970413f72f0SBarry Smith @*/
1971413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRegisterAll(void)
1972413f72f0SBarry Smith {
1973413f72f0SBarry Smith   PetscFunctionBegin;
1974413f72f0SBarry Smith   if (ISLocalToGlobalMappingRegisterAllCalled) PetscFunctionReturn(0);
1975413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAllCalled = PETSC_TRUE;
19769566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGBASIC, ISLocalToGlobalMappingCreate_Basic));
19779566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGHASH, ISLocalToGlobalMappingCreate_Hash));
1978413f72f0SBarry Smith   PetscFunctionReturn(0);
1979413f72f0SBarry Smith }
1980