xref: /petsc/src/vec/is/utils/isltog.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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;
51*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(pointIS, &numCells));
52*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject) pointIS, ISSTRIDE, &isStride));
53*5f80ce2aSJacob Faibussowitsch   if (isStride) CHKERRQ(ISStrideGetInfo(pointIS, pStart, &step));
549305a4c7SMatthew G. Knepley   *pEnd   = *pStart + numCells;
55*5f80ce2aSJacob Faibussowitsch   if (!isStride || step != 1) CHKERRQ(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;
88*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject) pointIS, ISSTRIDE, &isStride));
89*5f80ce2aSJacob Faibussowitsch   if (isStride) CHKERRQ(ISStrideGetInfo(pointIS, pStart, &step));
90*5f80ce2aSJacob Faibussowitsch   if (!isStride || step != 1) CHKERRQ(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) {
119*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISSetType(subpointIS, ISGENERAL));
120*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGeneralSetIndices(subpointIS, pEnd-pStart, &points[pStart], PETSC_USE_POINTER));
1219305a4c7SMatthew G. Knepley   } else {
122*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISSetType(subpointIS, ISSTRIDE));
123*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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)) {
154*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISLocalToGlobalMappingSetType(mapping,ISLOCALTOGLOBALMAPPINGHASH));
155413f72f0SBarry Smith     } else {
156*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISLocalToGlobalMappingSetType(mapping,ISLOCALTOGLOBALMAPPINGBASIC));
157413f72f0SBarry Smith     }
158413f72f0SBarry Smith   }
159*5f80ce2aSJacob Faibussowitsch   if (mapping->ops->globaltolocalmappingsetup) CHKERRQ((*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;
171*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNew(&map));
172*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
180*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
190*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNew(&map));
191*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHMapICreate(&map->globalht));
192413f72f0SBarry Smith   for (i=0; i<n; i++) {
193413f72f0SBarry Smith     if (idx[i] < 0) continue;
194*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscHMapISet(map->globalht,idx[i],i));
195413f72f0SBarry Smith   }
196413f72f0SBarry Smith   mapping->data = (void*)map;
197*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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);
207*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(map->globals));
208*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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);
218*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHMapIDestroy(&map->globalht));
219*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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);
279*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)ltog),ltog->bs,ltog->n,ltog->indices,PETSC_COPY_VALUES,nltog));
280*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingGetType(ltog,&l2gtype));
281*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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);
326*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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) {
352*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mapping),&viewer));
3533050cee2SBarry Smith   }
3540700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
3555a5d4f66SBarry Smith 
356*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mapping),&rank));
357*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
35832077d6dSBarry Smith   if (iascii) {
359*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectPrintClassNamePrefixType((PetscObject)mapping,viewer));
360*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPushSynchronized(viewer));
3615a5d4f66SBarry Smith     for (i=0; i<mapping->n; i++) {
362*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %" PetscInt_FMT " %" PetscInt_FMT "\n",rank,i,mapping->indices[i]));
3636831982aSBarry Smith     }
364*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerFlush(viewer));
365*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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 
399*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject)is,&comm));
400*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(is,&n));
401*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)is,ISBLOCK,&isblock));
4026006e8d2SBarry Smith   if (!isblock) {
403*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetIndices(is,&indices));
404*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISLocalToGlobalMappingCreate(comm,1,n,indices,PETSC_COPY_VALUES,mapping));
405*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISRestoreIndices(is,&indices));
4066006e8d2SBarry Smith   } else {
407*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetBlockSize(is,&bs));
408*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISBlockGetIndices(is,&indices));
409*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISLocalToGlobalMappingCreate(comm,bs,n/bs,indices,PETSC_COPY_VALUES,mapping));
410*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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   const PetscInt *ilocal;
439a4d96a55SJed Brown   MPI_Comm       comm;
440a4d96a55SJed Brown 
441a4d96a55SJed Brown   PetscFunctionBegin;
442a4d96a55SJed Brown   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
443a4d96a55SJed Brown   PetscValidPointer(mapping,3);
444a4d96a55SJed Brown 
445*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject)sf,&comm));
446*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL));
4479a535bafSVaclav Hapla   if (start == PETSC_DECIDE) {
4489a535bafSVaclav Hapla     start = 0;
449*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Exscan(&nroots,&start,1,MPIU_INT,MPI_SUM,comm));
4502c71b3e2SJacob Faibussowitsch   } else PetscCheckFalse(start < 0,comm, PETSC_ERR_ARG_OUTOFRANGE, "start must be nonnegative or PETSC_DECIDE");
451f6e5521dSKarl Rupp   if (ilocal) {
452f6e5521dSKarl Rupp     for (i=0,maxlocal=0; i<nleaves; i++) maxlocal = PetscMax(maxlocal,ilocal[i]+1);
453f6e5521dSKarl Rupp   }
454a4d96a55SJed Brown   else maxlocal = nleaves;
455*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nroots,&globals));
456*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(maxlocal,&ltog));
457a4d96a55SJed Brown   for (i=0; i<nroots; i++) globals[i] = start + i;
458a4d96a55SJed Brown   for (i=0; i<maxlocal; i++) ltog[i] = -1;
459*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastBegin(sf,MPIU_INT,globals,ltog,MPI_REPLACE));
460*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastEnd(sf,MPIU_INT,globals,ltog,MPI_REPLACE));
461*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingCreate(comm,1,maxlocal,ltog,PETSC_OWN_POINTER,mapping));
462*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(globals));
463a4d96a55SJed Brown   PetscFunctionReturn(0);
464a4d96a55SJed Brown }
465b46b645bSBarry Smith 
46663fa5c83Sstefano_zampini /*@
46763fa5c83Sstefano_zampini     ISLocalToGlobalMappingSetBlockSize - Sets the blocksize of the mapping
46863fa5c83Sstefano_zampini 
46963fa5c83Sstefano_zampini     Not collective
47063fa5c83Sstefano_zampini 
47163fa5c83Sstefano_zampini     Input Parameters:
472a2b725a8SWilliam Gropp +   mapping - mapping data structure
473a2b725a8SWilliam Gropp -   bs - the blocksize
47463fa5c83Sstefano_zampini 
47563fa5c83Sstefano_zampini     Level: advanced
47663fa5c83Sstefano_zampini 
47763fa5c83Sstefano_zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS()
47863fa5c83Sstefano_zampini @*/
47963fa5c83Sstefano_zampini PetscErrorCode  ISLocalToGlobalMappingSetBlockSize(ISLocalToGlobalMapping mapping,PetscInt bs)
48063fa5c83Sstefano_zampini {
481a59f3c4dSstefano_zampini   PetscInt       *nid;
482a59f3c4dSstefano_zampini   const PetscInt *oid;
483a59f3c4dSstefano_zampini   PetscInt       i,cn,on,obs,nn;
48463fa5c83Sstefano_zampini 
48563fa5c83Sstefano_zampini   PetscFunctionBegin;
48663fa5c83Sstefano_zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
4872c71b3e2SJacob Faibussowitsch   PetscCheckFalse(bs < 1,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid block size %" PetscInt_FMT,bs);
48863fa5c83Sstefano_zampini   if (bs == mapping->bs) PetscFunctionReturn(0);
48963fa5c83Sstefano_zampini   on  = mapping->n;
49063fa5c83Sstefano_zampini   obs = mapping->bs;
49163fa5c83Sstefano_zampini   oid = mapping->indices;
49263fa5c83Sstefano_zampini   nn  = (on*obs)/bs;
4932c71b3e2SJacob Faibussowitsch   PetscCheckFalse((on*obs)%bs,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);
494a59f3c4dSstefano_zampini 
495*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nn,&nid));
496*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingGetIndices(mapping,&oid));
497a59f3c4dSstefano_zampini   for (i=0;i<nn;i++) {
498a59f3c4dSstefano_zampini     PetscInt j;
499a59f3c4dSstefano_zampini     for (j=0,cn=0;j<bs-1;j++) {
500a59f3c4dSstefano_zampini       if (oid[i*bs+j] < 0) { cn++; continue; }
5012c71b3e2SJacob Faibussowitsch       PetscCheckFalse(oid[i*bs+j] != oid[i*bs+j+1]-1,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Block sizes %" PetscInt_FMT " and %" PetscInt_FMT " are incompatible with the block indices: non consecutive indices %" PetscInt_FMT " %" PetscInt_FMT,bs,obs,oid[i*bs+j],oid[i*bs+j+1]);
502a59f3c4dSstefano_zampini     }
503a59f3c4dSstefano_zampini     if (oid[i*bs+j] < 0) cn++;
5048b7cb0e6Sstefano_zampini     if (cn) {
5052c71b3e2SJacob Faibussowitsch       PetscCheckFalse(cn != bs,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Block sizes %" PetscInt_FMT " and %" PetscInt_FMT " are incompatible with the block indices: invalid number of negative entries in block %" PetscInt_FMT,bs,obs,cn);
506a59f3c4dSstefano_zampini       nid[i] = -1;
5078b7cb0e6Sstefano_zampini     } else {
508a59f3c4dSstefano_zampini       nid[i] = oid[i*bs]/bs;
50963fa5c83Sstefano_zampini     }
51063fa5c83Sstefano_zampini   }
511*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingRestoreIndices(mapping,&oid));
512a59f3c4dSstefano_zampini 
51363fa5c83Sstefano_zampini   mapping->n           = nn;
51463fa5c83Sstefano_zampini   mapping->bs          = bs;
515*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(mapping->indices));
51663fa5c83Sstefano_zampini   mapping->indices     = nid;
517c9345713Sstefano_zampini   mapping->globalstart = 0;
518c9345713Sstefano_zampini   mapping->globalend   = 0;
5191bd0b88eSStefano Zampini 
5201bd0b88eSStefano Zampini   /* reset the cached information */
521*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(mapping->info_procs));
522*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(mapping->info_numprocs));
5231bd0b88eSStefano Zampini   if (mapping->info_indices) {
5241bd0b88eSStefano Zampini     PetscInt i;
5251bd0b88eSStefano Zampini 
526*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree((mapping->info_indices)[0]));
5271bd0b88eSStefano Zampini     for (i=1; i<mapping->info_nproc; i++) {
528*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(mapping->info_indices[i]));
5291bd0b88eSStefano Zampini     }
530*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(mapping->info_indices));
5311bd0b88eSStefano Zampini   }
5321bd0b88eSStefano Zampini   mapping->info_cached = PETSC_FALSE;
5331bd0b88eSStefano Zampini 
534413f72f0SBarry Smith   if (mapping->ops->destroy) {
535*5f80ce2aSJacob Faibussowitsch     CHKERRQ((*mapping->ops->destroy)(mapping));
536413f72f0SBarry Smith   }
53763fa5c83Sstefano_zampini   PetscFunctionReturn(0);
53863fa5c83Sstefano_zampini }
53963fa5c83Sstefano_zampini 
54045b6f7e9SBarry Smith /*@
54145b6f7e9SBarry Smith     ISLocalToGlobalMappingGetBlockSize - Gets the blocksize of the mapping
54245b6f7e9SBarry Smith     ordering and a global parallel ordering.
54345b6f7e9SBarry Smith 
54445b6f7e9SBarry Smith     Not Collective
54545b6f7e9SBarry Smith 
54645b6f7e9SBarry Smith     Input Parameters:
54745b6f7e9SBarry Smith .   mapping - mapping data structure
54845b6f7e9SBarry Smith 
54945b6f7e9SBarry Smith     Output Parameter:
55045b6f7e9SBarry Smith .   bs - the blocksize
55145b6f7e9SBarry Smith 
55245b6f7e9SBarry Smith     Level: advanced
55345b6f7e9SBarry Smith 
55445b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS()
55545b6f7e9SBarry Smith @*/
55645b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping mapping,PetscInt *bs)
55745b6f7e9SBarry Smith {
55845b6f7e9SBarry Smith   PetscFunctionBegin;
559cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
56045b6f7e9SBarry Smith   *bs = mapping->bs;
56145b6f7e9SBarry Smith   PetscFunctionReturn(0);
56245b6f7e9SBarry Smith }
56345b6f7e9SBarry Smith 
564ba5bb76aSSatish Balay /*@
56590f02eecSBarry Smith     ISLocalToGlobalMappingCreate - Creates a mapping between a local (0 to n)
56690f02eecSBarry Smith     ordering and a global parallel ordering.
5672362add9SBarry Smith 
56889d82c54SBarry Smith     Not Collective, but communicator may have more than one process
569b9cd556bSLois Curfman McInnes 
5702362add9SBarry Smith     Input Parameters:
57189d82c54SBarry Smith +   comm - MPI communicator
572f0413b6fSBarry Smith .   bs - the block size
57328bc9809SBarry Smith .   n - the number of local elements divided by the block size, or equivalently the number of block indices
57428bc9809SBarry 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
575d5ad8652SBarry Smith -   mode - see PetscCopyMode
5762362add9SBarry Smith 
577a997ad1aSLois Curfman McInnes     Output Parameter:
57890f02eecSBarry Smith .   mapping - new mapping data structure
5792362add9SBarry Smith 
58095452b02SPatrick Sanan     Notes:
58195452b02SPatrick Sanan     There is one integer value in indices per block and it represents the actual indices bs*idx + j, where j=0,..,bs-1
582413f72f0SBarry Smith 
5839a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
584413f72f0SBarry 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.
585413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
586413f72f0SBarry Smith 
587a997ad1aSLois Curfman McInnes     Level: advanced
588a997ad1aSLois Curfman McInnes 
589413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingSetFromOptions(), ISLOCALTOGLOBALMAPPINGBASIC, ISLOCALTOGLOBALMAPPINGHASH
590413f72f0SBarry Smith           ISLocalToGlobalMappingSetType(), ISLocalToGlobalMappingType
5912362add9SBarry Smith @*/
59260c7cefcSBarry Smith PetscErrorCode  ISLocalToGlobalMappingCreate(MPI_Comm comm,PetscInt bs,PetscInt n,const PetscInt indices[],PetscCopyMode mode,ISLocalToGlobalMapping *mapping)
5932362add9SBarry Smith {
59432dcc486SBarry Smith   PetscInt       *in;
595b46b645bSBarry Smith 
596b46b645bSBarry Smith   PetscFunctionBegin;
597064a246eSJacob Faibussowitsch   if (n) PetscValidIntPointer(indices,4);
598064a246eSJacob Faibussowitsch   PetscValidPointer(mapping,6);
599b46b645bSBarry Smith 
6000298fd71SBarry Smith   *mapping = NULL;
601*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISInitializePackage());
6022362add9SBarry Smith 
603*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHeaderCreate(*mapping,IS_LTOGM_CLASSID,"ISLocalToGlobalMapping","Local to global mapping","IS",comm,ISLocalToGlobalMappingDestroy,ISLocalToGlobalMappingView));
604d4bb536fSBarry Smith   (*mapping)->n  = n;
605f0413b6fSBarry Smith   (*mapping)->bs = bs;
606d5ad8652SBarry Smith   if (mode == PETSC_COPY_VALUES) {
607*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(n,&in));
608*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(in,indices,n));
609d5ad8652SBarry Smith     (*mapping)->indices = in;
61071910c26SVaclav Hapla     (*mapping)->dealloc_indices = PETSC_TRUE;
611*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt)));
6126389a1a1SBarry Smith   } else if (mode == PETSC_OWN_POINTER) {
6136389a1a1SBarry Smith     (*mapping)->indices = (PetscInt*)indices;
61471910c26SVaclav Hapla     (*mapping)->dealloc_indices = PETSC_TRUE;
615*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt)));
61671910c26SVaclav Hapla   } else if (mode == PETSC_USE_POINTER) {
61771910c26SVaclav Hapla     (*mapping)->indices = (PetscInt*)indices;
6186389a1a1SBarry Smith   }
61998921bdaSJacob Faibussowitsch   else SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid mode %d", mode);
6203a40ed3dSBarry Smith   PetscFunctionReturn(0);
6212362add9SBarry Smith }
6222362add9SBarry Smith 
623413f72f0SBarry Smith PetscFunctionList ISLocalToGlobalMappingList = NULL;
624413f72f0SBarry Smith 
62590f02eecSBarry Smith /*@
6267e99dc12SLawrence Mitchell    ISLocalToGlobalMappingSetFromOptions - Set mapping options from the options database.
6277e99dc12SLawrence Mitchell 
6287e99dc12SLawrence Mitchell    Not collective
6297e99dc12SLawrence Mitchell 
6307e99dc12SLawrence Mitchell    Input Parameters:
6317e99dc12SLawrence Mitchell .  mapping - mapping data structure
6327e99dc12SLawrence Mitchell 
6337e99dc12SLawrence Mitchell    Level: advanced
6347e99dc12SLawrence Mitchell 
6357e99dc12SLawrence Mitchell @*/
6367e99dc12SLawrence Mitchell PetscErrorCode ISLocalToGlobalMappingSetFromOptions(ISLocalToGlobalMapping mapping)
6377e99dc12SLawrence Mitchell {
6387e99dc12SLawrence Mitchell   PetscErrorCode             ierr;
639413f72f0SBarry Smith   char                       type[256];
640413f72f0SBarry Smith   ISLocalToGlobalMappingType defaulttype = "Not set";
6417e99dc12SLawrence Mitchell   PetscBool                  flg;
6427e99dc12SLawrence Mitchell 
6437e99dc12SLawrence Mitchell   PetscFunctionBegin;
6447e99dc12SLawrence Mitchell   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
645*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingRegisterAll());
6467e99dc12SLawrence Mitchell   ierr = PetscObjectOptionsBegin((PetscObject)mapping);CHKERRQ(ierr);
647*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsFList("-islocaltoglobalmapping_type","ISLocalToGlobalMapping method","ISLocalToGlobalMappingSetType",ISLocalToGlobalMappingList,(char*)(((PetscObject)mapping)->type_name) ? ((PetscObject)mapping)->type_name : defaulttype,type,256,&flg));
648413f72f0SBarry Smith   if (flg) {
649*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISLocalToGlobalMappingSetType(mapping,type));
650413f72f0SBarry Smith   }
6517e99dc12SLawrence Mitchell   ierr = PetscOptionsEnd();CHKERRQ(ierr);
6527e99dc12SLawrence Mitchell   PetscFunctionReturn(0);
6537e99dc12SLawrence Mitchell }
6547e99dc12SLawrence Mitchell 
6557e99dc12SLawrence Mitchell /*@
65690f02eecSBarry Smith    ISLocalToGlobalMappingDestroy - Destroys a mapping between a local (0 to n)
65790f02eecSBarry Smith    ordering and a global parallel ordering.
65890f02eecSBarry Smith 
6590f5bd95cSBarry Smith    Note Collective
660b9cd556bSLois Curfman McInnes 
66190f02eecSBarry Smith    Input Parameters:
66290f02eecSBarry Smith .  mapping - mapping data structure
66390f02eecSBarry Smith 
664a997ad1aSLois Curfman McInnes    Level: advanced
665a997ad1aSLois Curfman McInnes 
6663acfe500SLois Curfman McInnes .seealso: ISLocalToGlobalMappingCreate()
66790f02eecSBarry Smith @*/
6686bf464f9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping *mapping)
66990f02eecSBarry Smith {
6703a40ed3dSBarry Smith   PetscFunctionBegin;
6716bf464f9SBarry Smith   if (!*mapping) PetscFunctionReturn(0);
6726bf464f9SBarry Smith   PetscValidHeaderSpecific((*mapping),IS_LTOGM_CLASSID,1);
6734c8fdceaSLisandro Dalcin   if (--((PetscObject)(*mapping))->refct > 0) {*mapping = NULL;PetscFunctionReturn(0);}
67471910c26SVaclav Hapla   if ((*mapping)->dealloc_indices) {
675*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree((*mapping)->indices));
67671910c26SVaclav Hapla   }
677*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree((*mapping)->info_procs));
678*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree((*mapping)->info_numprocs));
679268a049cSStefano Zampini   if ((*mapping)->info_indices) {
680268a049cSStefano Zampini     PetscInt i;
681268a049cSStefano Zampini 
682*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(((*mapping)->info_indices)[0]));
683268a049cSStefano Zampini     for (i=1; i<(*mapping)->info_nproc; i++) {
684*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(((*mapping)->info_indices)[i]));
685268a049cSStefano Zampini     }
686*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree((*mapping)->info_indices));
687268a049cSStefano Zampini   }
6881bd0b88eSStefano Zampini   if ((*mapping)->info_nodei) {
689*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(((*mapping)->info_nodei)[0]));
6901bd0b88eSStefano Zampini   }
691*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2((*mapping)->info_nodec,(*mapping)->info_nodei));
692413f72f0SBarry Smith   if ((*mapping)->ops->destroy) {
693*5f80ce2aSJacob Faibussowitsch     CHKERRQ((*(*mapping)->ops->destroy)(*mapping));
694413f72f0SBarry Smith   }
695*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHeaderDestroy(mapping));
6964c8fdceaSLisandro Dalcin   *mapping = NULL;
6973a40ed3dSBarry Smith   PetscFunctionReturn(0);
69890f02eecSBarry Smith }
69990f02eecSBarry Smith 
70090f02eecSBarry Smith /*@
7013acfe500SLois Curfman McInnes     ISLocalToGlobalMappingApplyIS - Creates from an IS in the local numbering
7023acfe500SLois Curfman McInnes     a new index set using the global numbering defined in an ISLocalToGlobalMapping
7033acfe500SLois Curfman McInnes     context.
70490f02eecSBarry Smith 
7054cb36875SStefano Zampini     Collective on is
706b9cd556bSLois Curfman McInnes 
70790f02eecSBarry Smith     Input Parameters:
708b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
709b9cd556bSLois Curfman McInnes -   is - index set in local numbering
71090f02eecSBarry Smith 
71190f02eecSBarry Smith     Output Parameters:
71290f02eecSBarry Smith .   newis - index set in global numbering
71390f02eecSBarry Smith 
7144cb36875SStefano Zampini     Notes:
7154cb36875SStefano Zampini     The output IS will have the same communicator of the input IS.
7164cb36875SStefano Zampini 
717a997ad1aSLois Curfman McInnes     Level: advanced
718a997ad1aSLois Curfman McInnes 
71990f02eecSBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),
720d4bb536fSBarry Smith           ISLocalToGlobalMappingDestroy(), ISGlobalToLocalMappingApply()
72190f02eecSBarry Smith @*/
7227087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping,IS is,IS *newis)
72390f02eecSBarry Smith {
724e24637baSBarry Smith   PetscInt       n,*idxout;
7255d0c19d7SBarry Smith   const PetscInt *idxin;
7263a40ed3dSBarry Smith 
7273a40ed3dSBarry Smith   PetscFunctionBegin;
7280700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
7290700a824SBarry Smith   PetscValidHeaderSpecific(is,IS_CLASSID,2);
7304482741eSBarry Smith   PetscValidPointer(newis,3);
73190f02eecSBarry Smith 
732*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(is,&n));
733*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(is,&idxin));
734*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n,&idxout));
735*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingApply(mapping,n,idxin,idxout));
736*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(is,&idxin));
737*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PetscObjectComm((PetscObject)is),n,idxout,PETSC_OWN_POINTER,newis));
7383a40ed3dSBarry Smith   PetscFunctionReturn(0);
73990f02eecSBarry Smith }
74090f02eecSBarry Smith 
741b89cb25eSSatish Balay /*@
7423acfe500SLois Curfman McInnes    ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering
7433acfe500SLois Curfman McInnes    and converts them to the global numbering.
74490f02eecSBarry Smith 
745b9cd556bSLois Curfman McInnes    Not collective
746b9cd556bSLois Curfman McInnes 
747bb25748dSBarry Smith    Input Parameters:
748b9cd556bSLois Curfman McInnes +  mapping - the local to global mapping context
749bb25748dSBarry Smith .  N - number of integers
750b9cd556bSLois Curfman McInnes -  in - input indices in local numbering
751bb25748dSBarry Smith 
752bb25748dSBarry Smith    Output Parameter:
753bb25748dSBarry Smith .  out - indices in global numbering
754bb25748dSBarry Smith 
755b9cd556bSLois Curfman McInnes    Notes:
756b9cd556bSLois Curfman McInnes    The in and out array parameters may be identical.
757d4bb536fSBarry Smith 
758a997ad1aSLois Curfman McInnes    Level: advanced
759a997ad1aSLois Curfman McInnes 
76045b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
7610752156aSBarry Smith           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
762d4bb536fSBarry Smith           AOPetscToApplication(), ISGlobalToLocalMappingApply()
763bb25748dSBarry Smith 
764afcb2eb5SJed Brown @*/
765afcb2eb5SJed Brown PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
766afcb2eb5SJed Brown {
767cbc1caf0SMatthew G. Knepley   PetscInt i,bs,Nmax;
76845b6f7e9SBarry Smith 
76945b6f7e9SBarry Smith   PetscFunctionBegin;
770cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
771cbc1caf0SMatthew G. Knepley   bs   = mapping->bs;
772cbc1caf0SMatthew G. Knepley   Nmax = bs*mapping->n;
77345b6f7e9SBarry Smith   if (bs == 1) {
774cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
77545b6f7e9SBarry Smith     for (i=0; i<N; i++) {
77645b6f7e9SBarry Smith       if (in[i] < 0) {
77745b6f7e9SBarry Smith         out[i] = in[i];
77845b6f7e9SBarry Smith         continue;
77945b6f7e9SBarry Smith       }
7802c71b3e2SJacob Faibussowitsch       PetscCheckFalse(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);
78145b6f7e9SBarry Smith       out[i] = idx[in[i]];
78245b6f7e9SBarry Smith     }
78345b6f7e9SBarry Smith   } else {
784cbc1caf0SMatthew G. Knepley     const PetscInt *idx = mapping->indices;
78545b6f7e9SBarry Smith     for (i=0; i<N; i++) {
78645b6f7e9SBarry Smith       if (in[i] < 0) {
78745b6f7e9SBarry Smith         out[i] = in[i];
78845b6f7e9SBarry Smith         continue;
78945b6f7e9SBarry Smith       }
7902c71b3e2SJacob Faibussowitsch       PetscCheckFalse(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);
79145b6f7e9SBarry Smith       out[i] = idx[in[i]/bs]*bs + (in[i] % bs);
79245b6f7e9SBarry Smith     }
79345b6f7e9SBarry Smith   }
79445b6f7e9SBarry Smith   PetscFunctionReturn(0);
79545b6f7e9SBarry Smith }
79645b6f7e9SBarry Smith 
79745b6f7e9SBarry Smith /*@
7986006e8d2SBarry Smith    ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local block numbering and converts them to the global block numbering
79945b6f7e9SBarry Smith 
80045b6f7e9SBarry Smith    Not collective
80145b6f7e9SBarry Smith 
80245b6f7e9SBarry Smith    Input Parameters:
80345b6f7e9SBarry Smith +  mapping - the local to global mapping context
80445b6f7e9SBarry Smith .  N - number of integers
8056006e8d2SBarry Smith -  in - input indices in local block numbering
80645b6f7e9SBarry Smith 
80745b6f7e9SBarry Smith    Output Parameter:
8086006e8d2SBarry Smith .  out - indices in global block numbering
80945b6f7e9SBarry Smith 
81045b6f7e9SBarry Smith    Notes:
81145b6f7e9SBarry Smith    The in and out array parameters may be identical.
81245b6f7e9SBarry Smith 
8136006e8d2SBarry Smith    Example:
8146006e8d2SBarry 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
8156006e8d2SBarry Smith      (the first block) would produce 0 and the mapping applied to 1 (the second block) would produce 3.
8166006e8d2SBarry Smith 
81745b6f7e9SBarry Smith    Level: advanced
81845b6f7e9SBarry Smith 
81945b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
82045b6f7e9SBarry Smith           ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
82145b6f7e9SBarry Smith           AOPetscToApplication(), ISGlobalToLocalMappingApply()
82245b6f7e9SBarry Smith 
82345b6f7e9SBarry Smith @*/
82445b6f7e9SBarry Smith PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
82545b6f7e9SBarry Smith {
8268a1f772fSStefano Zampini   PetscInt       i,Nmax;
8278a1f772fSStefano Zampini   const PetscInt *idx;
828d4bb536fSBarry Smith 
829a0d79125SStefano Zampini   PetscFunctionBegin;
830a0d79125SStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
8318a1f772fSStefano Zampini   Nmax = mapping->n;
8328a1f772fSStefano Zampini   idx = mapping->indices;
833afcb2eb5SJed Brown   for (i=0; i<N; i++) {
834afcb2eb5SJed Brown     if (in[i] < 0) {
835afcb2eb5SJed Brown       out[i] = in[i];
836afcb2eb5SJed Brown       continue;
837afcb2eb5SJed Brown     }
8382c71b3e2SJacob Faibussowitsch     PetscCheckFalse(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);
839afcb2eb5SJed Brown     out[i] = idx[in[i]];
840afcb2eb5SJed Brown   }
841afcb2eb5SJed Brown   PetscFunctionReturn(0);
842afcb2eb5SJed Brown }
843d4bb536fSBarry Smith 
8447e99dc12SLawrence Mitchell /*@
845a997ad1aSLois Curfman McInnes     ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers
846a997ad1aSLois Curfman McInnes     specified with a global numbering.
847d4bb536fSBarry Smith 
848b9cd556bSLois Curfman McInnes     Not collective
849b9cd556bSLois Curfman McInnes 
850d4bb536fSBarry Smith     Input Parameters:
851b9cd556bSLois Curfman McInnes +   mapping - mapping between local and global numbering
8520040bde1SJunchao Zhang .   type - IS_GTOLM_MASK - maps global indices with no local value to -1 in the output list (i.e., mask them)
853d4bb536fSBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
854d4bb536fSBarry Smith .   n - number of global indices to map
855b9cd556bSLois Curfman McInnes -   idx - global indices to map
856d4bb536fSBarry Smith 
857d4bb536fSBarry Smith     Output Parameters:
858b9cd556bSLois Curfman McInnes +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
859b9cd556bSLois Curfman McInnes -   idxout - local index of each global index, one must pass in an array long enough
860e182c471SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApply() with
8610298fd71SBarry Smith              idxout == NULL to determine the required length (returned in nout)
862e182c471SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApply()
863e182c471SBarry Smith              a second time to set the values.
864d4bb536fSBarry Smith 
865b9cd556bSLois Curfman McInnes     Notes:
8660298fd71SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
867d4bb536fSBarry Smith 
8689a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
869413f72f0SBarry 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.
870413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
8710f5bd95cSBarry Smith 
872a997ad1aSLois Curfman McInnes     Level: advanced
873a997ad1aSLois Curfman McInnes 
87432fd6b96SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
87532fd6b96SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
87632fd6b96SBarry Smith 
8779d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),
878413f72f0SBarry Smith           ISLocalToGlobalMappingDestroy()
879d4bb536fSBarry Smith @*/
880413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
881d4bb536fSBarry Smith {
8829d90f715SBarry Smith   PetscFunctionBegin;
8839d90f715SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
884413f72f0SBarry Smith   if (!mapping->data) {
885*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGlobalToLocalMappingSetUp(mapping));
8869d90f715SBarry Smith   }
887*5f80ce2aSJacob Faibussowitsch   CHKERRQ((*mapping->ops->globaltolocalmappingapply)(mapping,type,n,idx,nout,idxout));
8889d90f715SBarry Smith   PetscFunctionReturn(0);
8899d90f715SBarry Smith }
8909d90f715SBarry Smith 
891d4fe737eSStefano Zampini /*@
892d4fe737eSStefano Zampini     ISGlobalToLocalMappingApplyIS - Creates from an IS in the global numbering
893d4fe737eSStefano Zampini     a new index set using the local numbering defined in an ISLocalToGlobalMapping
894d4fe737eSStefano Zampini     context.
895d4fe737eSStefano Zampini 
896d4fe737eSStefano Zampini     Not collective
897d4fe737eSStefano Zampini 
898d4fe737eSStefano Zampini     Input Parameters:
899d4fe737eSStefano Zampini +   mapping - mapping between local and global numbering
9000040bde1SJunchao Zhang .   type - IS_GTOLM_MASK - maps global indices with no local value to -1 in the output list (i.e., mask them)
9012785b321SStefano Zampini            IS_GTOLM_DROP - drops the indices with no local value from the output list
902d4fe737eSStefano Zampini -   is - index set in global numbering
903d4fe737eSStefano Zampini 
904d4fe737eSStefano Zampini     Output Parameters:
905d4fe737eSStefano Zampini .   newis - index set in local numbering
906d4fe737eSStefano Zampini 
9074cb36875SStefano Zampini     Notes:
9084cb36875SStefano Zampini     The output IS will be sequential, as it encodes a purely local operation
9094cb36875SStefano Zampini 
910d4fe737eSStefano Zampini     Level: advanced
911d4fe737eSStefano Zampini 
912d4fe737eSStefano Zampini .seealso: ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(),
913d4fe737eSStefano Zampini           ISLocalToGlobalMappingDestroy()
914d4fe737eSStefano Zampini @*/
915413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,IS is,IS *newis)
916d4fe737eSStefano Zampini {
917d4fe737eSStefano Zampini   PetscInt       n,nout,*idxout;
918d4fe737eSStefano Zampini   const PetscInt *idxin;
919d4fe737eSStefano Zampini 
920d4fe737eSStefano Zampini   PetscFunctionBegin;
921d4fe737eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
922d4fe737eSStefano Zampini   PetscValidHeaderSpecific(is,IS_CLASSID,3);
923d4fe737eSStefano Zampini   PetscValidPointer(newis,4);
924d4fe737eSStefano Zampini 
925*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(is,&n));
926*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(is,&idxin));
927d4fe737eSStefano Zampini   if (type == IS_GTOLM_MASK) {
928*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(n,&idxout));
929d4fe737eSStefano Zampini   } else {
930*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,NULL));
931*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nout,&idxout));
932d4fe737eSStefano Zampini   }
933*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,idxout));
934*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(is,&idxin));
935*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,nout,idxout,PETSC_OWN_POINTER,newis));
936d4fe737eSStefano Zampini   PetscFunctionReturn(0);
937d4fe737eSStefano Zampini }
938d4fe737eSStefano Zampini 
9399d90f715SBarry Smith /*@
9409d90f715SBarry Smith     ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers
9419d90f715SBarry Smith     specified with a block global numbering.
9429d90f715SBarry Smith 
9439d90f715SBarry Smith     Not collective
9449d90f715SBarry Smith 
9459d90f715SBarry Smith     Input Parameters:
9469d90f715SBarry Smith +   mapping - mapping between local and global numbering
9470040bde1SJunchao Zhang .   type - IS_GTOLM_MASK - maps global indices with no local value to -1 in the output list (i.e., mask them)
9489d90f715SBarry Smith            IS_GTOLM_DROP - drops the indices with no local value from the output list
9499d90f715SBarry Smith .   n - number of global indices to map
9509d90f715SBarry Smith -   idx - global indices to map
9519d90f715SBarry Smith 
9529d90f715SBarry Smith     Output Parameters:
9539d90f715SBarry Smith +   nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n)
9549d90f715SBarry Smith -   idxout - local index of each global index, one must pass in an array long enough
9559d90f715SBarry Smith              to hold all the indices. You can call ISGlobalToLocalMappingApplyBlock() with
9569d90f715SBarry Smith              idxout == NULL to determine the required length (returned in nout)
9579d90f715SBarry Smith              and then allocate the required space and call ISGlobalToLocalMappingApplyBlock()
9589d90f715SBarry Smith              a second time to set the values.
9599d90f715SBarry Smith 
9609d90f715SBarry Smith     Notes:
9619d90f715SBarry Smith     Either nout or idxout may be NULL. idx and idxout may be identical.
9629d90f715SBarry Smith 
9639a7b7924SJed Brown     For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used;
964413f72f0SBarry 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.
965413f72f0SBarry Smith     Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used.
9669d90f715SBarry Smith 
9679d90f715SBarry Smith     Level: advanced
9689d90f715SBarry Smith 
9699d90f715SBarry Smith     Developer Note: The manual page states that idx and idxout may be identical but the calling
9709d90f715SBarry Smith        sequence declares idx as const so it cannot be the same as idxout.
9719d90f715SBarry Smith 
9729d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(),
973413f72f0SBarry Smith           ISLocalToGlobalMappingDestroy()
9749d90f715SBarry Smith @*/
975413f72f0SBarry Smith PetscErrorCode  ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,
9769d90f715SBarry Smith                                                  PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[])
9779d90f715SBarry Smith {
9783a40ed3dSBarry Smith   PetscFunctionBegin;
9790700a824SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
980413f72f0SBarry Smith   if (!mapping->data) {
981*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGlobalToLocalMappingSetUp(mapping));
982d4bb536fSBarry Smith   }
983*5f80ce2aSJacob Faibussowitsch   CHKERRQ((*mapping->ops->globaltolocalmappingapplyblock)(mapping,type,n,idx,nout,idxout));
9843a40ed3dSBarry Smith   PetscFunctionReturn(0);
985d4bb536fSBarry Smith }
98690f02eecSBarry Smith 
98789d82c54SBarry Smith /*@C
9886a818285SBarry Smith     ISLocalToGlobalMappingGetBlockInfo - Gets the neighbor information for each processor and
98989d82c54SBarry Smith      each index shared by more than one processor
99089d82c54SBarry Smith 
99189d82c54SBarry Smith     Collective on ISLocalToGlobalMapping
99289d82c54SBarry Smith 
993f899ff85SJose E. Roman     Input Parameter:
99489d82c54SBarry Smith .   mapping - the mapping from local to global indexing
99589d82c54SBarry Smith 
996d8d19677SJose E. Roman     Output Parameters:
99789d82c54SBarry Smith +   nproc - number of processors that are connected to this one
99889d82c54SBarry Smith .   proc - neighboring processors
99907b52d57SBarry Smith .   numproc - number of indices for each subdomain (processor)
10003463a7baSJed Brown -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
100189d82c54SBarry Smith 
100289d82c54SBarry Smith     Level: advanced
100389d82c54SBarry Smith 
10042cfcea29SBarry Smith     Fortran Usage:
10052cfcea29SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
10062cfcea29SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
10072cfcea29SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
10082cfcea29SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
10092cfcea29SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
10102cfcea29SBarry Smith 
101107b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
101207b52d57SBarry Smith           ISLocalToGlobalMappingRestoreInfo()
101389d82c54SBarry Smith @*/
10146a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
101589d82c54SBarry Smith {
1016268a049cSStefano Zampini   PetscFunctionBegin;
1017268a049cSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1018268a049cSStefano Zampini   if (mapping->info_cached) {
1019268a049cSStefano Zampini     *nproc    = mapping->info_nproc;
1020268a049cSStefano Zampini     *procs    = mapping->info_procs;
1021268a049cSStefano Zampini     *numprocs = mapping->info_numprocs;
1022268a049cSStefano Zampini     *indices  = mapping->info_indices;
1023268a049cSStefano Zampini   } else {
1024*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISLocalToGlobalMappingGetBlockInfo_Private(mapping,nproc,procs,numprocs,indices));
1025268a049cSStefano Zampini   }
1026268a049cSStefano Zampini   PetscFunctionReturn(0);
1027268a049cSStefano Zampini }
1028268a049cSStefano Zampini 
1029268a049cSStefano Zampini static PetscErrorCode  ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
1030268a049cSStefano Zampini {
103197f1f81fSBarry Smith   PetscMPIInt    size,rank,tag1,tag2,tag3,*len,*source,imdex;
103232dcc486SBarry Smith   PetscInt       i,n = mapping->n,Ng,ng,max = 0,*lindices = mapping->indices;
103332dcc486SBarry Smith   PetscInt       *nprocs,*owner,nsends,*sends,j,*starts,nmax,nrecvs,*recvs,proc;
103497f1f81fSBarry Smith   PetscInt       cnt,scale,*ownedsenders,*nownedsenders,rstart,nowned;
103532dcc486SBarry Smith   PetscInt       node,nownedm,nt,*sends2,nsends2,*starts2,*lens2,*dest,nrecvs2,*starts3,*recvs2,k,*bprocs,*tmp;
103632dcc486SBarry Smith   PetscInt       first_procs,first_numprocs,*first_indices;
103789d82c54SBarry Smith   MPI_Request    *recv_waits,*send_waits;
103830dcb7c9SBarry Smith   MPI_Status     recv_status,*send_status,*recv_statuses;
1039ce94432eSBarry Smith   MPI_Comm       comm;
1040ace3abfcSBarry Smith   PetscBool      debug = PETSC_FALSE;
104189d82c54SBarry Smith 
104289d82c54SBarry Smith   PetscFunctionBegin;
1043*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject)mapping,&comm));
1044*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm,&size));
1045*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm,&rank));
104624cf384cSBarry Smith   if (size == 1) {
104724cf384cSBarry Smith     *nproc         = 0;
10480298fd71SBarry Smith     *procs         = NULL;
1049*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscNew(numprocs));
10501e2105dcSBarry Smith     (*numprocs)[0] = 0;
1051*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscNew(indices));
10520298fd71SBarry Smith     (*indices)[0]  = NULL;
1053268a049cSStefano Zampini     /* save info for reuse */
1054268a049cSStefano Zampini     mapping->info_nproc = *nproc;
1055268a049cSStefano Zampini     mapping->info_procs = *procs;
1056268a049cSStefano Zampini     mapping->info_numprocs = *numprocs;
1057268a049cSStefano Zampini     mapping->info_indices = *indices;
1058268a049cSStefano Zampini     mapping->info_cached = PETSC_TRUE;
105924cf384cSBarry Smith     PetscFunctionReturn(0);
106024cf384cSBarry Smith   }
106124cf384cSBarry Smith 
1062*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(((PetscObject)mapping)->options,NULL,"-islocaltoglobalmappinggetinfo_debug",&debug,NULL));
106307b52d57SBarry Smith 
10643677ff5aSBarry Smith   /*
10656a818285SBarry Smith     Notes on ISLocalToGlobalMappingGetBlockInfo
10663677ff5aSBarry Smith 
10673677ff5aSBarry Smith     globally owned node - the nodes that have been assigned to this processor in global
10683677ff5aSBarry Smith            numbering, just for this routine.
10693677ff5aSBarry Smith 
10703677ff5aSBarry Smith     nontrivial globally owned node - node assigned to this processor that is on a subdomain
10713677ff5aSBarry Smith            boundary (i.e. is has more than one local owner)
10723677ff5aSBarry Smith 
10733677ff5aSBarry Smith     locally owned node - node that exists on this processors subdomain
10743677ff5aSBarry Smith 
10753677ff5aSBarry Smith     nontrivial locally owned node - node that is not in the interior (i.e. has more than one
10763677ff5aSBarry Smith            local subdomain
10773677ff5aSBarry Smith   */
1078*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetNewTag((PetscObject)mapping,&tag1));
1079*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetNewTag((PetscObject)mapping,&tag2));
1080*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetNewTag((PetscObject)mapping,&tag3));
108189d82c54SBarry Smith 
108289d82c54SBarry Smith   for (i=0; i<n; i++) {
108389d82c54SBarry Smith     if (lindices[i] > max) max = lindices[i];
108489d82c54SBarry Smith   }
1085*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPIU_Allreduce(&max,&Ng,1,MPIU_INT,MPI_MAX,comm));
108678058e43SBarry Smith   Ng++;
1087*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm,&size));
1088*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm,&rank));
1089bc8ff85bSBarry Smith   scale  = Ng/size + 1;
1090a2e34c3dSBarry Smith   ng     = scale; if (rank == size-1) ng = Ng - scale*(size-1); ng = PetscMax(1,ng);
1091caba0dd0SBarry Smith   rstart = scale*rank;
109289d82c54SBarry Smith 
109389d82c54SBarry Smith   /* determine ownership ranges of global indices */
1094*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(2*size,&nprocs));
1095*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArrayzero(nprocs,2*size));
109689d82c54SBarry Smith 
109789d82c54SBarry Smith   /* determine owners of each local node  */
1098*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n,&owner));
109989d82c54SBarry Smith   for (i=0; i<n; i++) {
11003677ff5aSBarry Smith     proc             = lindices[i]/scale; /* processor that globally owns this index */
110127c402fcSBarry Smith     nprocs[2*proc+1] = 1;                 /* processor globally owns at least one of ours */
11023677ff5aSBarry Smith     owner[i]         = proc;
110327c402fcSBarry Smith     nprocs[2*proc]++;                     /* count of how many that processor globally owns of ours */
110489d82c54SBarry Smith   }
110527c402fcSBarry Smith   nsends = 0; for (i=0; i<size; i++) nsends += nprocs[2*i+1];
1106*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInfo(mapping,"Number of global owners for my local data %" PetscInt_FMT "\n",nsends));
110789d82c54SBarry Smith 
110889d82c54SBarry Smith   /* inform other processors of number of messages and max length*/
1109*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMaxSum(comm,nprocs,&nmax,&nrecvs));
1110*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInfo(mapping,"Number of local owners for my global data %" PetscInt_FMT "\n",nrecvs));
111189d82c54SBarry Smith 
111289d82c54SBarry Smith   /* post receives for owned rows */
1113*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1((2*nrecvs+1)*(nmax+1),&recvs));
1114*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrecvs+1,&recv_waits));
111589d82c54SBarry Smith   for (i=0; i<nrecvs; i++) {
1116*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Irecv(recvs+2*nmax*i,2*nmax,MPIU_INT,MPI_ANY_SOURCE,tag1,comm,recv_waits+i));
111789d82c54SBarry Smith   }
111889d82c54SBarry Smith 
111989d82c54SBarry Smith   /* pack messages containing lists of local nodes to owners */
1120*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(2*n+1,&sends));
1121*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(size+1,&starts));
112289d82c54SBarry Smith   starts[0] = 0;
1123f6e5521dSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
112489d82c54SBarry Smith   for (i=0; i<n; i++) {
112589d82c54SBarry Smith     sends[starts[owner[i]]++] = lindices[i];
112630dcb7c9SBarry Smith     sends[starts[owner[i]]++] = i;
112789d82c54SBarry Smith   }
1128*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(owner));
112989d82c54SBarry Smith   starts[0] = 0;
1130f6e5521dSKarl Rupp   for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2];
113189d82c54SBarry Smith 
113289d82c54SBarry Smith   /* send the messages */
1133*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nsends+1,&send_waits));
1134*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nsends+1,&dest));
113589d82c54SBarry Smith   cnt = 0;
113689d82c54SBarry Smith   for (i=0; i<size; i++) {
113727c402fcSBarry Smith     if (nprocs[2*i]) {
1138*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Isend(sends+starts[i],2*nprocs[2*i],MPIU_INT,i,tag1,comm,send_waits+cnt));
113930dcb7c9SBarry Smith       dest[cnt] = i;
114089d82c54SBarry Smith       cnt++;
114189d82c54SBarry Smith     }
114289d82c54SBarry Smith   }
1143*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(starts));
114489d82c54SBarry Smith 
114589d82c54SBarry Smith   /* wait on receives */
1146*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrecvs+1,&source));
1147*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrecvs+1,&len));
114889d82c54SBarry Smith   cnt  = nrecvs;
1149*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(ng+1,&nownedsenders));
115089d82c54SBarry Smith   while (cnt) {
1151*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status));
115289d82c54SBarry Smith     /* unpack receives into our local space */
1153*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Get_count(&recv_status,MPIU_INT,&len[imdex]));
115489d82c54SBarry Smith     source[imdex] = recv_status.MPI_SOURCE;
115530dcb7c9SBarry Smith     len[imdex]    = len[imdex]/2;
1156caba0dd0SBarry Smith     /* count how many local owners for each of my global owned indices */
115730dcb7c9SBarry Smith     for (i=0; i<len[imdex]; i++) nownedsenders[recvs[2*imdex*nmax+2*i]-rstart]++;
115889d82c54SBarry Smith     cnt--;
115989d82c54SBarry Smith   }
1160*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(recv_waits));
116189d82c54SBarry Smith 
116230dcb7c9SBarry Smith   /* count how many globally owned indices are on an edge multiplied by how many processors own them. */
1163bc8ff85bSBarry Smith   nowned  = 0;
1164bc8ff85bSBarry Smith   nownedm = 0;
1165bc8ff85bSBarry Smith   for (i=0; i<ng; i++) {
1166bc8ff85bSBarry Smith     if (nownedsenders[i] > 1) {nownedm += nownedsenders[i]; nowned++;}
1167bc8ff85bSBarry Smith   }
1168bc8ff85bSBarry Smith 
1169bc8ff85bSBarry Smith   /* create single array to contain rank of all local owners of each globally owned index */
1170*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nownedm+1,&ownedsenders));
1171*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(ng+1,&starts));
1172bc8ff85bSBarry Smith   starts[0] = 0;
1173bc8ff85bSBarry Smith   for (i=1; i<ng; i++) {
1174bc8ff85bSBarry Smith     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
1175bc8ff85bSBarry Smith     else starts[i] = starts[i-1];
1176bc8ff85bSBarry Smith   }
1177bc8ff85bSBarry Smith 
117830dcb7c9SBarry Smith   /* for each nontrival globally owned node list all arriving processors */
1179bc8ff85bSBarry Smith   for (i=0; i<nrecvs; i++) {
1180bc8ff85bSBarry Smith     for (j=0; j<len[i]; j++) {
118130dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
1182f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) ownedsenders[starts[node]++] = source[i];
1183bc8ff85bSBarry Smith     }
1184bc8ff85bSBarry Smith   }
1185bc8ff85bSBarry Smith 
118607b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
118730dcb7c9SBarry Smith     starts[0] = 0;
118830dcb7c9SBarry Smith     for (i=1; i<ng; i++) {
118930dcb7c9SBarry Smith       if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
119030dcb7c9SBarry Smith       else starts[i] = starts[i-1];
119130dcb7c9SBarry Smith     }
119230dcb7c9SBarry Smith     for (i=0; i<ng; i++) {
119330dcb7c9SBarry Smith       if (nownedsenders[i] > 1) {
1194*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSynchronizedPrintf(comm,"[%d] global node %" PetscInt_FMT " local owner processors: ",rank,i+rstart));
119530dcb7c9SBarry Smith         for (j=0; j<nownedsenders[i]; j++) {
1196*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscSynchronizedPrintf(comm,"%" PetscInt_FMT " ",ownedsenders[starts[i]+j]));
119730dcb7c9SBarry Smith         }
1198*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSynchronizedPrintf(comm,"\n"));
119930dcb7c9SBarry Smith       }
120030dcb7c9SBarry Smith     }
1201*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSynchronizedFlush(comm,PETSC_STDOUT));
120207b52d57SBarry Smith   } /* -----------------------------------  */
120330dcb7c9SBarry Smith 
12043677ff5aSBarry Smith   /* wait on original sends */
12053a96401aSBarry Smith   if (nsends) {
1206*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nsends,&send_status));
1207*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Waitall(nsends,send_waits,send_status));
1208*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(send_status));
12093a96401aSBarry Smith   }
1210*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(send_waits));
1211*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(sends));
1212*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(nprocs));
12133677ff5aSBarry Smith 
12143677ff5aSBarry Smith   /* pack messages to send back to local owners */
121530dcb7c9SBarry Smith   starts[0] = 0;
121630dcb7c9SBarry Smith   for (i=1; i<ng; i++) {
121730dcb7c9SBarry Smith     if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1];
121830dcb7c9SBarry Smith     else starts[i] = starts[i-1];
121930dcb7c9SBarry Smith   }
122030dcb7c9SBarry Smith   nsends2 = nrecvs;
1221*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nsends2+1,&nprocs)); /* length of each message */
122230dcb7c9SBarry Smith   for (i=0; i<nrecvs; i++) {
122330dcb7c9SBarry Smith     nprocs[i] = 1;
122430dcb7c9SBarry Smith     for (j=0; j<len[i]; j++) {
122530dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
1226f6e5521dSKarl Rupp       if (nownedsenders[node] > 1) nprocs[i] += 2 + nownedsenders[node];
122730dcb7c9SBarry Smith     }
122830dcb7c9SBarry Smith   }
1229f6e5521dSKarl Rupp   nt = 0;
1230f6e5521dSKarl Rupp   for (i=0; i<nsends2; i++) nt += nprocs[i];
1231f6e5521dSKarl Rupp 
1232*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nt+1,&sends2));
1233*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nsends2+1,&starts2));
1234f6e5521dSKarl Rupp 
1235f6e5521dSKarl Rupp   starts2[0] = 0;
1236f6e5521dSKarl Rupp   for (i=1; i<nsends2; i++) starts2[i] = starts2[i-1] + nprocs[i-1];
123730dcb7c9SBarry Smith   /*
123830dcb7c9SBarry Smith      Each message is 1 + nprocs[i] long, and consists of
123930dcb7c9SBarry Smith        (0) the number of nodes being sent back
124030dcb7c9SBarry Smith        (1) the local node number,
124130dcb7c9SBarry Smith        (2) the number of processors sharing it,
124230dcb7c9SBarry Smith        (3) the processors sharing it
124330dcb7c9SBarry Smith   */
124430dcb7c9SBarry Smith   for (i=0; i<nsends2; i++) {
124530dcb7c9SBarry Smith     cnt = 1;
124630dcb7c9SBarry Smith     sends2[starts2[i]] = 0;
124730dcb7c9SBarry Smith     for (j=0; j<len[i]; j++) {
124830dcb7c9SBarry Smith       node = recvs[2*i*nmax+2*j]-rstart;
124930dcb7c9SBarry Smith       if (nownedsenders[node] > 1) {
125030dcb7c9SBarry Smith         sends2[starts2[i]]++;
125130dcb7c9SBarry Smith         sends2[starts2[i]+cnt++] = recvs[2*i*nmax+2*j+1];
125230dcb7c9SBarry Smith         sends2[starts2[i]+cnt++] = nownedsenders[node];
1253*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscArraycpy(&sends2[starts2[i]+cnt],&ownedsenders[starts[node]],nownedsenders[node]));
125430dcb7c9SBarry Smith         cnt += nownedsenders[node];
125530dcb7c9SBarry Smith       }
125630dcb7c9SBarry Smith     }
125730dcb7c9SBarry Smith   }
125830dcb7c9SBarry Smith 
125930dcb7c9SBarry Smith   /* receive the message lengths */
126030dcb7c9SBarry Smith   nrecvs2 = nsends;
1261*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrecvs2+1,&lens2));
1262*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrecvs2+1,&starts3));
1263*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrecvs2+1,&recv_waits));
126430dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
1265*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Irecv(&lens2[i],1,MPIU_INT,dest[i],tag2,comm,recv_waits+i));
126630dcb7c9SBarry Smith   }
1267d44834fbSBarry Smith 
12688a8e0b3aSBarry Smith   /* send the message lengths */
12698a8e0b3aSBarry Smith   for (i=0; i<nsends2; i++) {
1270*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Send(&nprocs[i],1,MPIU_INT,source[i],tag2,comm));
12718a8e0b3aSBarry Smith   }
12728a8e0b3aSBarry Smith 
1273d44834fbSBarry Smith   /* wait on receives of lens */
12740c468ba9SBarry Smith   if (nrecvs2) {
1275*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nrecvs2,&recv_statuses));
1276*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Waitall(nrecvs2,recv_waits,recv_statuses));
1277*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(recv_statuses));
12780c468ba9SBarry Smith   }
1279*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(recv_waits));
1280d44834fbSBarry Smith 
128130dcb7c9SBarry Smith   starts3[0] = 0;
1282d44834fbSBarry Smith   nt         = 0;
128330dcb7c9SBarry Smith   for (i=0; i<nrecvs2-1; i++) {
128430dcb7c9SBarry Smith     starts3[i+1] = starts3[i] + lens2[i];
1285d44834fbSBarry Smith     nt          += lens2[i];
128630dcb7c9SBarry Smith   }
128776466f69SStefano Zampini   if (nrecvs2) nt += lens2[nrecvs2-1];
1288d44834fbSBarry Smith 
1289*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nt+1,&recvs2));
1290*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrecvs2+1,&recv_waits));
129152b72c4aSBarry Smith   for (i=0; i<nrecvs2; i++) {
1292*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Irecv(recvs2+starts3[i],lens2[i],MPIU_INT,dest[i],tag3,comm,recv_waits+i));
129330dcb7c9SBarry Smith   }
129430dcb7c9SBarry Smith 
129530dcb7c9SBarry Smith   /* send the messages */
1296*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nsends2+1,&send_waits));
129730dcb7c9SBarry Smith   for (i=0; i<nsends2; i++) {
1298*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Isend(sends2+starts2[i],nprocs[i],MPIU_INT,source[i],tag3,comm,send_waits+i));
129930dcb7c9SBarry Smith   }
130030dcb7c9SBarry Smith 
130130dcb7c9SBarry Smith   /* wait on receives */
13020c468ba9SBarry Smith   if (nrecvs2) {
1303*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nrecvs2,&recv_statuses));
1304*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Waitall(nrecvs2,recv_waits,recv_statuses));
1305*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(recv_statuses));
13060c468ba9SBarry Smith   }
1307*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(recv_waits));
1308*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(nprocs));
130930dcb7c9SBarry Smith 
131007b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
131130dcb7c9SBarry Smith     cnt = 0;
131230dcb7c9SBarry Smith     for (i=0; i<nrecvs2; i++) {
131330dcb7c9SBarry Smith       nt = recvs2[cnt++];
131430dcb7c9SBarry Smith       for (j=0; j<nt; j++) {
1315*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSynchronizedPrintf(comm,"[%d] local node %" PetscInt_FMT " number of subdomains %" PetscInt_FMT ": ",rank,recvs2[cnt],recvs2[cnt+1]));
131630dcb7c9SBarry Smith         for (k=0; k<recvs2[cnt+1]; k++) {
1317*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscSynchronizedPrintf(comm,"%" PetscInt_FMT " ",recvs2[cnt+2+k]));
131830dcb7c9SBarry Smith         }
131930dcb7c9SBarry Smith         cnt += 2 + recvs2[cnt+1];
1320*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSynchronizedPrintf(comm,"\n"));
132130dcb7c9SBarry Smith       }
132230dcb7c9SBarry Smith     }
1323*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSynchronizedFlush(comm,PETSC_STDOUT));
132407b52d57SBarry Smith   } /* -----------------------------------  */
132530dcb7c9SBarry Smith 
132630dcb7c9SBarry Smith   /* count number subdomains for each local node */
1327*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(size,&nprocs));
132830dcb7c9SBarry Smith   cnt  = 0;
132930dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
133030dcb7c9SBarry Smith     nt = recvs2[cnt++];
133130dcb7c9SBarry Smith     for (j=0; j<nt; j++) {
1332f6e5521dSKarl Rupp       for (k=0; k<recvs2[cnt+1]; k++) nprocs[recvs2[cnt+2+k]]++;
133330dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt+1];
133430dcb7c9SBarry Smith     }
133530dcb7c9SBarry Smith   }
133630dcb7c9SBarry Smith   nt = 0; for (i=0; i<size; i++) nt += (nprocs[i] > 0);
133730dcb7c9SBarry Smith   *nproc    = nt;
1338*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nt+1,procs));
1339*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nt+1,numprocs));
1340*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nt+1,indices));
13410298fd71SBarry Smith   for (i=0;i<nt+1;i++) (*indices)[i]=NULL;
1342*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(size,&bprocs));
134330dcb7c9SBarry Smith   cnt  = 0;
134430dcb7c9SBarry Smith   for (i=0; i<size; i++) {
134530dcb7c9SBarry Smith     if (nprocs[i] > 0) {
134630dcb7c9SBarry Smith       bprocs[i]        = cnt;
134730dcb7c9SBarry Smith       (*procs)[cnt]    = i;
134830dcb7c9SBarry Smith       (*numprocs)[cnt] = nprocs[i];
1349*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(nprocs[i],&(*indices)[cnt]));
135030dcb7c9SBarry Smith       cnt++;
135130dcb7c9SBarry Smith     }
135230dcb7c9SBarry Smith   }
135330dcb7c9SBarry Smith 
135430dcb7c9SBarry Smith   /* make the list of subdomains for each nontrivial local node */
1355*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArrayzero(*numprocs,nt));
135630dcb7c9SBarry Smith   cnt  = 0;
135730dcb7c9SBarry Smith   for (i=0; i<nrecvs2; i++) {
135830dcb7c9SBarry Smith     nt = recvs2[cnt++];
135930dcb7c9SBarry Smith     for (j=0; j<nt; j++) {
1360f6e5521dSKarl Rupp       for (k=0; k<recvs2[cnt+1]; k++) (*indices)[bprocs[recvs2[cnt+2+k]]][(*numprocs)[bprocs[recvs2[cnt+2+k]]]++] = recvs2[cnt];
136130dcb7c9SBarry Smith       cnt += 2 + recvs2[cnt+1];
136230dcb7c9SBarry Smith     }
136330dcb7c9SBarry Smith   }
1364*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(bprocs));
1365*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(recvs2));
136630dcb7c9SBarry Smith 
136707b52d57SBarry Smith   /* sort the node indexing by their global numbers */
136807b52d57SBarry Smith   nt = *nproc;
136907b52d57SBarry Smith   for (i=0; i<nt; i++) {
1370*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1((*numprocs)[i],&tmp));
1371f6e5521dSKarl Rupp     for (j=0; j<(*numprocs)[i]; j++) tmp[j] = lindices[(*indices)[i][j]];
1372*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSortIntWithArray((*numprocs)[i],tmp,(*indices)[i]));
1373*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(tmp));
137407b52d57SBarry Smith   }
137507b52d57SBarry Smith 
137607b52d57SBarry Smith   if (debug) { /* -----------------------------------  */
137730dcb7c9SBarry Smith     nt = *nproc;
137830dcb7c9SBarry Smith     for (i=0; i<nt; i++) {
1379*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSynchronizedPrintf(comm,"[%d] subdomain %" PetscInt_FMT " number of indices %" PetscInt_FMT ": ",rank,(*procs)[i],(*numprocs)[i]));
138030dcb7c9SBarry Smith       for (j=0; j<(*numprocs)[i]; j++) {
1381*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSynchronizedPrintf(comm,"%" PetscInt_FMT " ",(*indices)[i][j]));
138230dcb7c9SBarry Smith       }
1383*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSynchronizedPrintf(comm,"\n"));
138430dcb7c9SBarry Smith     }
1385*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSynchronizedFlush(comm,PETSC_STDOUT));
138607b52d57SBarry Smith   } /* -----------------------------------  */
138730dcb7c9SBarry Smith 
138830dcb7c9SBarry Smith   /* wait on sends */
138930dcb7c9SBarry Smith   if (nsends2) {
1390*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nsends2,&send_status));
1391*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Waitall(nsends2,send_waits,send_status));
1392*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(send_status));
139330dcb7c9SBarry Smith   }
139430dcb7c9SBarry Smith 
1395*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(starts3));
1396*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(dest));
1397*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(send_waits));
13983677ff5aSBarry Smith 
1399*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(nownedsenders));
1400*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(ownedsenders));
1401*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(starts));
1402*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(starts2));
1403*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(lens2));
140489d82c54SBarry Smith 
1405*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(source));
1406*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(len));
1407*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(recvs));
1408*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(nprocs));
1409*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(sends2));
141024cf384cSBarry Smith 
141124cf384cSBarry Smith   /* put the information about myself as the first entry in the list */
141224cf384cSBarry Smith   first_procs    = (*procs)[0];
141324cf384cSBarry Smith   first_numprocs = (*numprocs)[0];
141424cf384cSBarry Smith   first_indices  = (*indices)[0];
141524cf384cSBarry Smith   for (i=0; i<*nproc; i++) {
141624cf384cSBarry Smith     if ((*procs)[i] == rank) {
141724cf384cSBarry Smith       (*procs)[0]    = (*procs)[i];
141824cf384cSBarry Smith       (*numprocs)[0] = (*numprocs)[i];
141924cf384cSBarry Smith       (*indices)[0]  = (*indices)[i];
142024cf384cSBarry Smith       (*procs)[i]    = first_procs;
142124cf384cSBarry Smith       (*numprocs)[i] = first_numprocs;
142224cf384cSBarry Smith       (*indices)[i]  = first_indices;
142324cf384cSBarry Smith       break;
142424cf384cSBarry Smith     }
142524cf384cSBarry Smith   }
1426268a049cSStefano Zampini 
1427268a049cSStefano Zampini   /* save info for reuse */
1428268a049cSStefano Zampini   mapping->info_nproc = *nproc;
1429268a049cSStefano Zampini   mapping->info_procs = *procs;
1430268a049cSStefano Zampini   mapping->info_numprocs = *numprocs;
1431268a049cSStefano Zampini   mapping->info_indices = *indices;
1432268a049cSStefano Zampini   mapping->info_cached = PETSC_TRUE;
143389d82c54SBarry Smith   PetscFunctionReturn(0);
143489d82c54SBarry Smith }
143589d82c54SBarry Smith 
14366a818285SBarry Smith /*@C
14376a818285SBarry Smith     ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by ISLocalToGlobalMappingGetBlockInfo()
14386a818285SBarry Smith 
14396a818285SBarry Smith     Collective on ISLocalToGlobalMapping
14406a818285SBarry Smith 
1441f899ff85SJose E. Roman     Input Parameter:
14426a818285SBarry Smith .   mapping - the mapping from local to global indexing
14436a818285SBarry Smith 
1444d8d19677SJose E. Roman     Output Parameters:
14456a818285SBarry Smith +   nproc - number of processors that are connected to this one
14466a818285SBarry Smith .   proc - neighboring processors
14476a818285SBarry Smith .   numproc - number of indices for each processor
14486a818285SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
14496a818285SBarry Smith 
14506a818285SBarry Smith     Level: advanced
14516a818285SBarry Smith 
14526a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
14536a818285SBarry Smith           ISLocalToGlobalMappingGetInfo()
14546a818285SBarry Smith @*/
14556a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
14566a818285SBarry Smith {
14576a818285SBarry Smith   PetscFunctionBegin;
1458cbc1caf0SMatthew G. Knepley   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1459268a049cSStefano Zampini   if (mapping->info_free) {
1460*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(*numprocs));
14616a818285SBarry Smith     if (*indices) {
1462268a049cSStefano Zampini       PetscInt i;
1463268a049cSStefano Zampini 
1464*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree((*indices)[0]));
14656a818285SBarry Smith       for (i=1; i<*nproc; i++) {
1466*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscFree((*indices)[i]));
14676a818285SBarry Smith       }
1468*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(*indices));
14696a818285SBarry Smith     }
1470268a049cSStefano Zampini   }
1471268a049cSStefano Zampini   *nproc    = 0;
1472268a049cSStefano Zampini   *procs    = NULL;
1473268a049cSStefano Zampini   *numprocs = NULL;
1474268a049cSStefano Zampini   *indices  = NULL;
14756a818285SBarry Smith   PetscFunctionReturn(0);
14766a818285SBarry Smith }
14776a818285SBarry Smith 
14786a818285SBarry Smith /*@C
14796a818285SBarry Smith     ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each processor and
14806a818285SBarry Smith      each index shared by more than one processor
14816a818285SBarry Smith 
14826a818285SBarry Smith     Collective on ISLocalToGlobalMapping
14836a818285SBarry Smith 
1484f899ff85SJose E. Roman     Input Parameter:
14856a818285SBarry Smith .   mapping - the mapping from local to global indexing
14866a818285SBarry Smith 
1487d8d19677SJose E. Roman     Output Parameters:
14886a818285SBarry Smith +   nproc - number of processors that are connected to this one
14896a818285SBarry Smith .   proc - neighboring processors
14906a818285SBarry Smith .   numproc - number of indices for each subdomain (processor)
14916a818285SBarry Smith -   indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering)
14926a818285SBarry Smith 
14936a818285SBarry Smith     Level: advanced
14946a818285SBarry Smith 
14951bd0b88eSStefano Zampini     Notes: The user needs to call ISLocalToGlobalMappingRestoreInfo when the data is no longer needed.
14961bd0b88eSStefano Zampini 
14976a818285SBarry Smith     Fortran Usage:
14986a818285SBarry Smith $        ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by
14996a818285SBarry Smith $        ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc],
15006a818285SBarry Smith           PetscInt indices[nproc][numprocmax],ierr)
15016a818285SBarry Smith         There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and
15026a818285SBarry Smith         indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough.
15036a818285SBarry Smith 
15046a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
15056a818285SBarry Smith           ISLocalToGlobalMappingRestoreInfo()
15066a818285SBarry Smith @*/
15076a818285SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
15086a818285SBarry Smith {
15098a1f772fSStefano Zampini   PetscInt       **bindices = NULL,*bnumprocs = NULL,bs,i,j,k;
15106a818285SBarry Smith 
15116a818285SBarry Smith   PetscFunctionBegin;
15126a818285SBarry Smith   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
15138a1f772fSStefano Zampini   bs = mapping->bs;
1514*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingGetBlockInfo(mapping,nproc,procs,&bnumprocs,&bindices));
1515268a049cSStefano Zampini   if (bs > 1) { /* we need to expand the cached info */
1516*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc1(*nproc,&*indices));
1517*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc1(*nproc,&*numprocs));
15186a818285SBarry Smith     for (i=0; i<*nproc; i++) {
1519*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(bs*bnumprocs[i],&(*indices)[i]));
1520268a049cSStefano Zampini       for (j=0; j<bnumprocs[i]; j++) {
15216a818285SBarry Smith         for (k=0; k<bs; k++) {
15226a818285SBarry Smith           (*indices)[i][j*bs+k] = bs*bindices[i][j] + k;
15236a818285SBarry Smith         }
15246a818285SBarry Smith       }
1525268a049cSStefano Zampini       (*numprocs)[i] = bnumprocs[i]*bs;
15266a818285SBarry Smith     }
1527268a049cSStefano Zampini     mapping->info_free = PETSC_TRUE;
1528268a049cSStefano Zampini   } else {
1529268a049cSStefano Zampini     *numprocs = bnumprocs;
1530268a049cSStefano Zampini     *indices  = bindices;
15316a818285SBarry Smith   }
15326a818285SBarry Smith   PetscFunctionReturn(0);
15336a818285SBarry Smith }
15346a818285SBarry Smith 
153507b52d57SBarry Smith /*@C
153607b52d57SBarry Smith     ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by ISLocalToGlobalMappingGetInfo()
153789d82c54SBarry Smith 
153807b52d57SBarry Smith     Collective on ISLocalToGlobalMapping
153907b52d57SBarry Smith 
1540f899ff85SJose E. Roman     Input Parameter:
154107b52d57SBarry Smith .   mapping - the mapping from local to global indexing
154207b52d57SBarry Smith 
1543d8d19677SJose E. Roman     Output Parameters:
154407b52d57SBarry Smith +   nproc - number of processors that are connected to this one
154507b52d57SBarry Smith .   proc - neighboring processors
154607b52d57SBarry Smith .   numproc - number of indices for each processor
154707b52d57SBarry Smith -   indices - indices of local nodes shared with neighbor (sorted by global numbering)
154807b52d57SBarry Smith 
154907b52d57SBarry Smith     Level: advanced
155007b52d57SBarry Smith 
155107b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
155207b52d57SBarry Smith           ISLocalToGlobalMappingGetInfo()
155307b52d57SBarry Smith @*/
15547087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[])
155507b52d57SBarry Smith {
155607b52d57SBarry Smith   PetscFunctionBegin;
1557*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingRestoreBlockInfo(mapping,nproc,procs,numprocs,indices));
155807b52d57SBarry Smith   PetscFunctionReturn(0);
155907b52d57SBarry Smith }
156086994e45SJed Brown 
156186994e45SJed Brown /*@C
15621bd0b88eSStefano Zampini     ISLocalToGlobalMappingGetNodeInfo - Gets the neighbor information for each node
15631bd0b88eSStefano Zampini 
15641bd0b88eSStefano Zampini     Collective on ISLocalToGlobalMapping
15651bd0b88eSStefano Zampini 
1566f899ff85SJose E. Roman     Input Parameter:
15671bd0b88eSStefano Zampini .   mapping - the mapping from local to global indexing
15681bd0b88eSStefano Zampini 
1569d8d19677SJose E. Roman     Output Parameters:
15701bd0b88eSStefano Zampini +   nnodes - number of local nodes (same ISLocalToGlobalMappingGetSize())
15711bd0b88eSStefano Zampini .   count - number of neighboring processors per node
15721bd0b88eSStefano Zampini -   indices - indices of processes sharing the node (sorted)
15731bd0b88eSStefano Zampini 
15741bd0b88eSStefano Zampini     Level: advanced
15751bd0b88eSStefano Zampini 
15761bd0b88eSStefano Zampini     Notes: The user needs to call ISLocalToGlobalMappingRestoreInfo when the data is no longer needed.
15771bd0b88eSStefano Zampini 
15781bd0b88eSStefano Zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
15791bd0b88eSStefano Zampini           ISLocalToGlobalMappingGetInfo(), ISLocalToGlobalMappingRestoreNodeInfo()
15801bd0b88eSStefano Zampini @*/
15811bd0b88eSStefano Zampini PetscErrorCode  ISLocalToGlobalMappingGetNodeInfo(ISLocalToGlobalMapping mapping,PetscInt *nnodes,PetscInt *count[],PetscInt **indices[])
15821bd0b88eSStefano Zampini {
15831bd0b88eSStefano Zampini   PetscInt       n;
15841bd0b88eSStefano Zampini 
15851bd0b88eSStefano Zampini   PetscFunctionBegin;
15861bd0b88eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
1587*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingGetSize(mapping,&n));
15881bd0b88eSStefano Zampini   if (!mapping->info_nodec) {
15891bd0b88eSStefano Zampini     PetscInt i,m,n_neigh,*neigh,*n_shared,**shared;
15901bd0b88eSStefano Zampini 
1591*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc2(n+1,&mapping->info_nodec,n,&mapping->info_nodei));
1592*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISLocalToGlobalMappingGetInfo(mapping,&n_neigh,&neigh,&n_shared,&shared));
1593071fcb05SBarry Smith     for (i=0;i<n;i++) { mapping->info_nodec[i] = 1;}
1594071fcb05SBarry Smith     m = n;
1595071fcb05SBarry Smith     mapping->info_nodec[n] = 0;
15961bd0b88eSStefano Zampini     for (i=1;i<n_neigh;i++) {
15971bd0b88eSStefano Zampini       PetscInt j;
15981bd0b88eSStefano Zampini 
15991bd0b88eSStefano Zampini       m += n_shared[i];
16001bd0b88eSStefano Zampini       for (j=0;j<n_shared[i];j++) mapping->info_nodec[shared[i][j]] += 1;
16011bd0b88eSStefano Zampini     }
1602*5f80ce2aSJacob Faibussowitsch     if (n) CHKERRQ(PetscMalloc1(m,&mapping->info_nodei[0]));
16031bd0b88eSStefano Zampini     for (i=1;i<n;i++) mapping->info_nodei[i] = mapping->info_nodei[i-1] + mapping->info_nodec[i-1];
1604*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArrayzero(mapping->info_nodec,n));
16051bd0b88eSStefano Zampini     for (i=0;i<n;i++) { mapping->info_nodec[i] = 1; mapping->info_nodei[i][0] = neigh[0]; }
16061bd0b88eSStefano Zampini     for (i=1;i<n_neigh;i++) {
16071bd0b88eSStefano Zampini       PetscInt j;
16081bd0b88eSStefano Zampini 
16091bd0b88eSStefano Zampini       for (j=0;j<n_shared[i];j++) {
16101bd0b88eSStefano Zampini         PetscInt k = shared[i][j];
16111bd0b88eSStefano Zampini 
16121bd0b88eSStefano Zampini         mapping->info_nodei[k][mapping->info_nodec[k]] = neigh[i];
16131bd0b88eSStefano Zampini         mapping->info_nodec[k] += 1;
16141bd0b88eSStefano Zampini       }
16151bd0b88eSStefano Zampini     }
1616*5f80ce2aSJacob Faibussowitsch     for (i=0;i<n;i++) CHKERRQ(PetscSortRemoveDupsInt(&mapping->info_nodec[i],mapping->info_nodei[i]));
1617*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISLocalToGlobalMappingRestoreInfo(mapping,&n_neigh,&neigh,&n_shared,&shared));
16181bd0b88eSStefano Zampini   }
16191bd0b88eSStefano Zampini   if (nnodes)  *nnodes  = n;
16201bd0b88eSStefano Zampini   if (count)   *count   = mapping->info_nodec;
16211bd0b88eSStefano Zampini   if (indices) *indices = mapping->info_nodei;
16221bd0b88eSStefano Zampini   PetscFunctionReturn(0);
16231bd0b88eSStefano Zampini }
16241bd0b88eSStefano Zampini 
16251bd0b88eSStefano Zampini /*@C
16261bd0b88eSStefano Zampini     ISLocalToGlobalMappingRestoreNodeInfo - Frees the memory allocated by ISLocalToGlobalMappingGetNodeInfo()
16271bd0b88eSStefano Zampini 
16281bd0b88eSStefano Zampini     Collective on ISLocalToGlobalMapping
16291bd0b88eSStefano Zampini 
1630f899ff85SJose E. Roman     Input Parameter:
16311bd0b88eSStefano Zampini .   mapping - the mapping from local to global indexing
16321bd0b88eSStefano Zampini 
1633d8d19677SJose E. Roman     Output Parameters:
16341bd0b88eSStefano Zampini +   nnodes - number of local nodes
16351bd0b88eSStefano Zampini .   count - number of neighboring processors per node
16361bd0b88eSStefano Zampini -   indices - indices of processes sharing the node (sorted)
16371bd0b88eSStefano Zampini 
16381bd0b88eSStefano Zampini     Level: advanced
16391bd0b88eSStefano Zampini 
16401bd0b88eSStefano Zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(),
16411bd0b88eSStefano Zampini           ISLocalToGlobalMappingGetInfo()
16421bd0b88eSStefano Zampini @*/
16431bd0b88eSStefano Zampini PetscErrorCode  ISLocalToGlobalMappingRestoreNodeInfo(ISLocalToGlobalMapping mapping,PetscInt *nnodes,PetscInt *count[],PetscInt **indices[])
16441bd0b88eSStefano Zampini {
16451bd0b88eSStefano Zampini   PetscFunctionBegin;
16461bd0b88eSStefano Zampini   PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1);
16471bd0b88eSStefano Zampini   if (nnodes)  *nnodes  = 0;
16481bd0b88eSStefano Zampini   if (count)   *count   = NULL;
16491bd0b88eSStefano Zampini   if (indices) *indices = NULL;
16501bd0b88eSStefano Zampini   PetscFunctionReturn(0);
16511bd0b88eSStefano Zampini }
16521bd0b88eSStefano Zampini 
16531bd0b88eSStefano Zampini /*@C
1654107e9a97SBarry Smith    ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped
165586994e45SJed Brown 
165686994e45SJed Brown    Not Collective
165786994e45SJed Brown 
16584165533cSJose E. Roman    Input Parameter:
165986994e45SJed Brown . ltog - local to global mapping
166086994e45SJed Brown 
16614165533cSJose E. Roman    Output Parameter:
1662565245c5SBarry Smith . array - array of indices, the length of this array may be obtained with ISLocalToGlobalMappingGetSize()
166386994e45SJed Brown 
166486994e45SJed Brown    Level: advanced
166586994e45SJed Brown 
166695452b02SPatrick Sanan    Notes:
166795452b02SPatrick Sanan     ISLocalToGlobalMappingGetSize() returns the length the this array
1668107e9a97SBarry Smith 
1669107e9a97SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreIndices(), ISLocalToGlobalMappingGetBlockIndices(), ISLocalToGlobalMappingRestoreBlockIndices()
167086994e45SJed Brown @*/
16717087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
167286994e45SJed Brown {
167386994e45SJed Brown   PetscFunctionBegin;
167486994e45SJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
167586994e45SJed Brown   PetscValidPointer(array,2);
167645b6f7e9SBarry Smith   if (ltog->bs == 1) {
167786994e45SJed Brown     *array = ltog->indices;
167845b6f7e9SBarry Smith   } else {
167945b6f7e9SBarry Smith     PetscInt       *jj,k,i,j,n = ltog->n, bs = ltog->bs;
168045b6f7e9SBarry Smith     const PetscInt *ii;
168145b6f7e9SBarry Smith 
1682*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(bs*n,&jj));
168345b6f7e9SBarry Smith     *array = jj;
168445b6f7e9SBarry Smith     k    = 0;
168545b6f7e9SBarry Smith     ii   = ltog->indices;
168645b6f7e9SBarry Smith     for (i=0; i<n; i++)
168745b6f7e9SBarry Smith       for (j=0; j<bs; j++)
168845b6f7e9SBarry Smith         jj[k++] = bs*ii[i] + j;
168945b6f7e9SBarry Smith   }
169086994e45SJed Brown   PetscFunctionReturn(0);
169186994e45SJed Brown }
169286994e45SJed Brown 
169386994e45SJed Brown /*@C
1694193a2b41SJulian Andrej    ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with ISLocalToGlobalMappingGetIndices()
169586994e45SJed Brown 
169686994e45SJed Brown    Not Collective
169786994e45SJed Brown 
16984165533cSJose E. Roman    Input Parameters:
169986994e45SJed Brown + ltog - local to global mapping
170086994e45SJed Brown - array - array of indices
170186994e45SJed Brown 
170286994e45SJed Brown    Level: advanced
170386994e45SJed Brown 
170486994e45SJed Brown .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices()
170586994e45SJed Brown @*/
17067087cfbeSBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
170786994e45SJed Brown {
170886994e45SJed Brown   PetscFunctionBegin;
170986994e45SJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
171086994e45SJed Brown   PetscValidPointer(array,2);
17112c71b3e2SJacob Faibussowitsch   PetscCheckFalse(ltog->bs == 1 && *array != ltog->indices,PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
171245b6f7e9SBarry Smith 
171345b6f7e9SBarry Smith   if (ltog->bs > 1) {
1714*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(*(void**)array));
171545b6f7e9SBarry Smith   }
171645b6f7e9SBarry Smith   PetscFunctionReturn(0);
171745b6f7e9SBarry Smith }
171845b6f7e9SBarry Smith 
171945b6f7e9SBarry Smith /*@C
172045b6f7e9SBarry Smith    ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block
172145b6f7e9SBarry Smith 
172245b6f7e9SBarry Smith    Not Collective
172345b6f7e9SBarry Smith 
17244165533cSJose E. Roman    Input Parameter:
172545b6f7e9SBarry Smith . ltog - local to global mapping
172645b6f7e9SBarry Smith 
17274165533cSJose E. Roman    Output Parameter:
172845b6f7e9SBarry Smith . array - array of indices
172945b6f7e9SBarry Smith 
173045b6f7e9SBarry Smith    Level: advanced
173145b6f7e9SBarry Smith 
173245b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreBlockIndices()
173345b6f7e9SBarry Smith @*/
173445b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
173545b6f7e9SBarry Smith {
173645b6f7e9SBarry Smith   PetscFunctionBegin;
173745b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
173845b6f7e9SBarry Smith   PetscValidPointer(array,2);
173945b6f7e9SBarry Smith   *array = ltog->indices;
174045b6f7e9SBarry Smith   PetscFunctionReturn(0);
174145b6f7e9SBarry Smith }
174245b6f7e9SBarry Smith 
174345b6f7e9SBarry Smith /*@C
174445b6f7e9SBarry Smith    ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with ISLocalToGlobalMappingGetBlockIndices()
174545b6f7e9SBarry Smith 
174645b6f7e9SBarry Smith    Not Collective
174745b6f7e9SBarry Smith 
17484165533cSJose E. Roman    Input Parameters:
174945b6f7e9SBarry Smith + ltog - local to global mapping
175045b6f7e9SBarry Smith - array - array of indices
175145b6f7e9SBarry Smith 
175245b6f7e9SBarry Smith    Level: advanced
175345b6f7e9SBarry Smith 
175445b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices()
175545b6f7e9SBarry Smith @*/
175645b6f7e9SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array)
175745b6f7e9SBarry Smith {
175845b6f7e9SBarry Smith   PetscFunctionBegin;
175945b6f7e9SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
176045b6f7e9SBarry Smith   PetscValidPointer(array,2);
17612c71b3e2SJacob Faibussowitsch   PetscCheckFalse(*array != ltog->indices,PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer");
17620298fd71SBarry Smith   *array = NULL;
176386994e45SJed Brown   PetscFunctionReturn(0);
176486994e45SJed Brown }
1765f7efa3c7SJed Brown 
1766f7efa3c7SJed Brown /*@C
1767f7efa3c7SJed Brown    ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings
1768f7efa3c7SJed Brown 
1769f7efa3c7SJed Brown    Not Collective
1770f7efa3c7SJed Brown 
17714165533cSJose E. Roman    Input Parameters:
1772f7efa3c7SJed Brown + comm - communicator for the new mapping, must contain the communicator of every mapping to concatenate
1773f7efa3c7SJed Brown . n - number of mappings to concatenate
1774f7efa3c7SJed Brown - ltogs - local to global mappings
1775f7efa3c7SJed Brown 
17764165533cSJose E. Roman    Output Parameter:
1777f7efa3c7SJed Brown . ltogcat - new mapping
1778f7efa3c7SJed Brown 
17799d90f715SBarry Smith    Note: this currently always returns a mapping with block size of 1
17809d90f715SBarry Smith 
17819d90f715SBarry Smith    Developer Note: If all the input mapping have the same block size we could easily handle that as a special case
17829d90f715SBarry Smith 
1783f7efa3c7SJed Brown    Level: advanced
1784f7efa3c7SJed Brown 
1785f7efa3c7SJed Brown .seealso: ISLocalToGlobalMappingCreate()
1786f7efa3c7SJed Brown @*/
1787f7efa3c7SJed Brown PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm,PetscInt n,const ISLocalToGlobalMapping ltogs[],ISLocalToGlobalMapping *ltogcat)
1788f7efa3c7SJed Brown {
1789f7efa3c7SJed Brown   PetscInt       i,cnt,m,*idx;
1790f7efa3c7SJed Brown 
1791f7efa3c7SJed Brown   PetscFunctionBegin;
17922c71b3e2SJacob Faibussowitsch   PetscCheckFalse(n < 0,comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have a non-negative number of mappings, given %" PetscInt_FMT,n);
1793f7efa3c7SJed Brown   if (n > 0) PetscValidPointer(ltogs,3);
1794f7efa3c7SJed Brown   for (i=0; i<n; i++) PetscValidHeaderSpecific(ltogs[i],IS_LTOGM_CLASSID,3);
1795f7efa3c7SJed Brown   PetscValidPointer(ltogcat,4);
1796f7efa3c7SJed Brown   for (cnt=0,i=0; i<n; i++) {
1797*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISLocalToGlobalMappingGetSize(ltogs[i],&m));
1798f7efa3c7SJed Brown     cnt += m;
1799f7efa3c7SJed Brown   }
1800*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(cnt,&idx));
1801f7efa3c7SJed Brown   for (cnt=0,i=0; i<n; i++) {
1802f7efa3c7SJed Brown     const PetscInt *subidx;
1803*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISLocalToGlobalMappingGetSize(ltogs[i],&m));
1804*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISLocalToGlobalMappingGetIndices(ltogs[i],&subidx));
1805*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(&idx[cnt],subidx,m));
1806*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISLocalToGlobalMappingRestoreIndices(ltogs[i],&subidx));
1807f7efa3c7SJed Brown     cnt += m;
1808f7efa3c7SJed Brown   }
1809*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingCreate(comm,1,cnt,idx,PETSC_OWN_POINTER,ltogcat));
1810f7efa3c7SJed Brown   PetscFunctionReturn(0);
1811f7efa3c7SJed Brown }
181204a59952SBarry Smith 
1813413f72f0SBarry Smith /*MC
1814413f72f0SBarry Smith       ISLOCALTOGLOBALMAPPINGBASIC - basic implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is
1815413f72f0SBarry Smith                                     used this is good for only small and moderate size problems.
1816413f72f0SBarry Smith 
1817413f72f0SBarry Smith    Options Database Keys:
1818a2b725a8SWilliam Gropp .   -islocaltoglobalmapping_type basic - select this method
1819413f72f0SBarry Smith 
1820413f72f0SBarry Smith    Level: beginner
1821413f72f0SBarry Smith 
1822413f72f0SBarry Smith .seealso:  ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType(), ISLOCALTOGLOBALMAPPINGHASH
1823413f72f0SBarry Smith M*/
1824413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Basic(ISLocalToGlobalMapping ltog)
1825413f72f0SBarry Smith {
1826413f72f0SBarry Smith   PetscFunctionBegin;
1827413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Basic;
1828413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Basic;
1829413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Basic;
1830413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Basic;
1831413f72f0SBarry Smith   PetscFunctionReturn(0);
1832413f72f0SBarry Smith }
1833413f72f0SBarry Smith 
1834413f72f0SBarry Smith /*MC
1835413f72f0SBarry Smith       ISLOCALTOGLOBALMAPPINGHASH - hash implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is
1836ed56e8eaSBarry Smith                                     used this is good for large memory problems.
1837413f72f0SBarry Smith 
1838413f72f0SBarry Smith    Options Database Keys:
1839a2b725a8SWilliam Gropp .   -islocaltoglobalmapping_type hash - select this method
1840413f72f0SBarry Smith 
184195452b02SPatrick Sanan    Notes:
184295452b02SPatrick Sanan     This is selected automatically for large problems if the user does not set the type.
1843ed56e8eaSBarry Smith 
1844413f72f0SBarry Smith    Level: beginner
1845413f72f0SBarry Smith 
1846413f72f0SBarry Smith .seealso:  ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType(), ISLOCALTOGLOBALMAPPINGHASH
1847413f72f0SBarry Smith M*/
1848413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Hash(ISLocalToGlobalMapping ltog)
1849413f72f0SBarry Smith {
1850413f72f0SBarry Smith   PetscFunctionBegin;
1851413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapply      = ISGlobalToLocalMappingApply_Hash;
1852413f72f0SBarry Smith   ltog->ops->globaltolocalmappingsetup      = ISGlobalToLocalMappingSetUp_Hash;
1853413f72f0SBarry Smith   ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Hash;
1854413f72f0SBarry Smith   ltog->ops->destroy                        = ISLocalToGlobalMappingDestroy_Hash;
1855413f72f0SBarry Smith   PetscFunctionReturn(0);
1856413f72f0SBarry Smith }
1857413f72f0SBarry Smith 
1858413f72f0SBarry Smith /*@C
1859413f72f0SBarry Smith     ISLocalToGlobalMappingRegister -  Adds a method for applying a global to local mapping with an ISLocalToGlobalMapping
1860413f72f0SBarry Smith 
1861413f72f0SBarry Smith    Not Collective
1862413f72f0SBarry Smith 
1863413f72f0SBarry Smith    Input Parameters:
1864413f72f0SBarry Smith +  sname - name of a new method
1865413f72f0SBarry Smith -  routine_create - routine to create method context
1866413f72f0SBarry Smith 
1867413f72f0SBarry Smith    Notes:
1868ed56e8eaSBarry Smith    ISLocalToGlobalMappingRegister() may be called multiple times to add several user-defined mappings.
1869413f72f0SBarry Smith 
1870413f72f0SBarry Smith    Sample usage:
1871413f72f0SBarry Smith .vb
1872413f72f0SBarry Smith    ISLocalToGlobalMappingRegister("my_mapper",MyCreate);
1873413f72f0SBarry Smith .ve
1874413f72f0SBarry Smith 
1875ed56e8eaSBarry Smith    Then, your mapping can be chosen with the procedural interface via
1876413f72f0SBarry Smith $     ISLocalToGlobalMappingSetType(ltog,"my_mapper")
1877413f72f0SBarry Smith    or at runtime via the option
1878ed56e8eaSBarry Smith $     -islocaltoglobalmapping_type my_mapper
1879413f72f0SBarry Smith 
1880413f72f0SBarry Smith    Level: advanced
1881413f72f0SBarry Smith 
1882413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingRegisterAll(), ISLocalToGlobalMappingRegisterDestroy(), ISLOCALTOGLOBALMAPPINGBASIC, ISLOCALTOGLOBALMAPPINGHASH
1883413f72f0SBarry Smith 
1884413f72f0SBarry Smith @*/
1885413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRegister(const char sname[],PetscErrorCode (*function)(ISLocalToGlobalMapping))
1886413f72f0SBarry Smith {
1887413f72f0SBarry Smith   PetscFunctionBegin;
1888*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISInitializePackage());
1889*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFunctionListAdd(&ISLocalToGlobalMappingList,sname,function));
1890413f72f0SBarry Smith   PetscFunctionReturn(0);
1891413f72f0SBarry Smith }
1892413f72f0SBarry Smith 
1893413f72f0SBarry Smith /*@C
1894ed56e8eaSBarry Smith    ISLocalToGlobalMappingSetType - Builds ISLocalToGlobalMapping for a particular global to local mapping approach.
1895413f72f0SBarry Smith 
1896413f72f0SBarry Smith    Logically Collective on ISLocalToGlobalMapping
1897413f72f0SBarry Smith 
1898413f72f0SBarry Smith    Input Parameters:
1899413f72f0SBarry Smith +  ltog - the ISLocalToGlobalMapping object
1900413f72f0SBarry Smith -  type - a known method
1901413f72f0SBarry Smith 
1902413f72f0SBarry Smith    Options Database Key:
1903ed56e8eaSBarry Smith .  -islocaltoglobalmapping_type  <method> - Sets the method; use -help for a list
1904413f72f0SBarry Smith     of available methods (for instance, basic or hash)
1905413f72f0SBarry Smith 
1906413f72f0SBarry Smith    Notes:
1907413f72f0SBarry Smith    See "petsc/include/petscis.h" for available methods
1908413f72f0SBarry Smith 
1909413f72f0SBarry Smith   Normally, it is best to use the ISLocalToGlobalMappingSetFromOptions() command and
1910413f72f0SBarry Smith   then set the ISLocalToGlobalMapping type from the options database rather than by using
1911413f72f0SBarry Smith   this routine.
1912413f72f0SBarry Smith 
1913413f72f0SBarry Smith   Level: intermediate
1914413f72f0SBarry Smith 
1915413f72f0SBarry Smith   Developer Note: ISLocalToGlobalMappingRegister() is used to add new types to ISLocalToGlobalMappingList from which they
1916413f72f0SBarry Smith   are accessed by ISLocalToGlobalMappingSetType().
1917413f72f0SBarry Smith 
1918a0d79125SStefano Zampini .seealso: ISLocalToGlobalMappingType, ISLocalToGlobalMappingRegister(), ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingGetType()
1919413f72f0SBarry Smith @*/
1920413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType type)
1921413f72f0SBarry Smith {
1922413f72f0SBarry Smith   PetscBool      match;
1923*5f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(ISLocalToGlobalMapping) = NULL;
1924413f72f0SBarry Smith 
1925413f72f0SBarry Smith   PetscFunctionBegin;
1926413f72f0SBarry Smith   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
1927a0d79125SStefano Zampini   if (type) PetscValidCharPointer(type,2);
1928413f72f0SBarry Smith 
1929*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)ltog,type,&match));
1930413f72f0SBarry Smith   if (match) PetscFunctionReturn(0);
1931413f72f0SBarry Smith 
1932a0d79125SStefano Zampini   /* L2G maps defer type setup at globaltolocal calls, allow passing NULL here */
1933a0d79125SStefano Zampini   if (type) {
1934*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFunctionListFind(ISLocalToGlobalMappingList,type,&r));
1935a0d79125SStefano Zampini     PetscCheck(r,PetscObjectComm((PetscObject)ltog),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested ISLocalToGlobalMapping type %s",type);
1936a0d79125SStefano Zampini   }
1937413f72f0SBarry Smith   /* Destroy the previous private LTOG context */
1938413f72f0SBarry Smith   if (ltog->ops->destroy) {
1939*5f80ce2aSJacob Faibussowitsch     CHKERRQ((*ltog->ops->destroy)(ltog));
1940413f72f0SBarry Smith     ltog->ops->destroy = NULL;
1941413f72f0SBarry Smith   }
1942*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectChangeTypeName((PetscObject)ltog,type));
1943*5f80ce2aSJacob Faibussowitsch   if (r) CHKERRQ((*r)(ltog));
1944a0d79125SStefano Zampini   PetscFunctionReturn(0);
1945a0d79125SStefano Zampini }
1946a0d79125SStefano Zampini 
1947a0d79125SStefano Zampini /*@C
1948a0d79125SStefano Zampini    ISLocalToGlobalMappingGetType - Get the type of the l2g map
1949a0d79125SStefano Zampini 
1950a0d79125SStefano Zampini    Not Collective
1951a0d79125SStefano Zampini 
1952a0d79125SStefano Zampini    Input Parameter:
1953a0d79125SStefano Zampini .  ltog - the ISLocalToGlobalMapping object
1954a0d79125SStefano Zampini 
1955a0d79125SStefano Zampini    Output Parameter:
1956a0d79125SStefano Zampini .  type - the type
1957a0d79125SStefano Zampini 
1958a0d79125SStefano Zampini .seealso: ISLocalToGlobalMappingType, ISLocalToGlobalMappingRegister(), ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType()
1959a0d79125SStefano Zampini @*/
1960a0d79125SStefano Zampini PetscErrorCode  ISLocalToGlobalMappingGetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType *type)
1961a0d79125SStefano Zampini {
1962a0d79125SStefano Zampini   PetscFunctionBegin;
1963a0d79125SStefano Zampini   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1);
1964a0d79125SStefano Zampini   PetscValidPointer(type,2);
1965a0d79125SStefano Zampini   *type = ((PetscObject)ltog)->type_name;
1966413f72f0SBarry Smith   PetscFunctionReturn(0);
1967413f72f0SBarry Smith }
1968413f72f0SBarry Smith 
1969413f72f0SBarry Smith PetscBool ISLocalToGlobalMappingRegisterAllCalled = PETSC_FALSE;
1970413f72f0SBarry Smith 
1971413f72f0SBarry Smith /*@C
1972413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAll - Registers all of the local to global mapping components in the IS package.
1973413f72f0SBarry Smith 
1974413f72f0SBarry Smith   Not Collective
1975413f72f0SBarry Smith 
1976413f72f0SBarry Smith   Level: advanced
1977413f72f0SBarry Smith 
1978413f72f0SBarry Smith .seealso:  ISRegister(),  ISLocalToGlobalRegister()
1979413f72f0SBarry Smith @*/
1980413f72f0SBarry Smith PetscErrorCode  ISLocalToGlobalMappingRegisterAll(void)
1981413f72f0SBarry Smith {
1982413f72f0SBarry Smith   PetscFunctionBegin;
1983413f72f0SBarry Smith   if (ISLocalToGlobalMappingRegisterAllCalled) PetscFunctionReturn(0);
1984413f72f0SBarry Smith   ISLocalToGlobalMappingRegisterAllCalled = PETSC_TRUE;
1985*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGBASIC, ISLocalToGlobalMappingCreate_Basic));
1986*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGHASH, ISLocalToGlobalMappingCreate_Hash));
1987413f72f0SBarry Smith   PetscFunctionReturn(0);
1988413f72f0SBarry Smith }
1989