xref: /petsc/src/vec/is/utils/isltog.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
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 
416528b96dSMatthew G. Knepley .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 
806528b96dSMatthew G. Knepley .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 
1136528b96dSMatthew G. Knepley .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   }
1599566063dSJacob Faibussowitsch   if (mapping->ops->globaltolocalmappingsetup) PetscCall((*mapping->ops->globaltolocalmappingsetup)(mapping));
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 
2716658fb44Sstefano_zampini .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 
2983b9aefa3SBarry Smith .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
320fe2efc57SMark .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 
3415a5d4f66SBarry Smith .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 
3867e99dc12SLawrence Mitchell .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 
4337e99dc12SLawrence Mitchell .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 
47363fa5c83Sstefano_zampini .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);
483*08401ef6SPierre 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;
489*08401ef6SPierre 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; }
497*08401ef6SPierre 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) {
501*08401ef6SPierre 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 
530413f72f0SBarry Smith   if (mapping->ops->destroy) {
5319566063dSJacob Faibussowitsch     PetscCall((*mapping->ops->destroy)(mapping));
532413f72f0SBarry Smith   }
53363fa5c83Sstefano_zampini   PetscFunctionReturn(0);
53463fa5c83Sstefano_zampini }
53563fa5c83Sstefano_zampini 
53645b6f7e9SBarry Smith /*@
53745b6f7e9SBarry Smith     ISLocalToGlobalMappingGetBlockSize - Gets the blocksize of the mapping
53845b6f7e9SBarry Smith     ordering and a global parallel ordering.
53945b6f7e9SBarry Smith 
54045b6f7e9SBarry Smith     Not Collective
54145b6f7e9SBarry Smith 
54245b6f7e9SBarry Smith     Input Parameters:
54345b6f7e9SBarry Smith .   mapping - mapping data structure
54445b6f7e9SBarry Smith 
54545b6f7e9SBarry Smith     Output Parameter:
54645b6f7e9SBarry Smith .   bs - the blocksize
54745b6f7e9SBarry Smith 
54845b6f7e9SBarry Smith     Level: advanced
54945b6f7e9SBarry Smith 
55045b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS()
55145b6f7e9SBarry Smith @*/
55245b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping mapping,PetscInt *bs)
55345b6f7e9SBarry Smith {
55445b6f7e9SBarry Smith   PetscFunctionBegin;
555cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
55645b6f7e9SBarry Smith   *bs = mapping->bs;
55745b6f7e9SBarry Smith   PetscFunctionReturn(0);
55845b6f7e9SBarry Smith }
55945b6f7e9SBarry Smith 
560ba5bb76aSSatish Balay /*@
56190f02eecSBarry Smith     ISLocalToGlobalMappingCreate - Creates a mapping between a local (0 to n)
56290f02eecSBarry Smith     ordering and a global parallel ordering.
5632362add9SBarry Smith 
56489d82c54SBarry Smith     Not Collective, but communicator may have more than one process
565b9cd556bSLois Curfman McInnes 
5662362add9SBarry Smith     Input Parameters:
56789d82c54SBarry Smith +   comm - MPI communicator
568f0413b6fSBarry Smith .   bs - the block size
56928bc9809SBarry Smith .   n - the number of local elements divided by the block size, or equivalently the number of block indices
57028bc9809SBarry Smith .   indices - the global index for each local element, these do not need to be in increasing order (sorted), these values should not be scaled (i.e. multiplied) by the blocksize bs
571d5ad8652SBarry Smith -   mode - see PetscCopyMode
5722362add9SBarry Smith 
573a997ad1aSLois Curfman McInnes     Output Parameter:
57490f02eecSBarry Smith .   mapping - new mapping data structure
5752362add9SBarry Smith 
57695452b02SPatrick Sanan     Notes:
57795452b02SPatrick Sanan     There is one integer value in indices per block and it represents the actual indices bs*idx + j, where j=0,..,bs-1
578413f72f0SBarry Smith 
5799a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
580413f72f0SBarry 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.
581413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
582413f72f0SBarry Smith 
583a997ad1aSLois Curfman McInnes     Level: advanced
584a997ad1aSLois Curfman McInnes 
585413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingSetFromOptions(), ISLOCALTOGLOBALMAPPINGBASIC, ISLOCALTOGLOBALMAPPINGHASH
586413f72f0SBarry Smith           ISLocalToGlobalMappingSetType(), ISLocalToGlobalMappingType
5872362add9SBarry Smith @*/
58860c7cefcSBarry Smith PetscErrorCode  ISLocalToGlobalMappingCreate(MPI_Comm comm,PetscInt bs,PetscInt n,const PetscInt indices[],PetscCopyMode mode,ISLocalToGlobalMapping *mapping)
5892362add9SBarry Smith {
59032dcc486SBarry Smith   PetscInt       *in;
591b46b645bSBarry Smith 
592b46b645bSBarry Smith   PetscFunctionBegin;
593064a246eSJacob Faibussowitsch   if (n) PetscValidIntPointer(indices,4);
594064a246eSJacob Faibussowitsch   PetscValidPointer(mapping,6);
595b46b645bSBarry Smith 
5960298fd71SBarry Smith   *mapping = NULL;
5979566063dSJacob Faibussowitsch   PetscCall(ISInitializePackage());
5982362add9SBarry Smith 
5999566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(*mapping,IS_LTOGM_CLASSID,"ISLocalToGlobalMapping","Local to global mapping","IS",comm,ISLocalToGlobalMappingDestroy,ISLocalToGlobalMappingView));
600d4bb536fSBarry Smith   (*mapping)->n  = n;
601f0413b6fSBarry Smith   (*mapping)->bs = bs;
602d5ad8652SBarry Smith   if (mode == PETSC_COPY_VALUES) {
6039566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n,&in));
6049566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(in,indices,n));
605d5ad8652SBarry Smith     (*mapping)->indices = in;
60671910c26SVaclav Hapla     (*mapping)->dealloc_indices = PETSC_TRUE;
6079566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt)));
6086389a1a1SBarry Smith   } else if (mode == PETSC_OWN_POINTER) {
6096389a1a1SBarry Smith     (*mapping)->indices = (PetscInt*)indices;
61071910c26SVaclav Hapla     (*mapping)->dealloc_indices = PETSC_TRUE;
6119566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt)));
61271910c26SVaclav Hapla   } else if (mode == PETSC_USE_POINTER) {
61371910c26SVaclav Hapla     (*mapping)->indices = (PetscInt*)indices;
6146389a1a1SBarry Smith   }
61598921bdaSJacob Faibussowitsch   else SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid mode %d", mode);
6163a40ed3dSBarry Smith   PetscFunctionReturn(0);
6172362add9SBarry Smith }
6182362add9SBarry Smith 
619413f72f0SBarry Smith PetscFunctionList ISLocalToGlobalMappingList = NULL;
620413f72f0SBarry Smith 
62190f02eecSBarry Smith /*@
6227e99dc12SLawrence Mitchell    ISLocalToGlobalMappingSetFromOptions - Set mapping options from the options database.
6237e99dc12SLawrence Mitchell 
6247e99dc12SLawrence Mitchell    Not collective
6257e99dc12SLawrence Mitchell 
6267e99dc12SLawrence Mitchell    Input Parameters:
6277e99dc12SLawrence Mitchell .  mapping - mapping data structure
6287e99dc12SLawrence Mitchell 
6297e99dc12SLawrence Mitchell    Level: advanced
6307e99dc12SLawrence Mitchell 
6317e99dc12SLawrence Mitchell @*/
6327e99dc12SLawrence Mitchell PetscErrorCode ISLocalToGlobalMappingSetFromOptions(ISLocalToGlobalMapping mapping)
6337e99dc12SLawrence Mitchell {
6347e99dc12SLawrence Mitchell   PetscErrorCode             ierr;
635413f72f0SBarry Smith   char                       type[256];
636413f72f0SBarry Smith   ISLocalToGlobalMappingType defaulttype = "Not set";
6377e99dc12SLawrence Mitchell   PetscBool                  flg;
6387e99dc12SLawrence Mitchell 
6397e99dc12SLawrence Mitchell   PetscFunctionBegin;
6407e99dc12SLawrence Mitchell   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
6419566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRegisterAll());
6429566063dSJacob Faibussowitsch   ierr = PetscObjectOptionsBegin((PetscObject)mapping);PetscCall(ierr);
6439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-islocaltoglobalmapping_type","ISLocalToGlobalMapping method","ISLocalToGlobalMappingSetType",ISLocalToGlobalMappingList,(char*)(((PetscObject)mapping)->type_name) ? ((PetscObject)mapping)->type_name : defaulttype,type,256,&flg));
644413f72f0SBarry Smith   if (flg) {
6459566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingSetType(mapping,type));
646413f72f0SBarry Smith   }
6479566063dSJacob Faibussowitsch   ierr = PetscOptionsEnd();PetscCall(ierr);
6487e99dc12SLawrence Mitchell   PetscFunctionReturn(0);
6497e99dc12SLawrence Mitchell }
6507e99dc12SLawrence Mitchell 
6517e99dc12SLawrence Mitchell /*@
65290f02eecSBarry Smith    ISLocalToGlobalMappingDestroy - Destroys a mapping between a local (0 to n)
65390f02eecSBarry Smith    ordering and a global parallel ordering.
65490f02eecSBarry Smith 
6550f5bd95cSBarry Smith    Note Collective
656b9cd556bSLois Curfman McInnes 
65790f02eecSBarry Smith    Input Parameters:
65890f02eecSBarry Smith .  mapping - mapping data structure
65990f02eecSBarry Smith 
660a997ad1aSLois Curfman McInnes    Level: advanced
661a997ad1aSLois Curfman McInnes 
6623acfe500SLois Curfman McInnes .seealso: ISLocalToGlobalMappingCreate()
66390f02eecSBarry Smith @*/
6646bf464f9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping *mapping)
66590f02eecSBarry Smith {
6663a40ed3dSBarry Smith   PetscFunctionBegin;
6676bf464f9SBarry Smith   if (!*mapping) PetscFunctionReturn(0);
6686bf464f9SBarry Smith   PetscValidHeaderSpecific((*mapping),IS_LTOGM_CLASSID,1);
6694c8fdceaSLisandro Dalcin   if (--((PetscObject)(*mapping))->refct > 0) {*mapping = NULL;PetscFunctionReturn(0);}
67071910c26SVaclav Hapla   if ((*mapping)->dealloc_indices) {
6719566063dSJacob Faibussowitsch     PetscCall(PetscFree((*mapping)->indices));
67271910c26SVaclav Hapla   }
6739566063dSJacob Faibussowitsch   PetscCall(PetscFree((*mapping)->info_procs));
6749566063dSJacob Faibussowitsch   PetscCall(PetscFree((*mapping)->info_numprocs));
675268a049cSStefano Zampini   if ((*mapping)->info_indices) {
676268a049cSStefano Zampini     PetscInt i;
677268a049cSStefano Zampini 
6789566063dSJacob Faibussowitsch     PetscCall(PetscFree(((*mapping)->info_indices)[0]));
679268a049cSStefano Zampini     for (i=1; i<(*mapping)->info_nproc; i++) {
6809566063dSJacob Faibussowitsch       PetscCall(PetscFree(((*mapping)->info_indices)[i]));
681268a049cSStefano Zampini     }
6829566063dSJacob Faibussowitsch     PetscCall(PetscFree((*mapping)->info_indices));
683268a049cSStefano Zampini   }
6841bd0b88eSStefano Zampini   if ((*mapping)->info_nodei) {
6859566063dSJacob Faibussowitsch     PetscCall(PetscFree(((*mapping)->info_nodei)[0]));
6861bd0b88eSStefano Zampini   }
6879566063dSJacob Faibussowitsch   PetscCall(PetscFree2((*mapping)->info_nodec,(*mapping)->info_nodei));
688413f72f0SBarry Smith   if ((*mapping)->ops->destroy) {
6899566063dSJacob Faibussowitsch     PetscCall((*(*mapping)->ops->destroy)(*mapping));
690413f72f0SBarry Smith   }
6919566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(mapping));
6924c8fdceaSLisandro Dalcin   *mapping = NULL;
6933a40ed3dSBarry Smith   PetscFunctionReturn(0);
69490f02eecSBarry Smith }
69590f02eecSBarry Smith 
69690f02eecSBarry Smith /*@
6973acfe500SLois Curfman McInnes     ISLocalToGlobalMappingApplyIS - Creates from an IS in the local numbering
6983acfe500SLois Curfman McInnes     a new index set using the global numbering defined in an ISLocalToGlobalMapping
6993acfe500SLois Curfman McInnes     context.
70090f02eecSBarry Smith 
7014cb36875SStefano Zampini     Collective on is
702b9cd556bSLois Curfman McInnes 
70390f02eecSBarry Smith     Input Parameters:
704b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
705b9cd556bSLois Curfman McInnes -   is - index set in local numbering
70690f02eecSBarry Smith 
70790f02eecSBarry Smith     Output Parameters:
70890f02eecSBarry Smith .   newis - index set in global numbering
70990f02eecSBarry Smith 
7104cb36875SStefano Zampini     Notes:
7114cb36875SStefano Zampini     The output IS will have the same communicator of the input IS.
7124cb36875SStefano Zampini 
713a997ad1aSLois Curfman McInnes     Level: advanced
714a997ad1aSLois Curfman McInnes 
71590f02eecSBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),
716d4bb536fSBarry Smith           ISLocalToGlobalMappingDestroy(), ISGlobalToLocalMappingApply()
71790f02eecSBarry Smith @*/
7187087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping,IS is,IS *newis)
71990f02eecSBarry Smith {
720e24637baSBarry Smith   PetscInt       n,*idxout;
7215d0c19d7SBarry Smith   const PetscInt *idxin;
7223a40ed3dSBarry Smith 
7233a40ed3dSBarry Smith   PetscFunctionBegin;
7240700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
7250700a824SBarry Smith   PetscValidHeaderSpecific(is,IS_CLASSID,2);
7264482741eSBarry Smith   PetscValidPointer(newis,3);
72790f02eecSBarry Smith 
7289566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is,&n));
7299566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is,&idxin));
7309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n,&idxout));
7319566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(mapping,n,idxin,idxout));
7329566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is,&idxin));
7339566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)is),n,idxout,PETSC_OWN_POINTER,newis));
7343a40ed3dSBarry Smith   PetscFunctionReturn(0);
73590f02eecSBarry Smith }
73690f02eecSBarry Smith 
737b89cb25eSSatish Balay /*@
7383acfe500SLois Curfman McInnes    ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering
7393acfe500SLois Curfman McInnes    and converts them to the global numbering.
74090f02eecSBarry Smith 
741b9cd556bSLois Curfman McInnes    Not collective
742b9cd556bSLois Curfman McInnes 
743bb25748dSBarry Smith    Input Parameters:
744b9cd556bSLois Curfman McInnes +  mapping - the local to global mapping context
745bb25748dSBarry Smith .  N - number of integers
746b9cd556bSLois Curfman McInnes -  in - input indices in local numbering
747bb25748dSBarry Smith 
748bb25748dSBarry Smith    Output Parameter:
749bb25748dSBarry Smith .  out - indices in global numbering
750bb25748dSBarry Smith 
751b9cd556bSLois Curfman McInnes    Notes:
752b9cd556bSLois Curfman McInnes    The in and out array parameters may be identical.
753d4bb536fSBarry Smith 
754a997ad1aSLois Curfman McInnes    Level: advanced
755a997ad1aSLois Curfman McInnes 
75645b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
7570752156aSBarry Smith           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
758d4bb536fSBarry Smith           AOPetscToApplication(), ISGlobalToLocalMappingApply()
759bb25748dSBarry Smith 
760afcb2eb5SJed Brown @*/
761afcb2eb5SJed Brown PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
762afcb2eb5SJed Brown {
763cbc1caf0SMatthew G. Knepley   PetscInt i,bs,Nmax;
76445b6f7e9SBarry Smith 
76545b6f7e9SBarry Smith   PetscFunctionBegin;
766cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
767cbc1caf0SMatthew G. Knepley   bs   = mapping->bs;
768cbc1caf0SMatthew G. Knepley   Nmax = bs*mapping->n;
76945b6f7e9SBarry Smith   if (bs == 1) {
770cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
77145b6f7e9SBarry Smith     for (i=0; i<N; i++) {
77245b6f7e9SBarry Smith       if (in[i] < 0) {
77345b6f7e9SBarry Smith         out[i] = in[i];
77445b6f7e9SBarry Smith         continue;
77545b6f7e9SBarry Smith       }
776*08401ef6SPierre 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);
77745b6f7e9SBarry Smith       out[i] = idx[in[i]];
77845b6f7e9SBarry Smith     }
77945b6f7e9SBarry Smith   } else {
780cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
78145b6f7e9SBarry Smith     for (i=0; i<N; i++) {
78245b6f7e9SBarry Smith       if (in[i] < 0) {
78345b6f7e9SBarry Smith         out[i] = in[i];
78445b6f7e9SBarry Smith         continue;
78545b6f7e9SBarry Smith       }
786*08401ef6SPierre 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);
78745b6f7e9SBarry Smith       out[i] = idx[in[i]/bs]*bs + (in[i] % bs);
78845b6f7e9SBarry Smith     }
78945b6f7e9SBarry Smith   }
79045b6f7e9SBarry Smith   PetscFunctionReturn(0);
79145b6f7e9SBarry Smith }
79245b6f7e9SBarry Smith 
79345b6f7e9SBarry Smith /*@
7946006e8d2SBarry Smith    ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local block numbering and converts them to the global block numbering
79545b6f7e9SBarry Smith 
79645b6f7e9SBarry Smith    Not collective
79745b6f7e9SBarry Smith 
79845b6f7e9SBarry Smith    Input Parameters:
79945b6f7e9SBarry Smith +  mapping - the local to global mapping context
80045b6f7e9SBarry Smith .  N - number of integers
8016006e8d2SBarry Smith -  in - input indices in local block numbering
80245b6f7e9SBarry Smith 
80345b6f7e9SBarry Smith    Output Parameter:
8046006e8d2SBarry Smith .  out - indices in global block numbering
80545b6f7e9SBarry Smith 
80645b6f7e9SBarry Smith    Notes:
80745b6f7e9SBarry Smith    The in and out array parameters may be identical.
80845b6f7e9SBarry Smith 
8096006e8d2SBarry Smith    Example:
8106006e8d2SBarry 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
8116006e8d2SBarry Smith      (the first block) would produce 0 and the mapping applied to 1 (the second block) would produce 3.
8126006e8d2SBarry Smith 
81345b6f7e9SBarry Smith    Level: advanced
81445b6f7e9SBarry Smith 
81545b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
81645b6f7e9SBarry Smith           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
81745b6f7e9SBarry Smith           AOPetscToApplication(), ISGlobalToLocalMappingApply()
81845b6f7e9SBarry Smith 
81945b6f7e9SBarry Smith @*/
82045b6f7e9SBarry Smith PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
82145b6f7e9SBarry Smith {
8228a1f772fSStefano Zampini   PetscInt       i,Nmax;
8238a1f772fSStefano Zampini   const PetscInt *idx;
824d4bb536fSBarry Smith 
825a0d79125SStefano Zampini   PetscFunctionBegin;
826a0d79125SStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
8278a1f772fSStefano Zampini   Nmax = mapping->n;
8288a1f772fSStefano Zampini   idx = mapping->indices;
829afcb2eb5SJed Brown   for (i=0; i<N; i++) {
830afcb2eb5SJed Brown     if (in[i] < 0) {
831afcb2eb5SJed Brown       out[i] = in[i];
832afcb2eb5SJed Brown       continue;
833afcb2eb5SJed Brown     }
834*08401ef6SPierre 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);
835afcb2eb5SJed Brown     out[i] = idx[in[i]];
836afcb2eb5SJed Brown   }
837afcb2eb5SJed Brown   PetscFunctionReturn(0);
838afcb2eb5SJed Brown }
839d4bb536fSBarry Smith 
8407e99dc12SLawrence Mitchell /*@
841a997ad1aSLois Curfman McInnes     ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers
842a997ad1aSLois Curfman McInnes     specified with a global numbering.
843d4bb536fSBarry Smith 
844b9cd556bSLois Curfman McInnes     Not collective
845b9cd556bSLois Curfman McInnes 
846d4bb536fSBarry Smith     Input Parameters:
847b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
8480040bde1SJunchao Zhang .   type - IS_GTOLM_MASK - maps global indices with no local value to -1 in the output list (i.e., mask them)
849d4bb536fSBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
850d4bb536fSBarry Smith .   n - number of global indices to map
851b9cd556bSLois Curfman McInnes -   idx - global indices to map
852d4bb536fSBarry Smith 
853d4bb536fSBarry Smith     Output Parameters:
854b9cd556bSLois Curfman McInnes +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
855b9cd556bSLois Curfman McInnes -   idxout - local index of each global index, one must pass in an array long enough
856e182c471SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApply() with
8570298fd71SBarry Smith              idxout == NULL to determine the required length (returned in nout)
858e182c471SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApply()
859e182c471SBarry Smith              a second time to set the values.
860d4bb536fSBarry Smith 
861b9cd556bSLois Curfman McInnes     Notes:
8620298fd71SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
863d4bb536fSBarry Smith 
8649a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
865413f72f0SBarry 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.
866413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
8670f5bd95cSBarry Smith 
868a997ad1aSLois Curfman McInnes     Level: advanced
869a997ad1aSLois Curfman McInnes 
87032fd6b96SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
87132fd6b96SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
87232fd6b96SBarry Smith 
8739d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),
874413f72f0SBarry Smith           ISLocalToGlobalMappingDestroy()
875d4bb536fSBarry Smith @*/
876413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
877d4bb536fSBarry Smith {
8789d90f715SBarry Smith   PetscFunctionBegin;
8799d90f715SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
880413f72f0SBarry Smith   if (!mapping->data) {
8819566063dSJacob Faibussowitsch     PetscCall(ISGlobalToLocalMappingSetUp(mapping));
8829d90f715SBarry Smith   }
8839566063dSJacob Faibussowitsch   PetscCall((*mapping->ops->globaltolocalmappingapply)(mapping,type,n,idx,nout,idxout));
8849d90f715SBarry Smith   PetscFunctionReturn(0);
8859d90f715SBarry Smith }
8869d90f715SBarry Smith 
887d4fe737eSStefano Zampini /*@
888d4fe737eSStefano Zampini     ISGlobalToLocalMappingApplyIS - Creates from an IS in the global numbering
889d4fe737eSStefano Zampini     a new index set using the local numbering defined in an ISLocalToGlobalMapping
890d4fe737eSStefano Zampini     context.
891d4fe737eSStefano Zampini 
892d4fe737eSStefano Zampini     Not collective
893d4fe737eSStefano Zampini 
894d4fe737eSStefano Zampini     Input Parameters:
895d4fe737eSStefano Zampini +   mapping - mapping between local and global numbering
8960040bde1SJunchao Zhang .   type - IS_GTOLM_MASK - maps global indices with no local value to -1 in the output list (i.e., mask them)
8972785b321SStefano Zampini            IS_GTOLM_DROP - drops the indices with no local value from the output list
898d4fe737eSStefano Zampini -   is - index set in global numbering
899d4fe737eSStefano Zampini 
900d4fe737eSStefano Zampini     Output Parameters:
901d4fe737eSStefano Zampini .   newis - index set in local numbering
902d4fe737eSStefano Zampini 
9034cb36875SStefano Zampini     Notes:
9044cb36875SStefano Zampini     The output IS will be sequential, as it encodes a purely local operation
9054cb36875SStefano Zampini 
906d4fe737eSStefano Zampini     Level: advanced
907d4fe737eSStefano Zampini 
908d4fe737eSStefano Zampini .seealso: ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(),
909d4fe737eSStefano Zampini           ISLocalToGlobalMappingDestroy()
910d4fe737eSStefano Zampini @*/
911413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,IS is,IS *newis)
912d4fe737eSStefano Zampini {
913d4fe737eSStefano Zampini   PetscInt       n,nout,*idxout;
914d4fe737eSStefano Zampini   const PetscInt *idxin;
915d4fe737eSStefano Zampini 
916d4fe737eSStefano Zampini   PetscFunctionBegin;
917d4fe737eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
918d4fe737eSStefano Zampini   PetscValidHeaderSpecific(is,IS_CLASSID,3);
919d4fe737eSStefano Zampini   PetscValidPointer(newis,4);
920d4fe737eSStefano Zampini 
9219566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is,&n));
9229566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is,&idxin));
923d4fe737eSStefano Zampini   if (type == IS_GTOLM_MASK) {
9249566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n,&idxout));
925d4fe737eSStefano Zampini   } else {
9269566063dSJacob Faibussowitsch     PetscCall(ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,NULL));
9279566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nout,&idxout));
928d4fe737eSStefano Zampini   }
9299566063dSJacob Faibussowitsch   PetscCall(ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,idxout));
9309566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is,&idxin));
9319566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,nout,idxout,PETSC_OWN_POINTER,newis));
932d4fe737eSStefano Zampini   PetscFunctionReturn(0);
933d4fe737eSStefano Zampini }
934d4fe737eSStefano Zampini 
9359d90f715SBarry Smith /*@
9369d90f715SBarry Smith     ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers
9379d90f715SBarry Smith     specified with a block global numbering.
9389d90f715SBarry Smith 
9399d90f715SBarry Smith     Not collective
9409d90f715SBarry Smith 
9419d90f715SBarry Smith     Input Parameters:
9429d90f715SBarry Smith +   mapping - mapping between local and global numbering
9430040bde1SJunchao Zhang .   type - IS_GTOLM_MASK - maps global indices with no local value to -1 in the output list (i.e., mask them)
9449d90f715SBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
9459d90f715SBarry Smith .   n - number of global indices to map
9469d90f715SBarry Smith -   idx - global indices to map
9479d90f715SBarry Smith 
9489d90f715SBarry Smith     Output Parameters:
9499d90f715SBarry Smith +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
9509d90f715SBarry Smith -   idxout - local index of each global index, one must pass in an array long enough
9519d90f715SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApplyBlock() with
9529d90f715SBarry Smith              idxout == NULL to determine the required length (returned in nout)
9539d90f715SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApplyBlock()
9549d90f715SBarry Smith              a second time to set the values.
9559d90f715SBarry Smith 
9569d90f715SBarry Smith     Notes:
9579d90f715SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
9589d90f715SBarry Smith 
9599a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
960413f72f0SBarry 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.
961413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
9629d90f715SBarry Smith 
9639d90f715SBarry Smith     Level: advanced
9649d90f715SBarry Smith 
9659d90f715SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
9669d90f715SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
9679d90f715SBarry Smith 
9689d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(),
969413f72f0SBarry Smith           ISLocalToGlobalMappingDestroy()
9709d90f715SBarry Smith @*/
971413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,
9729d90f715SBarry Smith                                                  PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
9739d90f715SBarry Smith {
9743a40ed3dSBarry Smith   PetscFunctionBegin;
9750700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
976413f72f0SBarry Smith   if (!mapping->data) {
9779566063dSJacob Faibussowitsch     PetscCall(ISGlobalToLocalMappingSetUp(mapping));
978d4bb536fSBarry Smith   }
9799566063dSJacob Faibussowitsch   PetscCall((*mapping->ops->globaltolocalmappingapplyblock)(mapping,type,n,idx,nout,idxout));
9803a40ed3dSBarry Smith   PetscFunctionReturn(0);
981d4bb536fSBarry Smith }
98290f02eecSBarry Smith 
98389d82c54SBarry Smith /*@C
9846a818285SBarry Smith     ISLocalToGlobalMappingGetBlockInfo - Gets the neighbor information for each processor and
98589d82c54SBarry Smith      each index shared by more than one processor
98689d82c54SBarry Smith 
98789d82c54SBarry Smith     Collective on ISLocalToGlobalMapping
98889d82c54SBarry Smith 
989f899ff85SJose E. Roman     Input Parameter:
99089d82c54SBarry Smith .   mapping - the mapping from local to global indexing
99189d82c54SBarry Smith 
992d8d19677SJose E. Roman     Output Parameters:
99389d82c54SBarry Smith +   nproc - number of processors that are connected to this one
99489d82c54SBarry Smith .   proc - neighboring processors
99507b52d57SBarry Smith .   numproc - number of indices for each subdomain (processor)
9963463a7baSJed Brown -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
99789d82c54SBarry Smith 
99889d82c54SBarry Smith     Level: advanced
99989d82c54SBarry Smith 
10002cfcea29SBarry Smith     Fortran Usage:
10012cfcea29SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
10022cfcea29SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
10032cfcea29SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
10042cfcea29SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
10052cfcea29SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
10062cfcea29SBarry Smith 
100707b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
100807b52d57SBarry Smith           ISLocalToGlobalMappingRestoreInfo()
100989d82c54SBarry Smith @*/
10106a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
101189d82c54SBarry Smith {
1012268a049cSStefano Zampini   PetscFunctionBegin;
1013268a049cSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1014268a049cSStefano Zampini   if (mapping->info_cached) {
1015268a049cSStefano Zampini     *nproc    = mapping->info_nproc;
1016268a049cSStefano Zampini     *procs    = mapping->info_procs;
1017268a049cSStefano Zampini     *numprocs = mapping->info_numprocs;
1018268a049cSStefano Zampini     *indices  = mapping->info_indices;
1019268a049cSStefano Zampini   } else {
10209566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetBlockInfo_Private(mapping,nproc,procs,numprocs,indices));
1021268a049cSStefano Zampini   }
1022268a049cSStefano Zampini   PetscFunctionReturn(0);
1023268a049cSStefano Zampini }
1024268a049cSStefano Zampini 
1025268a049cSStefano Zampini static PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
1026268a049cSStefano Zampini {
102797f1f81fSBarry Smith   PetscMPIInt    size,rank,tag1,tag2,tag3,*len,*source,imdex;
102832dcc486SBarry Smith   PetscInt       i,n = mapping->n,Ng,ng,max = 0,*lindices = mapping->indices;
102932dcc486SBarry Smith   PetscInt       *nprocs,*owner,nsends,*sends,j,*starts,nmax,nrecvs,*recvs,proc;
103097f1f81fSBarry Smith   PetscInt       cnt,scale,*ownedsenders,*nownedsenders,rstart,nowned;
103132dcc486SBarry Smith   PetscInt       node,nownedm,nt,*sends2,nsends2,*starts2,*lens2,*dest,nrecvs2,*starts3,*recvs2,k,*bprocs,*tmp;
103232dcc486SBarry Smith   PetscInt       first_procs,first_numprocs,*first_indices;
103389d82c54SBarry Smith   MPI_Request    *recv_waits,*send_waits;
103430dcb7c9SBarry Smith   MPI_Status     recv_status,*send_status,*recv_statuses;
1035ce94432eSBarry Smith   MPI_Comm       comm;
1036ace3abfcSBarry Smith   PetscBool      debug = PETSC_FALSE;
103789d82c54SBarry Smith 
103889d82c54SBarry Smith   PetscFunctionBegin;
10399566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mapping,&comm));
10409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
10419566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
104224cf384cSBarry Smith   if (size == 1) {
104324cf384cSBarry Smith     *nproc         = 0;
10440298fd71SBarry Smith     *procs         = NULL;
10459566063dSJacob Faibussowitsch     PetscCall(PetscNew(numprocs));
10461e2105dcSBarry Smith     (*numprocs)[0] = 0;
10479566063dSJacob Faibussowitsch     PetscCall(PetscNew(indices));
10480298fd71SBarry Smith     (*indices)[0]  = NULL;
1049268a049cSStefano Zampini     /* save info for reuse */
1050268a049cSStefano Zampini     mapping->info_nproc = *nproc;
1051268a049cSStefano Zampini     mapping->info_procs = *procs;
1052268a049cSStefano Zampini     mapping->info_numprocs = *numprocs;
1053268a049cSStefano Zampini     mapping->info_indices = *indices;
1054268a049cSStefano Zampini     mapping->info_cached = PETSC_TRUE;
105524cf384cSBarry Smith     PetscFunctionReturn(0);
105624cf384cSBarry Smith   }
105724cf384cSBarry Smith 
10589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)mapping)->options,NULL,"-islocaltoglobalmappinggetinfo_debug",&debug,NULL));
105907b52d57SBarry Smith 
10603677ff5aSBarry Smith   /*
10616a818285SBarry Smith     Notes on ISLocalToGlobalMappingGetBlockInfo
10623677ff5aSBarry Smith 
10633677ff5aSBarry Smith     globally owned node - the nodes that have been assigned to this processor in global
10643677ff5aSBarry Smith            numbering, just for this routine.
10653677ff5aSBarry Smith 
10663677ff5aSBarry Smith     nontrivial globally owned node - node assigned to this processor that is on a subdomain
10673677ff5aSBarry Smith            boundary (i.e. is has more than one local owner)
10683677ff5aSBarry Smith 
10693677ff5aSBarry Smith     locally owned node - node that exists on this processors subdomain
10703677ff5aSBarry Smith 
10713677ff5aSBarry Smith     nontrivial locally owned node - node that is not in the interior (i.e. has more than one
10723677ff5aSBarry Smith            local subdomain
10733677ff5aSBarry Smith   */
10749566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetNewTag((PetscObject)mapping,&tag1));
10759566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetNewTag((PetscObject)mapping,&tag2));
10769566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetNewTag((PetscObject)mapping,&tag3));
107789d82c54SBarry Smith 
107889d82c54SBarry Smith   for (i=0; i<n; i++) {
107989d82c54SBarry Smith     if (lindices[i] > max) max = lindices[i];
108089d82c54SBarry Smith   }
10811c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&max,&Ng,1,MPIU_INT,MPI_MAX,comm));
108278058e43SBarry Smith   Ng++;
10839566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
10849566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
1085bc8ff85bSBarry Smith   scale  = Ng/size + 1;
1086a2e34c3dSBarry Smith   ng     = scale; if (rank == size-1) ng = Ng - scale*(size-1); ng = PetscMax(1,ng);
1087caba0dd0SBarry Smith   rstart = scale*rank;
108889d82c54SBarry Smith 
108989d82c54SBarry Smith   /* determine ownership ranges of global indices */
10909566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2*size,&nprocs));
10919566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(nprocs,2*size));
109289d82c54SBarry Smith 
109389d82c54SBarry Smith   /* determine owners of each local node  */
10949566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n,&owner));
109589d82c54SBarry Smith   for (i=0; i<n; i++) {
10963677ff5aSBarry Smith     proc             = lindices[i]/scale; /* processor that globally owns this index */
109727c402fcSBarry Smith     nprocs[2*proc+1] = 1;                 /* processor globally owns at least one of ours */
10983677ff5aSBarry Smith     owner[i]         = proc;
109927c402fcSBarry Smith     nprocs[2*proc]++;                     /* count of how many that processor globally owns of ours */
110089d82c54SBarry Smith   }
110127c402fcSBarry Smith   nsends = 0; for (i=0; i<size; i++) nsends += nprocs[2*i+1];
11029566063dSJacob Faibussowitsch   PetscCall(PetscInfo(mapping,"Number of global owners for my local data %" PetscInt_FMT "\n",nsends));
110389d82c54SBarry Smith 
110489d82c54SBarry Smith   /* inform other processors of number of messages and max length*/
11059566063dSJacob Faibussowitsch   PetscCall(PetscMaxSum(comm,nprocs,&nmax,&nrecvs));
11069566063dSJacob Faibussowitsch   PetscCall(PetscInfo(mapping,"Number of local owners for my global data %" PetscInt_FMT "\n",nrecvs));
110789d82c54SBarry Smith 
110889d82c54SBarry Smith   /* post receives for owned rows */
11099566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((2*nrecvs+1)*(nmax+1),&recvs));
11109566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs+1,&recv_waits));
111189d82c54SBarry Smith   for (i=0; i<nrecvs; i++) {
11129566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Irecv(recvs+2*nmax*i,2*nmax,MPIU_INT,MPI_ANY_SOURCE,tag1,comm,recv_waits+i));
111389d82c54SBarry Smith   }
111489d82c54SBarry Smith 
111589d82c54SBarry Smith   /* pack messages containing lists of local nodes to owners */
11169566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2*n+1,&sends));
11179566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size+1,&starts));
111889d82c54SBarry Smith   starts[0] = 0;
1119f6e5521dSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
112089d82c54SBarry Smith   for (i=0; i<n; i++) {
112189d82c54SBarry Smith     sends[starts[owner[i]]++] = lindices[i];
112230dcb7c9SBarry Smith     sends[starts[owner[i]]++] = i;
112389d82c54SBarry Smith   }
11249566063dSJacob Faibussowitsch   PetscCall(PetscFree(owner));
112589d82c54SBarry Smith   starts[0] = 0;
1126f6e5521dSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
112789d82c54SBarry Smith 
112889d82c54SBarry Smith   /* send the messages */
11299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends+1,&send_waits));
11309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends+1,&dest));
113189d82c54SBarry Smith   cnt = 0;
113289d82c54SBarry Smith   for (i=0; i<size; i++) {
113327c402fcSBarry Smith     if (nprocs[2*i]) {
11349566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Isend(sends+starts[i],2*nprocs[2*i],MPIU_INT,i,tag1,comm,send_waits+cnt));
113530dcb7c9SBarry Smith       dest[cnt] = i;
113689d82c54SBarry Smith       cnt++;
113789d82c54SBarry Smith     }
113889d82c54SBarry Smith   }
11399566063dSJacob Faibussowitsch   PetscCall(PetscFree(starts));
114089d82c54SBarry Smith 
114189d82c54SBarry Smith   /* wait on receives */
11429566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs+1,&source));
11439566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs+1,&len));
114489d82c54SBarry Smith   cnt  = nrecvs;
11459566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(ng+1,&nownedsenders));
114689d82c54SBarry Smith   while (cnt) {
11479566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status));
114889d82c54SBarry Smith     /* unpack receives into our local space */
11499566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Get_count(&recv_status,MPIU_INT,&len[imdex]));
115089d82c54SBarry Smith     source[imdex] = recv_status.MPI_SOURCE;
115130dcb7c9SBarry Smith     len[imdex]    = len[imdex]/2;
1152caba0dd0SBarry Smith     /* count how many local owners for each of my global owned indices */
115330dcb7c9SBarry Smith     for (i=0; i<len[imdex]; i++) nownedsenders[recvs[2*imdex*nmax+2*i]-rstart]++;
115489d82c54SBarry Smith     cnt--;
115589d82c54SBarry Smith   }
11569566063dSJacob Faibussowitsch   PetscCall(PetscFree(recv_waits));
115789d82c54SBarry Smith 
115830dcb7c9SBarry Smith   /* count how many globally owned indices are on an edge multiplied by how many processors own them. */
1159bc8ff85bSBarry Smith   nowned  = 0;
1160bc8ff85bSBarry Smith   nownedm = 0;
1161bc8ff85bSBarry Smith   for (i=0; i<ng; i++) {
1162bc8ff85bSBarry Smith     if (nownedsenders[i] > 1) {nownedm += nownedsenders[i]; nowned++;}
1163bc8ff85bSBarry Smith   }
1164bc8ff85bSBarry Smith 
1165bc8ff85bSBarry Smith   /* create single array to contain rank of all local owners of each globally owned index */
11669566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nownedm+1,&ownedsenders));
11679566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ng+1,&starts));
1168bc8ff85bSBarry Smith   starts[0] = 0;
1169bc8ff85bSBarry Smith   for (i=1; i<ng; i++) {
1170bc8ff85bSBarry Smith     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
1171bc8ff85bSBarry Smith     else starts[i] = starts[i-1];
1172bc8ff85bSBarry Smith   }
1173bc8ff85bSBarry Smith 
117430dcb7c9SBarry Smith   /* for each nontrival globally owned node list all arriving processors */
1175bc8ff85bSBarry Smith   for (i=0; i<nrecvs; i++) {
1176bc8ff85bSBarry Smith     for (j=0; j<len[i]; j++) {
117730dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
1178f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) ownedsenders[starts[node]++] = source[i];
1179bc8ff85bSBarry Smith     }
1180bc8ff85bSBarry Smith   }
1181bc8ff85bSBarry Smith 
118207b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
118330dcb7c9SBarry Smith     starts[0] = 0;
118430dcb7c9SBarry Smith     for (i=1; i<ng; i++) {
118530dcb7c9SBarry Smith       if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
118630dcb7c9SBarry Smith       else starts[i] = starts[i-1];
118730dcb7c9SBarry Smith     }
118830dcb7c9SBarry Smith     for (i=0; i<ng; i++) {
118930dcb7c9SBarry Smith       if (nownedsenders[i] > 1) {
11909566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(comm,"[%d] global node %" PetscInt_FMT " local owner processors: ",rank,i+rstart));
119130dcb7c9SBarry Smith         for (j=0; j<nownedsenders[i]; j++) {
11929566063dSJacob Faibussowitsch           PetscCall(PetscSynchronizedPrintf(comm,"%" PetscInt_FMT " ",ownedsenders[starts[i]+j]));
119330dcb7c9SBarry Smith         }
11949566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(comm,"\n"));
119530dcb7c9SBarry Smith       }
119630dcb7c9SBarry Smith     }
11979566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm,PETSC_STDOUT));
119807b52d57SBarry Smith   } /* -----------------------------------  */
119930dcb7c9SBarry Smith 
12003677ff5aSBarry Smith   /* wait on original sends */
12013a96401aSBarry Smith   if (nsends) {
12029566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nsends,&send_status));
12039566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nsends,send_waits,send_status));
12049566063dSJacob Faibussowitsch     PetscCall(PetscFree(send_status));
12053a96401aSBarry Smith   }
12069566063dSJacob Faibussowitsch   PetscCall(PetscFree(send_waits));
12079566063dSJacob Faibussowitsch   PetscCall(PetscFree(sends));
12089566063dSJacob Faibussowitsch   PetscCall(PetscFree(nprocs));
12093677ff5aSBarry Smith 
12103677ff5aSBarry Smith   /* pack messages to send back to local owners */
121130dcb7c9SBarry Smith   starts[0] = 0;
121230dcb7c9SBarry Smith   for (i=1; i<ng; i++) {
121330dcb7c9SBarry Smith     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
121430dcb7c9SBarry Smith     else starts[i] = starts[i-1];
121530dcb7c9SBarry Smith   }
121630dcb7c9SBarry Smith   nsends2 = nrecvs;
12179566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends2+1,&nprocs)); /* length of each message */
121830dcb7c9SBarry Smith   for (i=0; i<nrecvs; i++) {
121930dcb7c9SBarry Smith     nprocs[i] = 1;
122030dcb7c9SBarry Smith     for (j=0; j<len[i]; j++) {
122130dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
1222f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) nprocs[i] += 2 + nownedsenders[node];
122330dcb7c9SBarry Smith     }
122430dcb7c9SBarry Smith   }
1225f6e5521dSKarl Rupp   nt = 0;
1226f6e5521dSKarl Rupp   for (i=0; i<nsends2; i++) nt += nprocs[i];
1227f6e5521dSKarl Rupp 
12289566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt+1,&sends2));
12299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends2+1,&starts2));
1230f6e5521dSKarl Rupp 
1231f6e5521dSKarl Rupp   starts2[0] = 0;
1232f6e5521dSKarl Rupp   for (i=1; i<nsends2; i++) starts2[i] = starts2[i-1] + nprocs[i-1];
123330dcb7c9SBarry Smith   /*
123430dcb7c9SBarry Smith      Each message is 1 + nprocs[i] long, and consists of
123530dcb7c9SBarry Smith        (0) the number of nodes being sent back
123630dcb7c9SBarry Smith        (1) the local node number,
123730dcb7c9SBarry Smith        (2) the number of processors sharing it,
123830dcb7c9SBarry Smith        (3) the processors sharing it
123930dcb7c9SBarry Smith   */
124030dcb7c9SBarry Smith   for (i=0; i<nsends2; i++) {
124130dcb7c9SBarry Smith     cnt = 1;
124230dcb7c9SBarry Smith     sends2[starts2[i]] = 0;
124330dcb7c9SBarry Smith     for (j=0; j<len[i]; j++) {
124430dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
124530dcb7c9SBarry Smith       if (nownedsenders[node] > 1) {
124630dcb7c9SBarry Smith         sends2[starts2[i]]++;
124730dcb7c9SBarry Smith         sends2[starts2[i]+cnt++] = recvs[2*i*nmax+2*j+1];
124830dcb7c9SBarry Smith         sends2[starts2[i]+cnt++] = nownedsenders[node];
12499566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(&sends2[starts2[i]+cnt],&ownedsenders[starts[node]],nownedsenders[node]));
125030dcb7c9SBarry Smith         cnt += nownedsenders[node];
125130dcb7c9SBarry Smith       }
125230dcb7c9SBarry Smith     }
125330dcb7c9SBarry Smith   }
125430dcb7c9SBarry Smith 
125530dcb7c9SBarry Smith   /* receive the message lengths */
125630dcb7c9SBarry Smith   nrecvs2 = nsends;
12579566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs2+1,&lens2));
12589566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs2+1,&starts3));
12599566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs2+1,&recv_waits));
126030dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
12619566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Irecv(&lens2[i],1,MPIU_INT,dest[i],tag2,comm,recv_waits+i));
126230dcb7c9SBarry Smith   }
1263d44834fbSBarry Smith 
12648a8e0b3aSBarry Smith   /* send the message lengths */
12658a8e0b3aSBarry Smith   for (i=0; i<nsends2; i++) {
12669566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send(&nprocs[i],1,MPIU_INT,source[i],tag2,comm));
12678a8e0b3aSBarry Smith   }
12688a8e0b3aSBarry Smith 
1269d44834fbSBarry Smith   /* wait on receives of lens */
12700c468ba9SBarry Smith   if (nrecvs2) {
12719566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrecvs2,&recv_statuses));
12729566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrecvs2,recv_waits,recv_statuses));
12739566063dSJacob Faibussowitsch     PetscCall(PetscFree(recv_statuses));
12740c468ba9SBarry Smith   }
12759566063dSJacob Faibussowitsch   PetscCall(PetscFree(recv_waits));
1276d44834fbSBarry Smith 
127730dcb7c9SBarry Smith   starts3[0] = 0;
1278d44834fbSBarry Smith   nt         = 0;
127930dcb7c9SBarry Smith   for (i=0; i<nrecvs2-1; i++) {
128030dcb7c9SBarry Smith     starts3[i+1] = starts3[i] + lens2[i];
1281d44834fbSBarry Smith     nt          += lens2[i];
128230dcb7c9SBarry Smith   }
128376466f69SStefano Zampini   if (nrecvs2) nt += lens2[nrecvs2-1];
1284d44834fbSBarry Smith 
12859566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt+1,&recvs2));
12869566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nrecvs2+1,&recv_waits));
128752b72c4aSBarry Smith   for (i=0; i<nrecvs2; i++) {
12889566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Irecv(recvs2+starts3[i],lens2[i],MPIU_INT,dest[i],tag3,comm,recv_waits+i));
128930dcb7c9SBarry Smith   }
129030dcb7c9SBarry Smith 
129130dcb7c9SBarry Smith   /* send the messages */
12929566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nsends2+1,&send_waits));
129330dcb7c9SBarry Smith   for (i=0; i<nsends2; i++) {
12949566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Isend(sends2+starts2[i],nprocs[i],MPIU_INT,source[i],tag3,comm,send_waits+i));
129530dcb7c9SBarry Smith   }
129630dcb7c9SBarry Smith 
129730dcb7c9SBarry Smith   /* wait on receives */
12980c468ba9SBarry Smith   if (nrecvs2) {
12999566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrecvs2,&recv_statuses));
13009566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nrecvs2,recv_waits,recv_statuses));
13019566063dSJacob Faibussowitsch     PetscCall(PetscFree(recv_statuses));
13020c468ba9SBarry Smith   }
13039566063dSJacob Faibussowitsch   PetscCall(PetscFree(recv_waits));
13049566063dSJacob Faibussowitsch   PetscCall(PetscFree(nprocs));
130530dcb7c9SBarry Smith 
130607b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
130730dcb7c9SBarry Smith     cnt = 0;
130830dcb7c9SBarry Smith     for (i=0; i<nrecvs2; i++) {
130930dcb7c9SBarry Smith       nt = recvs2[cnt++];
131030dcb7c9SBarry Smith       for (j=0; j<nt; j++) {
13119566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(comm,"[%d] local node %" PetscInt_FMT " number of subdomains %" PetscInt_FMT ": ",rank,recvs2[cnt],recvs2[cnt+1]));
131230dcb7c9SBarry Smith         for (k=0; k<recvs2[cnt+1]; k++) {
13139566063dSJacob Faibussowitsch           PetscCall(PetscSynchronizedPrintf(comm,"%" PetscInt_FMT " ",recvs2[cnt+2+k]));
131430dcb7c9SBarry Smith         }
131530dcb7c9SBarry Smith         cnt += 2 + recvs2[cnt+1];
13169566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(comm,"\n"));
131730dcb7c9SBarry Smith       }
131830dcb7c9SBarry Smith     }
13199566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm,PETSC_STDOUT));
132007b52d57SBarry Smith   } /* -----------------------------------  */
132130dcb7c9SBarry Smith 
132230dcb7c9SBarry Smith   /* count number subdomains for each local node */
13239566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(size,&nprocs));
132430dcb7c9SBarry Smith   cnt  = 0;
132530dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
132630dcb7c9SBarry Smith     nt = recvs2[cnt++];
132730dcb7c9SBarry Smith     for (j=0; j<nt; j++) {
1328f6e5521dSKarl Rupp       for (k=0; k<recvs2[cnt+1]; k++) nprocs[recvs2[cnt+2+k]]++;
132930dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt+1];
133030dcb7c9SBarry Smith     }
133130dcb7c9SBarry Smith   }
133230dcb7c9SBarry Smith   nt = 0; for (i=0; i<size; i++) nt += (nprocs[i] > 0);
133330dcb7c9SBarry Smith   *nproc    = nt;
13349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt+1,procs));
13359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt+1,numprocs));
13369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt+1,indices));
13370298fd71SBarry Smith   for (i=0;i<nt+1;i++) (*indices)[i]=NULL;
13389566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size,&bprocs));
133930dcb7c9SBarry Smith   cnt  = 0;
134030dcb7c9SBarry Smith   for (i=0; i<size; i++) {
134130dcb7c9SBarry Smith     if (nprocs[i] > 0) {
134230dcb7c9SBarry Smith       bprocs[i]        = cnt;
134330dcb7c9SBarry Smith       (*procs)[cnt]    = i;
134430dcb7c9SBarry Smith       (*numprocs)[cnt] = nprocs[i];
13459566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nprocs[i],&(*indices)[cnt]));
134630dcb7c9SBarry Smith       cnt++;
134730dcb7c9SBarry Smith     }
134830dcb7c9SBarry Smith   }
134930dcb7c9SBarry Smith 
135030dcb7c9SBarry Smith   /* make the list of subdomains for each nontrivial local node */
13519566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(*numprocs,nt));
135230dcb7c9SBarry Smith   cnt  = 0;
135330dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
135430dcb7c9SBarry Smith     nt = recvs2[cnt++];
135530dcb7c9SBarry Smith     for (j=0; j<nt; j++) {
1356f6e5521dSKarl Rupp       for (k=0; k<recvs2[cnt+1]; k++) (*indices)[bprocs[recvs2[cnt+2+k]]][(*numprocs)[bprocs[recvs2[cnt+2+k]]]++] = recvs2[cnt];
135730dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt+1];
135830dcb7c9SBarry Smith     }
135930dcb7c9SBarry Smith   }
13609566063dSJacob Faibussowitsch   PetscCall(PetscFree(bprocs));
13619566063dSJacob Faibussowitsch   PetscCall(PetscFree(recvs2));
136230dcb7c9SBarry Smith 
136307b52d57SBarry Smith   /* sort the node indexing by their global numbers */
136407b52d57SBarry Smith   nt = *nproc;
136507b52d57SBarry Smith   for (i=0; i<nt; i++) {
13669566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1((*numprocs)[i],&tmp));
1367f6e5521dSKarl Rupp     for (j=0; j<(*numprocs)[i]; j++) tmp[j] = lindices[(*indices)[i][j]];
13689566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithArray((*numprocs)[i],tmp,(*indices)[i]));
13699566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmp));
137007b52d57SBarry Smith   }
137107b52d57SBarry Smith 
137207b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
137330dcb7c9SBarry Smith     nt = *nproc;
137430dcb7c9SBarry Smith     for (i=0; i<nt; i++) {
13759566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedPrintf(comm,"[%d] subdomain %" PetscInt_FMT " number of indices %" PetscInt_FMT ": ",rank,(*procs)[i],(*numprocs)[i]));
137630dcb7c9SBarry Smith       for (j=0; j<(*numprocs)[i]; j++) {
13779566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(comm,"%" PetscInt_FMT " ",(*indices)[i][j]));
137830dcb7c9SBarry Smith       }
13799566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedPrintf(comm,"\n"));
138030dcb7c9SBarry Smith     }
13819566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm,PETSC_STDOUT));
138207b52d57SBarry Smith   } /* -----------------------------------  */
138330dcb7c9SBarry Smith 
138430dcb7c9SBarry Smith   /* wait on sends */
138530dcb7c9SBarry Smith   if (nsends2) {
13869566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nsends2,&send_status));
13879566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(nsends2,send_waits,send_status));
13889566063dSJacob Faibussowitsch     PetscCall(PetscFree(send_status));
138930dcb7c9SBarry Smith   }
139030dcb7c9SBarry Smith 
13919566063dSJacob Faibussowitsch   PetscCall(PetscFree(starts3));
13929566063dSJacob Faibussowitsch   PetscCall(PetscFree(dest));
13939566063dSJacob Faibussowitsch   PetscCall(PetscFree(send_waits));
13943677ff5aSBarry Smith 
13959566063dSJacob Faibussowitsch   PetscCall(PetscFree(nownedsenders));
13969566063dSJacob Faibussowitsch   PetscCall(PetscFree(ownedsenders));
13979566063dSJacob Faibussowitsch   PetscCall(PetscFree(starts));
13989566063dSJacob Faibussowitsch   PetscCall(PetscFree(starts2));
13999566063dSJacob Faibussowitsch   PetscCall(PetscFree(lens2));
140089d82c54SBarry Smith 
14019566063dSJacob Faibussowitsch   PetscCall(PetscFree(source));
14029566063dSJacob Faibussowitsch   PetscCall(PetscFree(len));
14039566063dSJacob Faibussowitsch   PetscCall(PetscFree(recvs));
14049566063dSJacob Faibussowitsch   PetscCall(PetscFree(nprocs));
14059566063dSJacob Faibussowitsch   PetscCall(PetscFree(sends2));
140624cf384cSBarry Smith 
140724cf384cSBarry Smith   /* put the information about myself as the first entry in the list */
140824cf384cSBarry Smith   first_procs    = (*procs)[0];
140924cf384cSBarry Smith   first_numprocs = (*numprocs)[0];
141024cf384cSBarry Smith   first_indices  = (*indices)[0];
141124cf384cSBarry Smith   for (i=0; i<*nproc; i++) {
141224cf384cSBarry Smith     if ((*procs)[i] == rank) {
141324cf384cSBarry Smith       (*procs)[0]    = (*procs)[i];
141424cf384cSBarry Smith       (*numprocs)[0] = (*numprocs)[i];
141524cf384cSBarry Smith       (*indices)[0]  = (*indices)[i];
141624cf384cSBarry Smith       (*procs)[i]    = first_procs;
141724cf384cSBarry Smith       (*numprocs)[i] = first_numprocs;
141824cf384cSBarry Smith       (*indices)[i]  = first_indices;
141924cf384cSBarry Smith       break;
142024cf384cSBarry Smith     }
142124cf384cSBarry Smith   }
1422268a049cSStefano Zampini 
1423268a049cSStefano Zampini   /* save info for reuse */
1424268a049cSStefano Zampini   mapping->info_nproc = *nproc;
1425268a049cSStefano Zampini   mapping->info_procs = *procs;
1426268a049cSStefano Zampini   mapping->info_numprocs = *numprocs;
1427268a049cSStefano Zampini   mapping->info_indices = *indices;
1428268a049cSStefano Zampini   mapping->info_cached = PETSC_TRUE;
142989d82c54SBarry Smith   PetscFunctionReturn(0);
143089d82c54SBarry Smith }
143189d82c54SBarry Smith 
14326a818285SBarry Smith /*@C
14336a818285SBarry Smith     ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by ISLocalToGlobalMappingGetBlockInfo()
14346a818285SBarry Smith 
14356a818285SBarry Smith     Collective on ISLocalToGlobalMapping
14366a818285SBarry Smith 
1437f899ff85SJose E. Roman     Input Parameter:
14386a818285SBarry Smith .   mapping - the mapping from local to global indexing
14396a818285SBarry Smith 
1440d8d19677SJose E. Roman     Output Parameters:
14416a818285SBarry Smith +   nproc - number of processors that are connected to this one
14426a818285SBarry Smith .   proc - neighboring processors
14436a818285SBarry Smith .   numproc - number of indices for each processor
14446a818285SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
14456a818285SBarry Smith 
14466a818285SBarry Smith     Level: advanced
14476a818285SBarry Smith 
14486a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
14496a818285SBarry Smith           ISLocalToGlobalMappingGetInfo()
14506a818285SBarry Smith @*/
14516a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
14526a818285SBarry Smith {
14536a818285SBarry Smith   PetscFunctionBegin;
1454cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1455268a049cSStefano Zampini   if (mapping->info_free) {
14569566063dSJacob Faibussowitsch     PetscCall(PetscFree(*numprocs));
14576a818285SBarry Smith     if (*indices) {
1458268a049cSStefano Zampini       PetscInt i;
1459268a049cSStefano Zampini 
14609566063dSJacob Faibussowitsch       PetscCall(PetscFree((*indices)[0]));
14616a818285SBarry Smith       for (i=1; i<*nproc; i++) {
14629566063dSJacob Faibussowitsch         PetscCall(PetscFree((*indices)[i]));
14636a818285SBarry Smith       }
14649566063dSJacob Faibussowitsch       PetscCall(PetscFree(*indices));
14656a818285SBarry Smith     }
1466268a049cSStefano Zampini   }
1467268a049cSStefano Zampini   *nproc    = 0;
1468268a049cSStefano Zampini   *procs    = NULL;
1469268a049cSStefano Zampini   *numprocs = NULL;
1470268a049cSStefano Zampini   *indices  = NULL;
14716a818285SBarry Smith   PetscFunctionReturn(0);
14726a818285SBarry Smith }
14736a818285SBarry Smith 
14746a818285SBarry Smith /*@C
14756a818285SBarry Smith     ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each processor and
14766a818285SBarry Smith      each index shared by more than one processor
14776a818285SBarry Smith 
14786a818285SBarry Smith     Collective on ISLocalToGlobalMapping
14796a818285SBarry Smith 
1480f899ff85SJose E. Roman     Input Parameter:
14816a818285SBarry Smith .   mapping - the mapping from local to global indexing
14826a818285SBarry Smith 
1483d8d19677SJose E. Roman     Output Parameters:
14846a818285SBarry Smith +   nproc - number of processors that are connected to this one
14856a818285SBarry Smith .   proc - neighboring processors
14866a818285SBarry Smith .   numproc - number of indices for each subdomain (processor)
14876a818285SBarry Smith -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
14886a818285SBarry Smith 
14896a818285SBarry Smith     Level: advanced
14906a818285SBarry Smith 
14911bd0b88eSStefano Zampini     Notes: The user needs to call ISLocalToGlobalMappingRestoreInfo when the data is no longer needed.
14921bd0b88eSStefano Zampini 
14936a818285SBarry Smith     Fortran Usage:
14946a818285SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
14956a818285SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
14966a818285SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
14976a818285SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
14986a818285SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
14996a818285SBarry Smith 
15006a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
15016a818285SBarry Smith           ISLocalToGlobalMappingRestoreInfo()
15026a818285SBarry Smith @*/
15036a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
15046a818285SBarry Smith {
15058a1f772fSStefano Zampini   PetscInt       **bindices = NULL,*bnumprocs = NULL,bs,i,j,k;
15066a818285SBarry Smith 
15076a818285SBarry Smith   PetscFunctionBegin;
15086a818285SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
15098a1f772fSStefano Zampini   bs = mapping->bs;
15109566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockInfo(mapping,nproc,procs,&bnumprocs,&bindices));
1511268a049cSStefano Zampini   if (bs > 1) { /* we need to expand the cached info */
15129566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(*nproc,&*indices));
15139566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(*nproc,&*numprocs));
15146a818285SBarry Smith     for (i=0; i<*nproc; i++) {
15159566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(bs*bnumprocs[i],&(*indices)[i]));
1516268a049cSStefano Zampini       for (j=0; j<bnumprocs[i]; j++) {
15176a818285SBarry Smith         for (k=0; k<bs; k++) {
15186a818285SBarry Smith           (*indices)[i][j*bs+k] = bs*bindices[i][j] + k;
15196a818285SBarry Smith         }
15206a818285SBarry Smith       }
1521268a049cSStefano Zampini       (*numprocs)[i] = bnumprocs[i]*bs;
15226a818285SBarry Smith     }
1523268a049cSStefano Zampini     mapping->info_free = PETSC_TRUE;
1524268a049cSStefano Zampini   } else {
1525268a049cSStefano Zampini     *numprocs = bnumprocs;
1526268a049cSStefano Zampini     *indices  = bindices;
15276a818285SBarry Smith   }
15286a818285SBarry Smith   PetscFunctionReturn(0);
15296a818285SBarry Smith }
15306a818285SBarry Smith 
153107b52d57SBarry Smith /*@C
153207b52d57SBarry Smith     ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by ISLocalToGlobalMappingGetInfo()
153389d82c54SBarry Smith 
153407b52d57SBarry Smith     Collective on ISLocalToGlobalMapping
153507b52d57SBarry Smith 
1536f899ff85SJose E. Roman     Input Parameter:
153707b52d57SBarry Smith .   mapping - the mapping from local to global indexing
153807b52d57SBarry Smith 
1539d8d19677SJose E. Roman     Output Parameters:
154007b52d57SBarry Smith +   nproc - number of processors that are connected to this one
154107b52d57SBarry Smith .   proc - neighboring processors
154207b52d57SBarry Smith .   numproc - number of indices for each processor
154307b52d57SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
154407b52d57SBarry Smith 
154507b52d57SBarry Smith     Level: advanced
154607b52d57SBarry Smith 
154707b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
154807b52d57SBarry Smith           ISLocalToGlobalMappingGetInfo()
154907b52d57SBarry Smith @*/
15507087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
155107b52d57SBarry Smith {
155207b52d57SBarry Smith   PetscFunctionBegin;
15539566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreBlockInfo(mapping,nproc,procs,numprocs,indices));
155407b52d57SBarry Smith   PetscFunctionReturn(0);
155507b52d57SBarry Smith }
155686994e45SJed Brown 
155786994e45SJed Brown /*@C
15581bd0b88eSStefano Zampini     ISLocalToGlobalMappingGetNodeInfo - Gets the neighbor information for each node
15591bd0b88eSStefano Zampini 
15601bd0b88eSStefano Zampini     Collective on ISLocalToGlobalMapping
15611bd0b88eSStefano Zampini 
1562f899ff85SJose E. Roman     Input Parameter:
15631bd0b88eSStefano Zampini .   mapping - the mapping from local to global indexing
15641bd0b88eSStefano Zampini 
1565d8d19677SJose E. Roman     Output Parameters:
15661bd0b88eSStefano Zampini +   nnodes - number of local nodes (same ISLocalToGlobalMappingGetSize())
15671bd0b88eSStefano Zampini .   count - number of neighboring processors per node
15681bd0b88eSStefano Zampini -   indices - indices of processes sharing the node (sorted)
15691bd0b88eSStefano Zampini 
15701bd0b88eSStefano Zampini     Level: advanced
15711bd0b88eSStefano Zampini 
15721bd0b88eSStefano Zampini     Notes: The user needs to call ISLocalToGlobalMappingRestoreInfo when the data is no longer needed.
15731bd0b88eSStefano Zampini 
15741bd0b88eSStefano Zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
15751bd0b88eSStefano Zampini           ISLocalToGlobalMappingGetInfo(), ISLocalToGlobalMappingRestoreNodeInfo()
15761bd0b88eSStefano Zampini @*/
15771bd0b88eSStefano Zampini PetscErrorCode  ISLocalToGlobalMappingGetNodeInfo(ISLocalToGlobalMapping mapping,PetscInt *nnodes,PetscInt *count[],PetscInt **indices[])
15781bd0b88eSStefano Zampini {
15791bd0b88eSStefano Zampini   PetscInt       n;
15801bd0b88eSStefano Zampini 
15811bd0b88eSStefano Zampini   PetscFunctionBegin;
15821bd0b88eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
15839566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetSize(mapping,&n));
15841bd0b88eSStefano Zampini   if (!mapping->info_nodec) {
15851bd0b88eSStefano Zampini     PetscInt i,m,n_neigh,*neigh,*n_shared,**shared;
15861bd0b88eSStefano Zampini 
15879566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(n+1,&mapping->info_nodec,n,&mapping->info_nodei));
15889566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetInfo(mapping,&n_neigh,&neigh,&n_shared,&shared));
1589071fcb05SBarry Smith     for (i=0;i<n;i++) { mapping->info_nodec[i] = 1;}
1590071fcb05SBarry Smith     m = n;
1591071fcb05SBarry Smith     mapping->info_nodec[n] = 0;
15921bd0b88eSStefano Zampini     for (i=1;i<n_neigh;i++) {
15931bd0b88eSStefano Zampini       PetscInt j;
15941bd0b88eSStefano Zampini 
15951bd0b88eSStefano Zampini       m += n_shared[i];
15961bd0b88eSStefano Zampini       for (j=0;j<n_shared[i];j++) mapping->info_nodec[shared[i][j]] += 1;
15971bd0b88eSStefano Zampini     }
15989566063dSJacob Faibussowitsch     if (n) PetscCall(PetscMalloc1(m,&mapping->info_nodei[0]));
15991bd0b88eSStefano Zampini     for (i=1;i<n;i++) mapping->info_nodei[i] = mapping->info_nodei[i-1] + mapping->info_nodec[i-1];
16009566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(mapping->info_nodec,n));
16011bd0b88eSStefano Zampini     for (i=0;i<n;i++) { mapping->info_nodec[i] = 1; mapping->info_nodei[i][0] = neigh[0]; }
16021bd0b88eSStefano Zampini     for (i=1;i<n_neigh;i++) {
16031bd0b88eSStefano Zampini       PetscInt j;
16041bd0b88eSStefano Zampini 
16051bd0b88eSStefano Zampini       for (j=0;j<n_shared[i];j++) {
16061bd0b88eSStefano Zampini         PetscInt k = shared[i][j];
16071bd0b88eSStefano Zampini 
16081bd0b88eSStefano Zampini         mapping->info_nodei[k][mapping->info_nodec[k]] = neigh[i];
16091bd0b88eSStefano Zampini         mapping->info_nodec[k] += 1;
16101bd0b88eSStefano Zampini       }
16111bd0b88eSStefano Zampini     }
16129566063dSJacob Faibussowitsch     for (i=0;i<n;i++) PetscCall(PetscSortRemoveDupsInt(&mapping->info_nodec[i],mapping->info_nodei[i]));
16139566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreInfo(mapping,&n_neigh,&neigh,&n_shared,&shared));
16141bd0b88eSStefano Zampini   }
16151bd0b88eSStefano Zampini   if (nnodes)  *nnodes  = n;
16161bd0b88eSStefano Zampini   if (count)   *count   = mapping->info_nodec;
16171bd0b88eSStefano Zampini   if (indices) *indices = mapping->info_nodei;
16181bd0b88eSStefano Zampini   PetscFunctionReturn(0);
16191bd0b88eSStefano Zampini }
16201bd0b88eSStefano Zampini 
16211bd0b88eSStefano Zampini /*@C
16221bd0b88eSStefano Zampini     ISLocalToGlobalMappingRestoreNodeInfo - Frees the memory allocated by ISLocalToGlobalMappingGetNodeInfo()
16231bd0b88eSStefano Zampini 
16241bd0b88eSStefano Zampini     Collective on ISLocalToGlobalMapping
16251bd0b88eSStefano Zampini 
1626f899ff85SJose E. Roman     Input Parameter:
16271bd0b88eSStefano Zampini .   mapping - the mapping from local to global indexing
16281bd0b88eSStefano Zampini 
1629d8d19677SJose E. Roman     Output Parameters:
16301bd0b88eSStefano Zampini +   nnodes - number of local nodes
16311bd0b88eSStefano Zampini .   count - number of neighboring processors per node
16321bd0b88eSStefano Zampini -   indices - indices of processes sharing the node (sorted)
16331bd0b88eSStefano Zampini 
16341bd0b88eSStefano Zampini     Level: advanced
16351bd0b88eSStefano Zampini 
16361bd0b88eSStefano Zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
16371bd0b88eSStefano Zampini           ISLocalToGlobalMappingGetInfo()
16381bd0b88eSStefano Zampini @*/
16391bd0b88eSStefano Zampini PetscErrorCode  ISLocalToGlobalMappingRestoreNodeInfo(ISLocalToGlobalMapping mapping,PetscInt *nnodes,PetscInt *count[],PetscInt **indices[])
16401bd0b88eSStefano Zampini {
16411bd0b88eSStefano Zampini   PetscFunctionBegin;
16421bd0b88eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
16431bd0b88eSStefano Zampini   if (nnodes)  *nnodes  = 0;
16441bd0b88eSStefano Zampini   if (count)   *count   = NULL;
16451bd0b88eSStefano Zampini   if (indices) *indices = NULL;
16461bd0b88eSStefano Zampini   PetscFunctionReturn(0);
16471bd0b88eSStefano Zampini }
16481bd0b88eSStefano Zampini 
16491bd0b88eSStefano Zampini /*@C
1650107e9a97SBarry Smith    ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped
165186994e45SJed Brown 
165286994e45SJed Brown    Not Collective
165386994e45SJed Brown 
16544165533cSJose E. Roman    Input Parameter:
165586994e45SJed Brown . ltog - local to global mapping
165686994e45SJed Brown 
16574165533cSJose E. Roman    Output Parameter:
1658565245c5SBarry Smith . array - array of indices, the length of this array may be obtained with ISLocalToGlobalMappingGetSize()
165986994e45SJed Brown 
166086994e45SJed Brown    Level: advanced
166186994e45SJed Brown 
166295452b02SPatrick Sanan    Notes:
166395452b02SPatrick Sanan     ISLocalToGlobalMappingGetSize() returns the length the this array
1664107e9a97SBarry Smith 
1665107e9a97SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreIndices(), ISLocalToGlobalMappingGetBlockIndices(), ISLocalToGlobalMappingRestoreBlockIndices()
166686994e45SJed Brown @*/
16677087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
166886994e45SJed Brown {
166986994e45SJed Brown   PetscFunctionBegin;
167086994e45SJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
167186994e45SJed Brown   PetscValidPointer(array,2);
167245b6f7e9SBarry Smith   if (ltog->bs == 1) {
167386994e45SJed Brown     *array = ltog->indices;
167445b6f7e9SBarry Smith   } else {
167545b6f7e9SBarry Smith     PetscInt       *jj,k,i,j,n = ltog->n, bs = ltog->bs;
167645b6f7e9SBarry Smith     const PetscInt *ii;
167745b6f7e9SBarry Smith 
16789566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(bs*n,&jj));
167945b6f7e9SBarry Smith     *array = jj;
168045b6f7e9SBarry Smith     k    = 0;
168145b6f7e9SBarry Smith     ii   = ltog->indices;
168245b6f7e9SBarry Smith     for (i=0; i<n; i++)
168345b6f7e9SBarry Smith       for (j=0; j<bs; j++)
168445b6f7e9SBarry Smith         jj[k++] = bs*ii[i] + j;
168545b6f7e9SBarry Smith   }
168686994e45SJed Brown   PetscFunctionReturn(0);
168786994e45SJed Brown }
168886994e45SJed Brown 
168986994e45SJed Brown /*@C
1690193a2b41SJulian Andrej    ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with ISLocalToGlobalMappingGetIndices()
169186994e45SJed Brown 
169286994e45SJed Brown    Not Collective
169386994e45SJed Brown 
16944165533cSJose E. Roman    Input Parameters:
169586994e45SJed Brown + ltog - local to global mapping
169686994e45SJed Brown - array - array of indices
169786994e45SJed Brown 
169886994e45SJed Brown    Level: advanced
169986994e45SJed Brown 
170086994e45SJed Brown .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices()
170186994e45SJed Brown @*/
17027087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
170386994e45SJed Brown {
170486994e45SJed Brown   PetscFunctionBegin;
170586994e45SJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
170686994e45SJed Brown   PetscValidPointer(array,2);
17072c71b3e2SJacob Faibussowitsch   PetscCheckFalse(ltog->bs == 1 && *array != ltog->indices,PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
170845b6f7e9SBarry Smith 
170945b6f7e9SBarry Smith   if (ltog->bs > 1) {
17109566063dSJacob Faibussowitsch     PetscCall(PetscFree(*(void**)array));
171145b6f7e9SBarry Smith   }
171245b6f7e9SBarry Smith   PetscFunctionReturn(0);
171345b6f7e9SBarry Smith }
171445b6f7e9SBarry Smith 
171545b6f7e9SBarry Smith /*@C
171645b6f7e9SBarry Smith    ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block
171745b6f7e9SBarry Smith 
171845b6f7e9SBarry Smith    Not Collective
171945b6f7e9SBarry Smith 
17204165533cSJose E. Roman    Input Parameter:
172145b6f7e9SBarry Smith . ltog - local to global mapping
172245b6f7e9SBarry Smith 
17234165533cSJose E. Roman    Output Parameter:
172445b6f7e9SBarry Smith . array - array of indices
172545b6f7e9SBarry Smith 
172645b6f7e9SBarry Smith    Level: advanced
172745b6f7e9SBarry Smith 
172845b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreBlockIndices()
172945b6f7e9SBarry Smith @*/
173045b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
173145b6f7e9SBarry Smith {
173245b6f7e9SBarry Smith   PetscFunctionBegin;
173345b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
173445b6f7e9SBarry Smith   PetscValidPointer(array,2);
173545b6f7e9SBarry Smith   *array = ltog->indices;
173645b6f7e9SBarry Smith   PetscFunctionReturn(0);
173745b6f7e9SBarry Smith }
173845b6f7e9SBarry Smith 
173945b6f7e9SBarry Smith /*@C
174045b6f7e9SBarry Smith    ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with ISLocalToGlobalMappingGetBlockIndices()
174145b6f7e9SBarry Smith 
174245b6f7e9SBarry Smith    Not Collective
174345b6f7e9SBarry Smith 
17444165533cSJose E. Roman    Input Parameters:
174545b6f7e9SBarry Smith + ltog - local to global mapping
174645b6f7e9SBarry Smith - array - array of indices
174745b6f7e9SBarry Smith 
174845b6f7e9SBarry Smith    Level: advanced
174945b6f7e9SBarry Smith 
175045b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices()
175145b6f7e9SBarry Smith @*/
175245b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
175345b6f7e9SBarry Smith {
175445b6f7e9SBarry Smith   PetscFunctionBegin;
175545b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
175645b6f7e9SBarry Smith   PetscValidPointer(array,2);
1757*08401ef6SPierre Jolivet   PetscCheck(*array == ltog->indices,PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
17580298fd71SBarry Smith   *array = NULL;
175986994e45SJed Brown   PetscFunctionReturn(0);
176086994e45SJed Brown }
1761f7efa3c7SJed Brown 
1762f7efa3c7SJed Brown /*@C
1763f7efa3c7SJed Brown    ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings
1764f7efa3c7SJed Brown 
1765f7efa3c7SJed Brown    Not Collective
1766f7efa3c7SJed Brown 
17674165533cSJose E. Roman    Input Parameters:
1768f7efa3c7SJed Brown + comm - communicator for the new mapping, must contain the communicator of every mapping to concatenate
1769f7efa3c7SJed Brown . n - number of mappings to concatenate
1770f7efa3c7SJed Brown - ltogs - local to global mappings
1771f7efa3c7SJed Brown 
17724165533cSJose E. Roman    Output Parameter:
1773f7efa3c7SJed Brown . ltogcat - new mapping
1774f7efa3c7SJed Brown 
17759d90f715SBarry Smith    Note: this currently always returns a mapping with block size of 1
17769d90f715SBarry Smith 
17779d90f715SBarry Smith    Developer Note: If all the input mapping have the same block size we could easily handle that as a special case
17789d90f715SBarry Smith 
1779f7efa3c7SJed Brown    Level: advanced
1780f7efa3c7SJed Brown 
1781f7efa3c7SJed Brown .seealso: ISLocalToGlobalMappingCreate()
1782f7efa3c7SJed Brown @*/
1783f7efa3c7SJed Brown PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm,PetscInt n,const ISLocalToGlobalMapping ltogs[],ISLocalToGlobalMapping *ltogcat)
1784f7efa3c7SJed Brown {
1785f7efa3c7SJed Brown   PetscInt       i,cnt,m,*idx;
1786f7efa3c7SJed Brown 
1787f7efa3c7SJed Brown   PetscFunctionBegin;
1788*08401ef6SPierre Jolivet   PetscCheck(n >= 0,comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have a non-negative number of mappings, given %" PetscInt_FMT,n);
1789f7efa3c7SJed Brown   if (n > 0) PetscValidPointer(ltogs,3);
1790f7efa3c7SJed Brown   for (i=0; i<n; i++) PetscValidHeaderSpecific(ltogs[i],IS_LTOGM_CLASSID,3);
1791f7efa3c7SJed Brown   PetscValidPointer(ltogcat,4);
1792f7efa3c7SJed Brown   for (cnt=0,i=0; i<n; i++) {
17939566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(ltogs[i],&m));
1794f7efa3c7SJed Brown     cnt += m;
1795f7efa3c7SJed Brown   }
17969566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(cnt,&idx));
1797f7efa3c7SJed Brown   for (cnt=0,i=0; i<n; i++) {
1798f7efa3c7SJed Brown     const PetscInt *subidx;
17999566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(ltogs[i],&m));
18009566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetIndices(ltogs[i],&subidx));
18019566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(&idx[cnt],subidx,m));
18029566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreIndices(ltogs[i],&subidx));
1803f7efa3c7SJed Brown     cnt += m;
1804f7efa3c7SJed Brown   }
18059566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(comm,1,cnt,idx,PETSC_OWN_POINTER,ltogcat));
1806f7efa3c7SJed Brown   PetscFunctionReturn(0);
1807f7efa3c7SJed Brown }
180804a59952SBarry Smith 
1809413f72f0SBarry Smith /*MC
1810413f72f0SBarry Smith       ISLOCALTOGLOBALMAPPINGBASIC - basic implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is
1811413f72f0SBarry Smith                                     used this is good for only small and moderate size problems.
1812413f72f0SBarry Smith 
1813413f72f0SBarry Smith    Options Database Keys:
1814a2b725a8SWilliam Gropp .   -islocaltoglobalmapping_type basic - select this method
1815413f72f0SBarry Smith 
1816413f72f0SBarry Smith    Level: beginner
1817413f72f0SBarry Smith 
1818413f72f0SBarry Smith .seealso:  ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType(), ISLOCALTOGLOBALMAPPINGHASH
1819413f72f0SBarry Smith M*/
1820413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Basic(ISLocalToGlobalMapping ltog)
1821413f72f0SBarry Smith {
1822413f72f0SBarry Smith   PetscFunctionBegin;
1823413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Basic;
1824413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Basic;
1825413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Basic;
1826413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Basic;
1827413f72f0SBarry Smith   PetscFunctionReturn(0);
1828413f72f0SBarry Smith }
1829413f72f0SBarry Smith 
1830413f72f0SBarry Smith /*MC
1831413f72f0SBarry Smith       ISLOCALTOGLOBALMAPPINGHASH - hash implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is
1832ed56e8eaSBarry Smith                                     used this is good for large memory problems.
1833413f72f0SBarry Smith 
1834413f72f0SBarry Smith    Options Database Keys:
1835a2b725a8SWilliam Gropp .   -islocaltoglobalmapping_type hash - select this method
1836413f72f0SBarry Smith 
183795452b02SPatrick Sanan    Notes:
183895452b02SPatrick Sanan     This is selected automatically for large problems if the user does not set the type.
1839ed56e8eaSBarry Smith 
1840413f72f0SBarry Smith    Level: beginner
1841413f72f0SBarry Smith 
1842413f72f0SBarry Smith .seealso:  ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType(), ISLOCALTOGLOBALMAPPINGHASH
1843413f72f0SBarry Smith M*/
1844413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Hash(ISLocalToGlobalMapping ltog)
1845413f72f0SBarry Smith {
1846413f72f0SBarry Smith   PetscFunctionBegin;
1847413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Hash;
1848413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Hash;
1849413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Hash;
1850413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Hash;
1851413f72f0SBarry Smith   PetscFunctionReturn(0);
1852413f72f0SBarry Smith }
1853413f72f0SBarry Smith 
1854413f72f0SBarry Smith /*@C
1855413f72f0SBarry Smith     ISLocalToGlobalMappingRegister -  Adds a method for applying a global to local mapping with an ISLocalToGlobalMapping
1856413f72f0SBarry Smith 
1857413f72f0SBarry Smith    Not Collective
1858413f72f0SBarry Smith 
1859413f72f0SBarry Smith    Input Parameters:
1860413f72f0SBarry Smith +  sname - name of a new method
1861413f72f0SBarry Smith -  routine_create - routine to create method context
1862413f72f0SBarry Smith 
1863413f72f0SBarry Smith    Notes:
1864ed56e8eaSBarry Smith    ISLocalToGlobalMappingRegister() may be called multiple times to add several user-defined mappings.
1865413f72f0SBarry Smith 
1866413f72f0SBarry Smith    Sample usage:
1867413f72f0SBarry Smith .vb
1868413f72f0SBarry Smith    ISLocalToGlobalMappingRegister("my_mapper",MyCreate);
1869413f72f0SBarry Smith .ve
1870413f72f0SBarry Smith 
1871ed56e8eaSBarry Smith    Then, your mapping can be chosen with the procedural interface via
1872413f72f0SBarry Smith $     ISLocalToGlobalMappingSetType(ltog,"my_mapper")
1873413f72f0SBarry Smith    or at runtime via the option
1874ed56e8eaSBarry Smith $     -islocaltoglobalmapping_type my_mapper
1875413f72f0SBarry Smith 
1876413f72f0SBarry Smith    Level: advanced
1877413f72f0SBarry Smith 
1878413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingRegisterAll(), ISLocalToGlobalMappingRegisterDestroy(), ISLOCALTOGLOBALMAPPINGBASIC, ISLOCALTOGLOBALMAPPINGHASH
1879413f72f0SBarry Smith 
1880413f72f0SBarry Smith @*/
1881413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRegister(const char sname[],PetscErrorCode (*function)(ISLocalToGlobalMapping))
1882413f72f0SBarry Smith {
1883413f72f0SBarry Smith   PetscFunctionBegin;
18849566063dSJacob Faibussowitsch   PetscCall(ISInitializePackage());
18859566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&ISLocalToGlobalMappingList,sname,function));
1886413f72f0SBarry Smith   PetscFunctionReturn(0);
1887413f72f0SBarry Smith }
1888413f72f0SBarry Smith 
1889413f72f0SBarry Smith /*@C
1890ed56e8eaSBarry Smith    ISLocalToGlobalMappingSetType - Builds ISLocalToGlobalMapping for a particular global to local mapping approach.
1891413f72f0SBarry Smith 
1892413f72f0SBarry Smith    Logically Collective on ISLocalToGlobalMapping
1893413f72f0SBarry Smith 
1894413f72f0SBarry Smith    Input Parameters:
1895413f72f0SBarry Smith +  ltog - the ISLocalToGlobalMapping object
1896413f72f0SBarry Smith -  type - a known method
1897413f72f0SBarry Smith 
1898413f72f0SBarry Smith    Options Database Key:
1899ed56e8eaSBarry Smith .  -islocaltoglobalmapping_type  <method> - Sets the method; use -help for a list
1900413f72f0SBarry Smith     of available methods (for instance, basic or hash)
1901413f72f0SBarry Smith 
1902413f72f0SBarry Smith    Notes:
1903413f72f0SBarry Smith    See "petsc/include/petscis.h" for available methods
1904413f72f0SBarry Smith 
1905413f72f0SBarry Smith   Normally, it is best to use the ISLocalToGlobalMappingSetFromOptions() command and
1906413f72f0SBarry Smith   then set the ISLocalToGlobalMapping type from the options database rather than by using
1907413f72f0SBarry Smith   this routine.
1908413f72f0SBarry Smith 
1909413f72f0SBarry Smith   Level: intermediate
1910413f72f0SBarry Smith 
1911413f72f0SBarry Smith   Developer Note: ISLocalToGlobalMappingRegister() is used to add new types to ISLocalToGlobalMappingList from which they
1912413f72f0SBarry Smith   are accessed by ISLocalToGlobalMappingSetType().
1913413f72f0SBarry Smith 
1914a0d79125SStefano Zampini .seealso: ISLocalToGlobalMappingType, ISLocalToGlobalMappingRegister(), ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingGetType()
1915413f72f0SBarry Smith @*/
1916413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType type)
1917413f72f0SBarry Smith {
1918413f72f0SBarry Smith   PetscBool      match;
19195f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(ISLocalToGlobalMapping) = NULL;
1920413f72f0SBarry Smith 
1921413f72f0SBarry Smith   PetscFunctionBegin;
1922413f72f0SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
1923a0d79125SStefano Zampini   if (type) PetscValidCharPointer(type,2);
1924413f72f0SBarry Smith 
19259566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ltog,type,&match));
1926413f72f0SBarry Smith   if (match) PetscFunctionReturn(0);
1927413f72f0SBarry Smith 
1928a0d79125SStefano Zampini   /* L2G maps defer type setup at globaltolocal calls, allow passing NULL here */
1929a0d79125SStefano Zampini   if (type) {
19309566063dSJacob Faibussowitsch     PetscCall(PetscFunctionListFind(ISLocalToGlobalMappingList,type,&r));
1931a0d79125SStefano Zampini     PetscCheck(r,PetscObjectComm((PetscObject)ltog),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested ISLocalToGlobalMapping type %s",type);
1932a0d79125SStefano Zampini   }
1933413f72f0SBarry Smith   /* Destroy the previous private LTOG context */
1934413f72f0SBarry Smith   if (ltog->ops->destroy) {
19359566063dSJacob Faibussowitsch     PetscCall((*ltog->ops->destroy)(ltog));
1936413f72f0SBarry Smith     ltog->ops->destroy = NULL;
1937413f72f0SBarry Smith   }
19389566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)ltog,type));
19399566063dSJacob Faibussowitsch   if (r) PetscCall((*r)(ltog));
1940a0d79125SStefano Zampini   PetscFunctionReturn(0);
1941a0d79125SStefano Zampini }
1942a0d79125SStefano Zampini 
1943a0d79125SStefano Zampini /*@C
1944a0d79125SStefano Zampini    ISLocalToGlobalMappingGetType - Get the type of the l2g map
1945a0d79125SStefano Zampini 
1946a0d79125SStefano Zampini    Not Collective
1947a0d79125SStefano Zampini 
1948a0d79125SStefano Zampini    Input Parameter:
1949a0d79125SStefano Zampini .  ltog - the ISLocalToGlobalMapping object
1950a0d79125SStefano Zampini 
1951a0d79125SStefano Zampini    Output Parameter:
1952a0d79125SStefano Zampini .  type - the type
1953a0d79125SStefano Zampini 
195449762cbcSSatish Balay    Level: intermediate
195549762cbcSSatish Balay 
1956a0d79125SStefano Zampini .seealso: ISLocalToGlobalMappingType, ISLocalToGlobalMappingRegister(), ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType()
1957a0d79125SStefano Zampini @*/
1958a0d79125SStefano Zampini PetscErrorCode  ISLocalToGlobalMappingGetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType *type)
1959a0d79125SStefano Zampini {
1960a0d79125SStefano Zampini   PetscFunctionBegin;
1961a0d79125SStefano Zampini   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
1962a0d79125SStefano Zampini   PetscValidPointer(type,2);
1963a0d79125SStefano Zampini   *type = ((PetscObject)ltog)->type_name;
1964413f72f0SBarry Smith   PetscFunctionReturn(0);
1965413f72f0SBarry Smith }
1966413f72f0SBarry Smith 
1967413f72f0SBarry Smith PetscBool ISLocalToGlobalMappingRegisterAllCalled = PETSC_FALSE;
1968413f72f0SBarry Smith 
1969413f72f0SBarry Smith /*@C
1970413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAll - Registers all of the local to global mapping components in the IS package.
1971413f72f0SBarry Smith 
1972413f72f0SBarry Smith   Not Collective
1973413f72f0SBarry Smith 
1974413f72f0SBarry Smith   Level: advanced
1975413f72f0SBarry Smith 
1976413f72f0SBarry Smith .seealso:  ISRegister(),  ISLocalToGlobalRegister()
1977413f72f0SBarry Smith @*/
1978413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRegisterAll(void)
1979413f72f0SBarry Smith {
1980413f72f0SBarry Smith   PetscFunctionBegin;
1981413f72f0SBarry Smith   if (ISLocalToGlobalMappingRegisterAllCalled) PetscFunctionReturn(0);
1982413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAllCalled = PETSC_TRUE;
19839566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGBASIC, ISLocalToGlobalMappingCreate_Basic));
19849566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGHASH, ISLocalToGlobalMappingCreate_Hash));
1985413f72f0SBarry Smith   PetscFunctionReturn(0);
1986413f72f0SBarry Smith }
1987