12362add9SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/isimpl.h> /*I "petscis.h" I*/ 3e8f14785SLisandro Dalcin #include <petsc/private/hashmapi.h> 40c312b8eSJed Brown #include <petscsf.h> 5665c2dedSJed Brown #include <petscviewer.h> 62362add9SBarry Smith 77087cfbeSBarry Smith PetscClassId IS_LTOGM_CLASSID; 8268a049cSStefano Zampini static PetscErrorCode ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping, PetscInt *, PetscInt **, PetscInt **, PetscInt ***); 98e58c17dSMatthew Knepley 10413f72f0SBarry Smith typedef struct { 11413f72f0SBarry Smith PetscInt *globals; 12413f72f0SBarry Smith } ISLocalToGlobalMapping_Basic; 13413f72f0SBarry Smith 14413f72f0SBarry Smith typedef struct { 15e8f14785SLisandro Dalcin PetscHMapI globalht; 16413f72f0SBarry Smith } ISLocalToGlobalMapping_Hash; 17413f72f0SBarry Smith 186528b96dSMatthew G. Knepley /*@C 196528b96dSMatthew G. Knepley ISGetPointRange - Returns a description of the points in an IS suitable for traversal 20413f72f0SBarry Smith 216528b96dSMatthew G. Knepley Not collective 226528b96dSMatthew G. Knepley 236528b96dSMatthew G. Knepley Input Parameter: 246528b96dSMatthew G. Knepley . pointIS - The IS object 256528b96dSMatthew G. Knepley 266528b96dSMatthew G. Knepley Output Parameters: 276528b96dSMatthew G. Knepley + pStart - The first index, see notes 286528b96dSMatthew G. Knepley . pEnd - One past the last index, see notes 296528b96dSMatthew G. Knepley - points - The indices, see notes 306528b96dSMatthew G. Knepley 316528b96dSMatthew G. Knepley Notes: 326528b96dSMatthew G. Knepley If the IS contains contiguous indices in an ISSTRIDE, then the indices are contained in [pStart, pEnd) and points = NULL. Otherwise, pStart = 0, pEnd = numIndices, and points is an array of the indices. This supports the following pattern 336528b96dSMatthew G. Knepley $ ISGetPointRange(is, &pStart, &pEnd, &points); 346528b96dSMatthew G. Knepley $ for (p = pStart; p < pEnd; ++p) { 356528b96dSMatthew G. Knepley $ const PetscInt point = points ? points[p] : p; 366528b96dSMatthew G. Knepley $ } 376528b96dSMatthew G. Knepley $ ISRestorePointRange(is, &pstart, &pEnd, &points); 386528b96dSMatthew G. Knepley 396528b96dSMatthew G. Knepley Level: intermediate 406528b96dSMatthew G. Knepley 41db781477SPatrick Sanan .seealso: `ISRestorePointRange()`, `ISGetPointSubrange()`, `ISGetIndices()`, `ISCreateStride()` 426528b96dSMatthew G. Knepley @*/ 43*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGetPointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points) 44*d71ae5a4SJacob Faibussowitsch { 459305a4c7SMatthew G. Knepley PetscInt numCells, step = 1; 469305a4c7SMatthew G. Knepley PetscBool isStride; 479305a4c7SMatthew G. Knepley 489305a4c7SMatthew G. Knepley PetscFunctionBeginHot; 499305a4c7SMatthew G. Knepley *pStart = 0; 509305a4c7SMatthew G. Knepley *points = NULL; 519566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(pointIS, &numCells)); 529566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pointIS, ISSTRIDE, &isStride)); 539566063dSJacob Faibussowitsch if (isStride) PetscCall(ISStrideGetInfo(pointIS, pStart, &step)); 549305a4c7SMatthew G. Knepley *pEnd = *pStart + numCells; 559566063dSJacob Faibussowitsch if (!isStride || step != 1) PetscCall(ISGetIndices(pointIS, points)); 569305a4c7SMatthew G. Knepley PetscFunctionReturn(0); 579305a4c7SMatthew G. Knepley } 589305a4c7SMatthew G. Knepley 596528b96dSMatthew G. Knepley /*@C 606528b96dSMatthew G. Knepley ISRestorePointRange - Destroys the traversal description 616528b96dSMatthew G. Knepley 626528b96dSMatthew G. Knepley Not collective 636528b96dSMatthew G. Knepley 646528b96dSMatthew G. Knepley Input Parameters: 656528b96dSMatthew G. Knepley + pointIS - The IS object 666528b96dSMatthew G. Knepley . pStart - The first index, from ISGetPointRange() 676528b96dSMatthew G. Knepley . pEnd - One past the last index, from ISGetPointRange() 686528b96dSMatthew G. Knepley - points - The indices, from ISGetPointRange() 696528b96dSMatthew G. Knepley 706528b96dSMatthew G. Knepley Notes: 716528b96dSMatthew G. Knepley If the IS contains contiguous indices in an ISSTRIDE, then the indices are contained in [pStart, pEnd) and points = NULL. Otherwise, pStart = 0, pEnd = numIndices, and points is an array of the indices. This supports the following pattern 726528b96dSMatthew G. Knepley $ ISGetPointRange(is, &pStart, &pEnd, &points); 736528b96dSMatthew G. Knepley $ for (p = pStart; p < pEnd; ++p) { 746528b96dSMatthew G. Knepley $ const PetscInt point = points ? points[p] : p; 756528b96dSMatthew G. Knepley $ } 766528b96dSMatthew G. Knepley $ ISRestorePointRange(is, &pstart, &pEnd, &points); 776528b96dSMatthew G. Knepley 786528b96dSMatthew G. Knepley Level: intermediate 796528b96dSMatthew G. Knepley 80db781477SPatrick Sanan .seealso: `ISGetPointRange()`, `ISGetPointSubrange()`, `ISGetIndices()`, `ISCreateStride()` 816528b96dSMatthew G. Knepley @*/ 82*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISRestorePointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points) 83*d71ae5a4SJacob Faibussowitsch { 849305a4c7SMatthew G. Knepley PetscInt step = 1; 859305a4c7SMatthew G. Knepley PetscBool isStride; 869305a4c7SMatthew G. Knepley 879305a4c7SMatthew G. Knepley PetscFunctionBeginHot; 889566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pointIS, ISSTRIDE, &isStride)); 899566063dSJacob Faibussowitsch if (isStride) PetscCall(ISStrideGetInfo(pointIS, pStart, &step)); 909566063dSJacob Faibussowitsch if (!isStride || step != 1) PetscCall(ISGetIndices(pointIS, points)); 919305a4c7SMatthew G. Knepley PetscFunctionReturn(0); 929305a4c7SMatthew G. Knepley } 939305a4c7SMatthew G. Knepley 946528b96dSMatthew G. Knepley /*@C 956528b96dSMatthew G. Knepley ISGetPointSubrange - Configures the input IS to be a subrange for the traversal information given 966528b96dSMatthew G. Knepley 976528b96dSMatthew G. Knepley Not collective 986528b96dSMatthew G. Knepley 996528b96dSMatthew G. Knepley Input Parameters: 1006528b96dSMatthew G. Knepley + subpointIS - The IS object to be configured 1016528b96dSMatthew G. Knepley . pStar t - The first index of the subrange 1026528b96dSMatthew G. Knepley . pEnd - One past the last index for the subrange 1036528b96dSMatthew G. Knepley - points - The indices for the entire range, from ISGetPointRange() 1046528b96dSMatthew G. Knepley 1056528b96dSMatthew G. Knepley Output Parameters: 1066528b96dSMatthew G. Knepley . subpointIS - The IS object now configured to be a subrange 1076528b96dSMatthew G. Knepley 1086528b96dSMatthew G. Knepley Notes: 1096528b96dSMatthew G. Knepley The input IS will now respond properly to calls to ISGetPointRange() and return the subrange. 1106528b96dSMatthew G. Knepley 1116528b96dSMatthew G. Knepley Level: intermediate 1126528b96dSMatthew G. Knepley 113db781477SPatrick Sanan .seealso: `ISGetPointRange()`, `ISRestorePointRange()`, `ISGetIndices()`, `ISCreateStride()` 1146528b96dSMatthew G. Knepley @*/ 115*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGetPointSubrange(IS subpointIS, PetscInt pStart, PetscInt pEnd, const PetscInt *points) 116*d71ae5a4SJacob Faibussowitsch { 1179305a4c7SMatthew G. Knepley PetscFunctionBeginHot; 1189305a4c7SMatthew G. Knepley if (points) { 1199566063dSJacob Faibussowitsch PetscCall(ISSetType(subpointIS, ISGENERAL)); 1209566063dSJacob Faibussowitsch PetscCall(ISGeneralSetIndices(subpointIS, pEnd - pStart, &points[pStart], PETSC_USE_POINTER)); 1219305a4c7SMatthew G. Knepley } else { 1229566063dSJacob Faibussowitsch PetscCall(ISSetType(subpointIS, ISSTRIDE)); 1239566063dSJacob Faibussowitsch PetscCall(ISStrideSetStride(subpointIS, pEnd - pStart, pStart, 1)); 1249305a4c7SMatthew G. Knepley } 1259305a4c7SMatthew G. Knepley PetscFunctionReturn(0); 1269305a4c7SMatthew G. Knepley } 1279305a4c7SMatthew G. Knepley 128413f72f0SBarry Smith /* -----------------------------------------------------------------------------------------*/ 129413f72f0SBarry Smith 130413f72f0SBarry Smith /* 131413f72f0SBarry Smith Creates the global mapping information in the ISLocalToGlobalMapping structure 132413f72f0SBarry Smith 133413f72f0SBarry Smith If the user has not selected how to handle the global to local mapping then use HASH for "large" problems 134413f72f0SBarry Smith */ 135*d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISGlobalToLocalMappingSetUp(ISLocalToGlobalMapping mapping) 136*d71ae5a4SJacob Faibussowitsch { 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 } 1499371c9d4SSatish Balay if (start > end) { 1509371c9d4SSatish Balay start = 0; 1519371c9d4SSatish Balay end = -1; 1529371c9d4SSatish Balay } 153413f72f0SBarry Smith mapping->globalstart = start; 154413f72f0SBarry Smith mapping->globalend = end; 155413f72f0SBarry Smith if (!((PetscObject)mapping)->type_name) { 156413f72f0SBarry Smith if ((end - start) > PetscMax(4 * n, 1000000)) { 1579566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingSetType(mapping, ISLOCALTOGLOBALMAPPINGHASH)); 158413f72f0SBarry Smith } else { 1599566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingSetType(mapping, ISLOCALTOGLOBALMAPPINGBASIC)); 160413f72f0SBarry Smith } 161413f72f0SBarry Smith } 162dbbe0bcdSBarry Smith PetscTryTypeMethod(mapping, globaltolocalmappingsetup); 163413f72f0SBarry Smith PetscFunctionReturn(0); 164413f72f0SBarry Smith } 165413f72f0SBarry Smith 166*d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISGlobalToLocalMappingSetUp_Basic(ISLocalToGlobalMapping mapping) 167*d71ae5a4SJacob Faibussowitsch { 168413f72f0SBarry Smith PetscInt i, *idx = mapping->indices, n = mapping->n, end, start, *globals; 169413f72f0SBarry Smith ISLocalToGlobalMapping_Basic *map; 170413f72f0SBarry Smith 171413f72f0SBarry Smith PetscFunctionBegin; 172413f72f0SBarry Smith start = mapping->globalstart; 173413f72f0SBarry Smith end = mapping->globalend; 1749566063dSJacob Faibussowitsch PetscCall(PetscNew(&map)); 1759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(end - start + 2, &globals)); 176413f72f0SBarry Smith map->globals = globals; 177413f72f0SBarry Smith for (i = 0; i < end - start + 1; i++) globals[i] = -1; 178413f72f0SBarry Smith for (i = 0; i < n; i++) { 179413f72f0SBarry Smith if (idx[i] < 0) continue; 180413f72f0SBarry Smith globals[idx[i] - start] = i; 181413f72f0SBarry Smith } 182413f72f0SBarry Smith mapping->data = (void *)map; 183413f72f0SBarry Smith PetscFunctionReturn(0); 184413f72f0SBarry Smith } 185413f72f0SBarry Smith 186*d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISGlobalToLocalMappingSetUp_Hash(ISLocalToGlobalMapping mapping) 187*d71ae5a4SJacob Faibussowitsch { 188413f72f0SBarry Smith PetscInt i, *idx = mapping->indices, n = mapping->n; 189413f72f0SBarry Smith ISLocalToGlobalMapping_Hash *map; 190413f72f0SBarry Smith 191413f72f0SBarry Smith PetscFunctionBegin; 1929566063dSJacob Faibussowitsch PetscCall(PetscNew(&map)); 1939566063dSJacob Faibussowitsch PetscCall(PetscHMapICreate(&map->globalht)); 194413f72f0SBarry Smith for (i = 0; i < n; i++) { 195413f72f0SBarry Smith if (idx[i] < 0) continue; 1969566063dSJacob Faibussowitsch PetscCall(PetscHMapISet(map->globalht, idx[i], i)); 197413f72f0SBarry Smith } 198413f72f0SBarry Smith mapping->data = (void *)map; 199413f72f0SBarry Smith PetscFunctionReturn(0); 200413f72f0SBarry Smith } 201413f72f0SBarry Smith 202*d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISLocalToGlobalMappingDestroy_Basic(ISLocalToGlobalMapping mapping) 203*d71ae5a4SJacob Faibussowitsch { 204413f72f0SBarry Smith ISLocalToGlobalMapping_Basic *map = (ISLocalToGlobalMapping_Basic *)mapping->data; 205413f72f0SBarry Smith 206413f72f0SBarry Smith PetscFunctionBegin; 207413f72f0SBarry Smith if (!map) PetscFunctionReturn(0); 2089566063dSJacob Faibussowitsch PetscCall(PetscFree(map->globals)); 2099566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->data)); 210413f72f0SBarry Smith PetscFunctionReturn(0); 211413f72f0SBarry Smith } 212413f72f0SBarry Smith 213*d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISLocalToGlobalMappingDestroy_Hash(ISLocalToGlobalMapping mapping) 214*d71ae5a4SJacob Faibussowitsch { 215413f72f0SBarry Smith ISLocalToGlobalMapping_Hash *map = (ISLocalToGlobalMapping_Hash *)mapping->data; 216413f72f0SBarry Smith 217413f72f0SBarry Smith PetscFunctionBegin; 218413f72f0SBarry Smith if (!map) PetscFunctionReturn(0); 2199566063dSJacob Faibussowitsch PetscCall(PetscHMapIDestroy(&map->globalht)); 2209566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->data)); 221413f72f0SBarry Smith PetscFunctionReturn(0); 222413f72f0SBarry Smith } 223413f72f0SBarry Smith 224413f72f0SBarry Smith #define GTOLTYPE _Basic 225413f72f0SBarry Smith #define GTOLNAME _Basic 226541bf97eSAdrian Croucher #define GTOLBS mapping->bs 2279371c9d4SSatish Balay #define GTOL(g, local) \ 2289371c9d4SSatish Balay do { \ 229413f72f0SBarry Smith local = map->globals[g / bs - start]; \ 2300040bde1SJunchao Zhang if (local >= 0) local = bs * local + (g % bs); \ 231413f72f0SBarry Smith } while (0) 232541bf97eSAdrian Croucher 233413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h> 234413f72f0SBarry Smith 235413f72f0SBarry Smith #define GTOLTYPE _Basic 236413f72f0SBarry Smith #define GTOLNAME Block_Basic 237541bf97eSAdrian Croucher #define GTOLBS 1 2389371c9d4SSatish Balay #define GTOL(g, local) \ 239*d71ae5a4SJacob Faibussowitsch do { \ 240*d71ae5a4SJacob Faibussowitsch local = map->globals[g - start]; \ 241*d71ae5a4SJacob Faibussowitsch } while (0) 242413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h> 243413f72f0SBarry Smith 244413f72f0SBarry Smith #define GTOLTYPE _Hash 245413f72f0SBarry Smith #define GTOLNAME _Hash 246541bf97eSAdrian Croucher #define GTOLBS mapping->bs 2479371c9d4SSatish Balay #define GTOL(g, local) \ 2489371c9d4SSatish Balay do { \ 249e8f14785SLisandro Dalcin (void)PetscHMapIGet(map->globalht, g / bs, &local); \ 2500040bde1SJunchao Zhang if (local >= 0) local = bs * local + (g % bs); \ 251413f72f0SBarry Smith } while (0) 252413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h> 253413f72f0SBarry Smith 254413f72f0SBarry Smith #define GTOLTYPE _Hash 255413f72f0SBarry Smith #define GTOLNAME Block_Hash 256541bf97eSAdrian Croucher #define GTOLBS 1 2579371c9d4SSatish Balay #define GTOL(g, local) \ 258*d71ae5a4SJacob Faibussowitsch do { \ 259*d71ae5a4SJacob Faibussowitsch (void)PetscHMapIGet(map->globalht, g, &local); \ 260*d71ae5a4SJacob Faibussowitsch } while (0) 261413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h> 262413f72f0SBarry Smith 2636658fb44Sstefano_zampini /*@ 2646658fb44Sstefano_zampini ISLocalToGlobalMappingDuplicate - Duplicates the local to global mapping object 2656658fb44Sstefano_zampini 2666658fb44Sstefano_zampini Not Collective 2676658fb44Sstefano_zampini 2686658fb44Sstefano_zampini Input Parameter: 2696658fb44Sstefano_zampini . ltog - local to global mapping 2706658fb44Sstefano_zampini 2716658fb44Sstefano_zampini Output Parameter: 2726658fb44Sstefano_zampini . nltog - the duplicated local to global mapping 2736658fb44Sstefano_zampini 2746658fb44Sstefano_zampini Level: advanced 2756658fb44Sstefano_zampini 276db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()` 2776658fb44Sstefano_zampini @*/ 278*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingDuplicate(ISLocalToGlobalMapping ltog, ISLocalToGlobalMapping *nltog) 279*d71ae5a4SJacob Faibussowitsch { 280a0d79125SStefano Zampini ISLocalToGlobalMappingType l2gtype; 2816658fb44Sstefano_zampini 2826658fb44Sstefano_zampini PetscFunctionBegin; 2836658fb44Sstefano_zampini PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1); 2849566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)ltog), ltog->bs, ltog->n, ltog->indices, PETSC_COPY_VALUES, nltog)); 2859566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetType(ltog, &l2gtype)); 2869566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingSetType(*nltog, l2gtype)); 2876658fb44Sstefano_zampini PetscFunctionReturn(0); 2886658fb44Sstefano_zampini } 2896658fb44Sstefano_zampini 290565245c5SBarry Smith /*@ 291107e9a97SBarry Smith ISLocalToGlobalMappingGetSize - Gets the local size of a local to global mapping 2923b9aefa3SBarry Smith 2933b9aefa3SBarry Smith Not Collective 2943b9aefa3SBarry Smith 2953b9aefa3SBarry Smith Input Parameter: 2963b9aefa3SBarry Smith . ltog - local to global mapping 2973b9aefa3SBarry Smith 2983b9aefa3SBarry Smith Output Parameter: 299107e9a97SBarry Smith . n - the number of entries in the local mapping, ISLocalToGlobalMappingGetIndices() returns an array of this length 3003b9aefa3SBarry Smith 3013b9aefa3SBarry Smith Level: advanced 3023b9aefa3SBarry Smith 303db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()` 3043b9aefa3SBarry Smith @*/ 305*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping mapping, PetscInt *n) 306*d71ae5a4SJacob Faibussowitsch { 3073b9aefa3SBarry Smith PetscFunctionBegin; 3080700a824SBarry Smith PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 3094482741eSBarry Smith PetscValidIntPointer(n, 2); 310107e9a97SBarry Smith *n = mapping->bs * mapping->n; 3113b9aefa3SBarry Smith PetscFunctionReturn(0); 3123b9aefa3SBarry Smith } 3133b9aefa3SBarry Smith 3145a5d4f66SBarry Smith /*@C 315fe2efc57SMark ISLocalToGlobalMappingViewFromOptions - View from Options 316fe2efc57SMark 317fe2efc57SMark Collective on ISLocalToGlobalMapping 318fe2efc57SMark 319fe2efc57SMark Input Parameters: 320fe2efc57SMark + A - the local to global mapping object 321736c3998SJose E. Roman . obj - Optional object 322736c3998SJose E. Roman - name - command line option 323fe2efc57SMark 324fe2efc57SMark Level: intermediate 325db781477SPatrick Sanan .seealso: `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingView`, `PetscObjectViewFromOptions()`, `ISLocalToGlobalMappingCreate()` 326fe2efc57SMark @*/ 327*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingViewFromOptions(ISLocalToGlobalMapping A, PetscObject obj, const char name[]) 328*d71ae5a4SJacob Faibussowitsch { 329fe2efc57SMark PetscFunctionBegin; 330fe2efc57SMark PetscValidHeaderSpecific(A, IS_LTOGM_CLASSID, 1); 3319566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name)); 332fe2efc57SMark PetscFunctionReturn(0); 333fe2efc57SMark } 334fe2efc57SMark 335fe2efc57SMark /*@C 3365a5d4f66SBarry Smith ISLocalToGlobalMappingView - View a local to global mapping 3375a5d4f66SBarry Smith 338b9cd556bSLois Curfman McInnes Not Collective 339b9cd556bSLois Curfman McInnes 3405a5d4f66SBarry Smith Input Parameters: 3413b9aefa3SBarry Smith + ltog - local to global mapping 3423b9aefa3SBarry Smith - viewer - viewer 3435a5d4f66SBarry Smith 344a997ad1aSLois Curfman McInnes Level: advanced 345a997ad1aSLois Curfman McInnes 346db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()` 3475a5d4f66SBarry Smith @*/ 348*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingView(ISLocalToGlobalMapping mapping, PetscViewer viewer) 349*d71ae5a4SJacob Faibussowitsch { 35032dcc486SBarry Smith PetscInt i; 35132dcc486SBarry Smith PetscMPIInt rank; 352ace3abfcSBarry Smith PetscBool iascii; 3535a5d4f66SBarry Smith 3545a5d4f66SBarry Smith PetscFunctionBegin; 3550700a824SBarry Smith PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 35648a46eb9SPierre Jolivet if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mapping), &viewer)); 3570700a824SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 3585a5d4f66SBarry Smith 3599566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mapping), &rank)); 3609566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 36132077d6dSBarry Smith if (iascii) { 3629566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)mapping, viewer)); 3639566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 36448a46eb9SPierre Jolivet for (i = 0; i < mapping->n; i++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %" PetscInt_FMT " %" PetscInt_FMT "\n", rank, i, mapping->indices[i])); 3659566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 3669566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 3671575c14dSBarry Smith } 3685a5d4f66SBarry Smith PetscFunctionReturn(0); 3695a5d4f66SBarry Smith } 3705a5d4f66SBarry Smith 3711f428162SBarry Smith /*@ 3722bdab257SBarry Smith ISLocalToGlobalMappingCreateIS - Creates a mapping between a local (0 to n) 3732bdab257SBarry Smith ordering and a global parallel ordering. 3742bdab257SBarry Smith 3750f5bd95cSBarry Smith Not collective 376b9cd556bSLois Curfman McInnes 377a997ad1aSLois Curfman McInnes Input Parameter: 3788c03b21aSDmitry Karpeev . is - index set containing the global numbers for each local number 3792bdab257SBarry Smith 380a997ad1aSLois Curfman McInnes Output Parameter: 3812bdab257SBarry Smith . mapping - new mapping data structure 3822bdab257SBarry Smith 38395452b02SPatrick Sanan Notes: 38495452b02SPatrick Sanan the block size of the IS determines the block size of the mapping 385a997ad1aSLois Curfman McInnes Level: advanced 386a997ad1aSLois Curfman McInnes 387db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetFromOptions()` 3882bdab257SBarry Smith @*/ 389*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingCreateIS(IS is, ISLocalToGlobalMapping *mapping) 390*d71ae5a4SJacob Faibussowitsch { 3913bbf0e92SBarry Smith PetscInt n, bs; 3925d0c19d7SBarry Smith const PetscInt *indices; 3932bdab257SBarry Smith MPI_Comm comm; 3943bbf0e92SBarry Smith PetscBool isblock; 3953a40ed3dSBarry Smith 3963a40ed3dSBarry Smith PetscFunctionBegin; 3970700a824SBarry Smith PetscValidHeaderSpecific(is, IS_CLASSID, 1); 3984482741eSBarry Smith PetscValidPointer(mapping, 2); 3992bdab257SBarry Smith 4009566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)is, &comm)); 4019566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &n)); 4029566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)is, ISBLOCK, &isblock)); 4036006e8d2SBarry Smith if (!isblock) { 4049566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &indices)); 4059566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm, 1, n, indices, PETSC_COPY_VALUES, mapping)); 4069566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &indices)); 4076006e8d2SBarry Smith } else { 4089566063dSJacob Faibussowitsch PetscCall(ISGetBlockSize(is, &bs)); 4099566063dSJacob Faibussowitsch PetscCall(ISBlockGetIndices(is, &indices)); 4109566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm, bs, n / bs, indices, PETSC_COPY_VALUES, mapping)); 4119566063dSJacob Faibussowitsch PetscCall(ISBlockRestoreIndices(is, &indices)); 4126006e8d2SBarry Smith } 4133a40ed3dSBarry Smith PetscFunctionReturn(0); 4142bdab257SBarry Smith } 4155a5d4f66SBarry Smith 416a4d96a55SJed Brown /*@C 417a4d96a55SJed Brown ISLocalToGlobalMappingCreateSF - Creates a mapping between a local (0 to n) 418a4d96a55SJed Brown ordering and a global parallel ordering. 419a4d96a55SJed Brown 420a4d96a55SJed Brown Collective 421a4d96a55SJed Brown 422d8d19677SJose E. Roman Input Parameters: 423a4d96a55SJed Brown + sf - star forest mapping contiguous local indices to (rank, offset) 4249a535bafSVaclav Hapla - start - first global index on this process, or PETSC_DECIDE to compute contiguous global numbering automatically 425a4d96a55SJed Brown 426a4d96a55SJed Brown Output Parameter: 427a4d96a55SJed Brown . mapping - new mapping data structure 428a4d96a55SJed Brown 429a4d96a55SJed Brown Level: advanced 430a4d96a55SJed Brown 4319a535bafSVaclav Hapla Notes: 4329a535bafSVaclav Hapla If any processor calls this with start = PETSC_DECIDE then all processors must, otherwise the program will hang. 4339a535bafSVaclav Hapla 434db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingSetFromOptions()` 435a4d96a55SJed Brown @*/ 436*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF sf, PetscInt start, ISLocalToGlobalMapping *mapping) 437*d71ae5a4SJacob Faibussowitsch { 438a4d96a55SJed Brown PetscInt i, maxlocal, nroots, nleaves, *globals, *ltog; 439a4d96a55SJed Brown MPI_Comm comm; 440a4d96a55SJed Brown 441a4d96a55SJed Brown PetscFunctionBegin; 442a4d96a55SJed Brown PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 443a4d96a55SJed Brown PetscValidPointer(mapping, 3); 4449566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 44541f4c31fSVaclav Hapla PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL)); 4469a535bafSVaclav Hapla if (start == PETSC_DECIDE) { 4479a535bafSVaclav Hapla start = 0; 4489566063dSJacob Faibussowitsch PetscCallMPI(MPI_Exscan(&nroots, &start, 1, MPIU_INT, MPI_SUM, comm)); 44941f4c31fSVaclav Hapla } else PetscCheck(start >= 0, comm, PETSC_ERR_ARG_OUTOFRANGE, "start must be nonnegative or PETSC_DECIDE"); 45041f4c31fSVaclav Hapla PetscCall(PetscSFGetLeafRange(sf, NULL, &maxlocal)); 45141f4c31fSVaclav Hapla ++maxlocal; 4529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nroots, &globals)); 4539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxlocal, <og)); 454a4d96a55SJed Brown for (i = 0; i < nroots; i++) globals[i] = start + i; 455a4d96a55SJed Brown for (i = 0; i < maxlocal; i++) ltog[i] = -1; 4569566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf, MPIU_INT, globals, ltog, MPI_REPLACE)); 4579566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf, MPIU_INT, globals, ltog, MPI_REPLACE)); 4589566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm, 1, maxlocal, ltog, PETSC_OWN_POINTER, mapping)); 4599566063dSJacob Faibussowitsch PetscCall(PetscFree(globals)); 460a4d96a55SJed Brown PetscFunctionReturn(0); 461a4d96a55SJed Brown } 462b46b645bSBarry Smith 46363fa5c83Sstefano_zampini /*@ 46463fa5c83Sstefano_zampini ISLocalToGlobalMappingSetBlockSize - Sets the blocksize of the mapping 46563fa5c83Sstefano_zampini 46663fa5c83Sstefano_zampini Not collective 46763fa5c83Sstefano_zampini 46863fa5c83Sstefano_zampini Input Parameters: 469a2b725a8SWilliam Gropp + mapping - mapping data structure 470a2b725a8SWilliam Gropp - bs - the blocksize 47163fa5c83Sstefano_zampini 47263fa5c83Sstefano_zampini Level: advanced 47363fa5c83Sstefano_zampini 474db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()` 47563fa5c83Sstefano_zampini @*/ 476*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingSetBlockSize(ISLocalToGlobalMapping mapping, PetscInt bs) 477*d71ae5a4SJacob Faibussowitsch { 478a59f3c4dSstefano_zampini PetscInt *nid; 479a59f3c4dSstefano_zampini const PetscInt *oid; 480a59f3c4dSstefano_zampini PetscInt i, cn, on, obs, nn; 48163fa5c83Sstefano_zampini 48263fa5c83Sstefano_zampini PetscFunctionBegin; 48363fa5c83Sstefano_zampini PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 48408401ef6SPierre Jolivet PetscCheck(bs >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid block size %" PetscInt_FMT, bs); 48563fa5c83Sstefano_zampini if (bs == mapping->bs) PetscFunctionReturn(0); 48663fa5c83Sstefano_zampini on = mapping->n; 48763fa5c83Sstefano_zampini obs = mapping->bs; 48863fa5c83Sstefano_zampini oid = mapping->indices; 48963fa5c83Sstefano_zampini nn = (on * obs) / bs; 49008401ef6SPierre Jolivet PetscCheck((on * obs) % bs == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Block size %" PetscInt_FMT " is inconsistent with block size %" PetscInt_FMT " and number of block indices %" PetscInt_FMT, bs, obs, on); 491a59f3c4dSstefano_zampini 4929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nn, &nid)); 4939566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(mapping, &oid)); 494a59f3c4dSstefano_zampini for (i = 0; i < nn; i++) { 495a59f3c4dSstefano_zampini PetscInt j; 496a59f3c4dSstefano_zampini for (j = 0, cn = 0; j < bs - 1; j++) { 4979371c9d4SSatish Balay if (oid[i * bs + j] < 0) { 4989371c9d4SSatish Balay cn++; 4999371c9d4SSatish Balay continue; 5009371c9d4SSatish Balay } 50108401ef6SPierre Jolivet PetscCheck(oid[i * bs + j] == oid[i * bs + j + 1] - 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Block sizes %" PetscInt_FMT " and %" PetscInt_FMT " are incompatible with the block indices: non consecutive indices %" PetscInt_FMT " %" PetscInt_FMT, bs, obs, oid[i * bs + j], oid[i * bs + j + 1]); 502a59f3c4dSstefano_zampini } 503a59f3c4dSstefano_zampini if (oid[i * bs + j] < 0) cn++; 5048b7cb0e6Sstefano_zampini if (cn) { 50508401ef6SPierre Jolivet PetscCheck(cn == bs, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Block sizes %" PetscInt_FMT " and %" PetscInt_FMT " are incompatible with the block indices: invalid number of negative entries in block %" PetscInt_FMT, bs, obs, cn); 506a59f3c4dSstefano_zampini nid[i] = -1; 5078b7cb0e6Sstefano_zampini } else { 508a59f3c4dSstefano_zampini nid[i] = oid[i * bs] / bs; 50963fa5c83Sstefano_zampini } 51063fa5c83Sstefano_zampini } 5119566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(mapping, &oid)); 512a59f3c4dSstefano_zampini 51363fa5c83Sstefano_zampini mapping->n = nn; 51463fa5c83Sstefano_zampini mapping->bs = bs; 5159566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->indices)); 51663fa5c83Sstefano_zampini mapping->indices = nid; 517c9345713Sstefano_zampini mapping->globalstart = 0; 518c9345713Sstefano_zampini mapping->globalend = 0; 5191bd0b88eSStefano Zampini 5201bd0b88eSStefano Zampini /* reset the cached information */ 5219566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->info_procs)); 5229566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->info_numprocs)); 5231bd0b88eSStefano Zampini if (mapping->info_indices) { 5241bd0b88eSStefano Zampini PetscInt i; 5251bd0b88eSStefano Zampini 5269566063dSJacob Faibussowitsch PetscCall(PetscFree((mapping->info_indices)[0])); 52748a46eb9SPierre Jolivet for (i = 1; i < mapping->info_nproc; i++) PetscCall(PetscFree(mapping->info_indices[i])); 5289566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->info_indices)); 5291bd0b88eSStefano Zampini } 5301bd0b88eSStefano Zampini mapping->info_cached = PETSC_FALSE; 5311bd0b88eSStefano Zampini 532dbbe0bcdSBarry Smith PetscTryTypeMethod(mapping, destroy); 53363fa5c83Sstefano_zampini PetscFunctionReturn(0); 53463fa5c83Sstefano_zampini } 53563fa5c83Sstefano_zampini 53645b6f7e9SBarry Smith /*@ 53745b6f7e9SBarry Smith ISLocalToGlobalMappingGetBlockSize - Gets the blocksize of the mapping 53845b6f7e9SBarry Smith ordering and a global parallel ordering. 53945b6f7e9SBarry Smith 54045b6f7e9SBarry Smith Not Collective 54145b6f7e9SBarry Smith 54245b6f7e9SBarry Smith Input Parameters: 54345b6f7e9SBarry Smith . mapping - mapping data structure 54445b6f7e9SBarry Smith 54545b6f7e9SBarry Smith Output Parameter: 54645b6f7e9SBarry Smith . bs - the blocksize 54745b6f7e9SBarry Smith 54845b6f7e9SBarry Smith Level: advanced 54945b6f7e9SBarry Smith 550db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()` 55145b6f7e9SBarry Smith @*/ 552*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping mapping, PetscInt *bs) 553*d71ae5a4SJacob Faibussowitsch { 55445b6f7e9SBarry Smith PetscFunctionBegin; 555cbc1caf0SMatthew G. Knepley PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 55645b6f7e9SBarry Smith *bs = mapping->bs; 55745b6f7e9SBarry Smith PetscFunctionReturn(0); 55845b6f7e9SBarry Smith } 55945b6f7e9SBarry Smith 560ba5bb76aSSatish Balay /*@ 56190f02eecSBarry Smith ISLocalToGlobalMappingCreate - Creates a mapping between a local (0 to n) 56290f02eecSBarry Smith ordering and a global parallel ordering. 5632362add9SBarry Smith 56489d82c54SBarry Smith Not Collective, but communicator may have more than one process 565b9cd556bSLois Curfman McInnes 5662362add9SBarry Smith Input Parameters: 56789d82c54SBarry Smith + comm - MPI communicator 568f0413b6fSBarry Smith . bs - the block size 56928bc9809SBarry Smith . n - the number of local elements divided by the block size, or equivalently the number of block indices 57028bc9809SBarry Smith . indices - the global index for each local element, these do not need to be in increasing order (sorted), these values should not be scaled (i.e. multiplied) by the blocksize bs 571d5ad8652SBarry Smith - mode - see PetscCopyMode 5722362add9SBarry Smith 573a997ad1aSLois Curfman McInnes Output Parameter: 57490f02eecSBarry Smith . mapping - new mapping data structure 5752362add9SBarry Smith 57695452b02SPatrick Sanan Notes: 57795452b02SPatrick Sanan There is one integer value in indices per block and it represents the actual indices bs*idx + j, where j=0,..,bs-1 578413f72f0SBarry Smith 5799a7b7924SJed Brown For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used; 580413f72f0SBarry Smith this uses more memory but is faster; this approach is not scalable for extremely large mappings. For large problems ISLOCALTOGLOBALMAPPINGHASH is used, this is scalable. 581413f72f0SBarry Smith Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used. 582413f72f0SBarry Smith 583a997ad1aSLois Curfman McInnes Level: advanced 584a997ad1aSLois Curfman McInnes 585db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingSetFromOptions()`, `ISLOCALTOGLOBALMAPPINGBASIC`, `ISLOCALTOGLOBALMAPPINGHASH` 586db781477SPatrick Sanan `ISLocalToGlobalMappingSetType()`, `ISLocalToGlobalMappingType` 5872362add9SBarry Smith @*/ 588*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingCreate(MPI_Comm comm, PetscInt bs, PetscInt n, const PetscInt indices[], PetscCopyMode mode, ISLocalToGlobalMapping *mapping) 589*d71ae5a4SJacob Faibussowitsch { 59032dcc486SBarry Smith PetscInt *in; 591b46b645bSBarry Smith 592b46b645bSBarry Smith PetscFunctionBegin; 593064a246eSJacob Faibussowitsch if (n) PetscValidIntPointer(indices, 4); 594064a246eSJacob Faibussowitsch PetscValidPointer(mapping, 6); 595b46b645bSBarry Smith 5960298fd71SBarry Smith *mapping = NULL; 5979566063dSJacob Faibussowitsch PetscCall(ISInitializePackage()); 5982362add9SBarry Smith 5999566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(*mapping, IS_LTOGM_CLASSID, "ISLocalToGlobalMapping", "Local to global mapping", "IS", comm, ISLocalToGlobalMappingDestroy, ISLocalToGlobalMappingView)); 600d4bb536fSBarry Smith (*mapping)->n = n; 601f0413b6fSBarry Smith (*mapping)->bs = bs; 602d5ad8652SBarry Smith if (mode == PETSC_COPY_VALUES) { 6039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &in)); 6049566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(in, indices, n)); 605d5ad8652SBarry Smith (*mapping)->indices = in; 60671910c26SVaclav Hapla (*mapping)->dealloc_indices = PETSC_TRUE; 6076389a1a1SBarry Smith } else if (mode == PETSC_OWN_POINTER) { 6086389a1a1SBarry Smith (*mapping)->indices = (PetscInt *)indices; 60971910c26SVaclav Hapla (*mapping)->dealloc_indices = PETSC_TRUE; 61071910c26SVaclav Hapla } else if (mode == PETSC_USE_POINTER) { 61171910c26SVaclav Hapla (*mapping)->indices = (PetscInt *)indices; 6129371c9d4SSatish Balay } else SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid mode %d", mode); 6133a40ed3dSBarry Smith PetscFunctionReturn(0); 6142362add9SBarry Smith } 6152362add9SBarry Smith 616413f72f0SBarry Smith PetscFunctionList ISLocalToGlobalMappingList = NULL; 617413f72f0SBarry Smith 61890f02eecSBarry Smith /*@ 6197e99dc12SLawrence Mitchell ISLocalToGlobalMappingSetFromOptions - Set mapping options from the options database. 6207e99dc12SLawrence Mitchell 6217e99dc12SLawrence Mitchell Not collective 6227e99dc12SLawrence Mitchell 6237e99dc12SLawrence Mitchell Input Parameters: 6247e99dc12SLawrence Mitchell . mapping - mapping data structure 6257e99dc12SLawrence Mitchell 6267e99dc12SLawrence Mitchell Level: advanced 6277e99dc12SLawrence Mitchell 6287e99dc12SLawrence Mitchell @*/ 629*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingSetFromOptions(ISLocalToGlobalMapping mapping) 630*d71ae5a4SJacob Faibussowitsch { 631413f72f0SBarry Smith char type[256]; 632413f72f0SBarry Smith ISLocalToGlobalMappingType defaulttype = "Not set"; 6337e99dc12SLawrence Mitchell PetscBool flg; 6347e99dc12SLawrence Mitchell 6357e99dc12SLawrence Mitchell PetscFunctionBegin; 6367e99dc12SLawrence Mitchell PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 6379566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRegisterAll()); 638d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)mapping); 6399566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-islocaltoglobalmapping_type", "ISLocalToGlobalMapping method", "ISLocalToGlobalMappingSetType", ISLocalToGlobalMappingList, (char *)(((PetscObject)mapping)->type_name) ? ((PetscObject)mapping)->type_name : defaulttype, type, 256, &flg)); 6401baa6e33SBarry Smith if (flg) PetscCall(ISLocalToGlobalMappingSetType(mapping, type)); 641d0609cedSBarry Smith PetscOptionsEnd(); 6427e99dc12SLawrence Mitchell PetscFunctionReturn(0); 6437e99dc12SLawrence Mitchell } 6447e99dc12SLawrence Mitchell 6457e99dc12SLawrence Mitchell /*@ 64690f02eecSBarry Smith ISLocalToGlobalMappingDestroy - Destroys a mapping between a local (0 to n) 64790f02eecSBarry Smith ordering and a global parallel ordering. 64890f02eecSBarry Smith 6490f5bd95cSBarry Smith Note Collective 650b9cd556bSLois Curfman McInnes 65190f02eecSBarry Smith Input Parameters: 65290f02eecSBarry Smith . mapping - mapping data structure 65390f02eecSBarry Smith 654a997ad1aSLois Curfman McInnes Level: advanced 655a997ad1aSLois Curfman McInnes 656db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()` 65790f02eecSBarry Smith @*/ 658*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping *mapping) 659*d71ae5a4SJacob Faibussowitsch { 6603a40ed3dSBarry Smith PetscFunctionBegin; 6616bf464f9SBarry Smith if (!*mapping) PetscFunctionReturn(0); 6626bf464f9SBarry Smith PetscValidHeaderSpecific((*mapping), IS_LTOGM_CLASSID, 1); 6639371c9d4SSatish Balay if (--((PetscObject)(*mapping))->refct > 0) { 6649371c9d4SSatish Balay *mapping = NULL; 6659371c9d4SSatish Balay PetscFunctionReturn(0); 66671910c26SVaclav Hapla } 66748a46eb9SPierre Jolivet if ((*mapping)->dealloc_indices) PetscCall(PetscFree((*mapping)->indices)); 6689566063dSJacob Faibussowitsch PetscCall(PetscFree((*mapping)->info_procs)); 6699566063dSJacob Faibussowitsch PetscCall(PetscFree((*mapping)->info_numprocs)); 670268a049cSStefano Zampini if ((*mapping)->info_indices) { 671268a049cSStefano Zampini PetscInt i; 672268a049cSStefano Zampini 6739566063dSJacob Faibussowitsch PetscCall(PetscFree(((*mapping)->info_indices)[0])); 67448a46eb9SPierre Jolivet for (i = 1; i < (*mapping)->info_nproc; i++) PetscCall(PetscFree(((*mapping)->info_indices)[i])); 6759566063dSJacob Faibussowitsch PetscCall(PetscFree((*mapping)->info_indices)); 676268a049cSStefano Zampini } 67748a46eb9SPierre Jolivet if ((*mapping)->info_nodei) PetscCall(PetscFree(((*mapping)->info_nodei)[0])); 6789566063dSJacob Faibussowitsch PetscCall(PetscFree2((*mapping)->info_nodec, (*mapping)->info_nodei)); 67948a46eb9SPierre Jolivet if ((*mapping)->ops->destroy) PetscCall((*(*mapping)->ops->destroy)(*mapping)); 6809566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(mapping)); 6814c8fdceaSLisandro Dalcin *mapping = NULL; 6823a40ed3dSBarry Smith PetscFunctionReturn(0); 68390f02eecSBarry Smith } 68490f02eecSBarry Smith 68590f02eecSBarry Smith /*@ 6863acfe500SLois Curfman McInnes ISLocalToGlobalMappingApplyIS - Creates from an IS in the local numbering 6873acfe500SLois Curfman McInnes a new index set using the global numbering defined in an ISLocalToGlobalMapping 6883acfe500SLois Curfman McInnes context. 68990f02eecSBarry Smith 6904cb36875SStefano Zampini Collective on is 691b9cd556bSLois Curfman McInnes 69290f02eecSBarry Smith Input Parameters: 693b9cd556bSLois Curfman McInnes + mapping - mapping between local and global numbering 694b9cd556bSLois Curfman McInnes - is - index set in local numbering 69590f02eecSBarry Smith 69690f02eecSBarry Smith Output Parameters: 69790f02eecSBarry Smith . newis - index set in global numbering 69890f02eecSBarry Smith 6994cb36875SStefano Zampini Notes: 7004cb36875SStefano Zampini The output IS will have the same communicator of the input IS. 7014cb36875SStefano Zampini 702a997ad1aSLois Curfman McInnes Level: advanced 703a997ad1aSLois Curfman McInnes 704db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingCreate()`, 705db781477SPatrick Sanan `ISLocalToGlobalMappingDestroy()`, `ISGlobalToLocalMappingApply()` 70690f02eecSBarry Smith @*/ 707*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping, IS is, IS *newis) 708*d71ae5a4SJacob Faibussowitsch { 709e24637baSBarry Smith PetscInt n, *idxout; 7105d0c19d7SBarry Smith const PetscInt *idxin; 7113a40ed3dSBarry Smith 7123a40ed3dSBarry Smith PetscFunctionBegin; 7130700a824SBarry Smith PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 7140700a824SBarry Smith PetscValidHeaderSpecific(is, IS_CLASSID, 2); 7154482741eSBarry Smith PetscValidPointer(newis, 3); 71690f02eecSBarry Smith 7179566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &n)); 7189566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &idxin)); 7199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &idxout)); 7209566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(mapping, n, idxin, idxout)); 7219566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &idxin)); 7229566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)is), n, idxout, PETSC_OWN_POINTER, newis)); 7233a40ed3dSBarry Smith PetscFunctionReturn(0); 72490f02eecSBarry Smith } 72590f02eecSBarry Smith 726b89cb25eSSatish Balay /*@ 7273acfe500SLois Curfman McInnes ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering 7283acfe500SLois Curfman McInnes and converts them to the global numbering. 72990f02eecSBarry Smith 730b9cd556bSLois Curfman McInnes Not collective 731b9cd556bSLois Curfman McInnes 732bb25748dSBarry Smith Input Parameters: 733b9cd556bSLois Curfman McInnes + mapping - the local to global mapping context 734bb25748dSBarry Smith . N - number of integers 735b9cd556bSLois Curfman McInnes - in - input indices in local numbering 736bb25748dSBarry Smith 737bb25748dSBarry Smith Output Parameter: 738bb25748dSBarry Smith . out - indices in global numbering 739bb25748dSBarry Smith 740b9cd556bSLois Curfman McInnes Notes: 741b9cd556bSLois Curfman McInnes The in and out array parameters may be identical. 742d4bb536fSBarry Smith 743a997ad1aSLois Curfman McInnes Level: advanced 744a997ad1aSLois Curfman McInnes 745c2e3fba1SPatrick Sanan .seealso: `ISLocalToGlobalMappingApplyBlock()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingDestroy()`, 746c2e3fba1SPatrick Sanan `ISLocalToGlobalMappingApplyIS()`, `AOCreateBasic()`, `AOApplicationToPetsc()`, 747db781477SPatrick Sanan `AOPetscToApplication()`, `ISGlobalToLocalMappingApply()` 748bb25748dSBarry Smith 749afcb2eb5SJed Brown @*/ 750*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping, PetscInt N, const PetscInt in[], PetscInt out[]) 751*d71ae5a4SJacob Faibussowitsch { 752cbc1caf0SMatthew G. Knepley PetscInt i, bs, Nmax; 75345b6f7e9SBarry Smith 75445b6f7e9SBarry Smith PetscFunctionBegin; 755cbc1caf0SMatthew G. Knepley PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 756cbc1caf0SMatthew G. Knepley bs = mapping->bs; 757cbc1caf0SMatthew G. Knepley Nmax = bs * mapping->n; 75845b6f7e9SBarry Smith if (bs == 1) { 759cbc1caf0SMatthew G. Knepley const PetscInt *idx = mapping->indices; 76045b6f7e9SBarry Smith for (i = 0; i < N; i++) { 76145b6f7e9SBarry Smith if (in[i] < 0) { 76245b6f7e9SBarry Smith out[i] = in[i]; 76345b6f7e9SBarry Smith continue; 76445b6f7e9SBarry Smith } 76508401ef6SPierre Jolivet PetscCheck(in[i] < Nmax, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local index %" PetscInt_FMT " too large %" PetscInt_FMT " (max) at %" PetscInt_FMT, in[i], Nmax - 1, i); 76645b6f7e9SBarry Smith out[i] = idx[in[i]]; 76745b6f7e9SBarry Smith } 76845b6f7e9SBarry Smith } else { 769cbc1caf0SMatthew G. Knepley const PetscInt *idx = mapping->indices; 77045b6f7e9SBarry Smith for (i = 0; i < N; i++) { 77145b6f7e9SBarry Smith if (in[i] < 0) { 77245b6f7e9SBarry Smith out[i] = in[i]; 77345b6f7e9SBarry Smith continue; 77445b6f7e9SBarry Smith } 77508401ef6SPierre Jolivet PetscCheck(in[i] < Nmax, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local index %" PetscInt_FMT " too large %" PetscInt_FMT " (max) at %" PetscInt_FMT, in[i], Nmax - 1, i); 77645b6f7e9SBarry Smith out[i] = idx[in[i] / bs] * bs + (in[i] % bs); 77745b6f7e9SBarry Smith } 77845b6f7e9SBarry Smith } 77945b6f7e9SBarry Smith PetscFunctionReturn(0); 78045b6f7e9SBarry Smith } 78145b6f7e9SBarry Smith 78245b6f7e9SBarry Smith /*@ 7836006e8d2SBarry Smith ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local block numbering and converts them to the global block numbering 78445b6f7e9SBarry Smith 78545b6f7e9SBarry Smith Not collective 78645b6f7e9SBarry Smith 78745b6f7e9SBarry Smith Input Parameters: 78845b6f7e9SBarry Smith + mapping - the local to global mapping context 78945b6f7e9SBarry Smith . N - number of integers 7906006e8d2SBarry Smith - in - input indices in local block numbering 79145b6f7e9SBarry Smith 79245b6f7e9SBarry Smith Output Parameter: 7936006e8d2SBarry Smith . out - indices in global block numbering 79445b6f7e9SBarry Smith 79545b6f7e9SBarry Smith Notes: 79645b6f7e9SBarry Smith The in and out array parameters may be identical. 79745b6f7e9SBarry Smith 7986006e8d2SBarry Smith Example: 7996006e8d2SBarry 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 8006006e8d2SBarry Smith (the first block) would produce 0 and the mapping applied to 1 (the second block) would produce 3. 8016006e8d2SBarry Smith 80245b6f7e9SBarry Smith Level: advanced 80345b6f7e9SBarry Smith 804c2e3fba1SPatrick Sanan .seealso: `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingDestroy()`, 805c2e3fba1SPatrick Sanan `ISLocalToGlobalMappingApplyIS()`, `AOCreateBasic()`, `AOApplicationToPetsc()`, 806db781477SPatrick Sanan `AOPetscToApplication()`, `ISGlobalToLocalMappingApply()` 80745b6f7e9SBarry Smith 80845b6f7e9SBarry Smith @*/ 809*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping, PetscInt N, const PetscInt in[], PetscInt out[]) 810*d71ae5a4SJacob Faibussowitsch { 8118a1f772fSStefano Zampini PetscInt i, Nmax; 8128a1f772fSStefano Zampini const PetscInt *idx; 813d4bb536fSBarry Smith 814a0d79125SStefano Zampini PetscFunctionBegin; 815a0d79125SStefano Zampini PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 8168a1f772fSStefano Zampini Nmax = mapping->n; 8178a1f772fSStefano Zampini idx = mapping->indices; 818afcb2eb5SJed Brown for (i = 0; i < N; i++) { 819afcb2eb5SJed Brown if (in[i] < 0) { 820afcb2eb5SJed Brown out[i] = in[i]; 821afcb2eb5SJed Brown continue; 822afcb2eb5SJed Brown } 82308401ef6SPierre Jolivet PetscCheck(in[i] < Nmax, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local block index %" PetscInt_FMT " too large %" PetscInt_FMT " (max) at %" PetscInt_FMT, in[i], Nmax - 1, i); 824afcb2eb5SJed Brown out[i] = idx[in[i]]; 825afcb2eb5SJed Brown } 826afcb2eb5SJed Brown PetscFunctionReturn(0); 827afcb2eb5SJed Brown } 828d4bb536fSBarry Smith 8297e99dc12SLawrence Mitchell /*@ 830a997ad1aSLois Curfman McInnes ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers 831a997ad1aSLois Curfman McInnes specified with a global numbering. 832d4bb536fSBarry Smith 833b9cd556bSLois Curfman McInnes Not collective 834b9cd556bSLois Curfman McInnes 835d4bb536fSBarry Smith Input Parameters: 836b9cd556bSLois Curfman McInnes + mapping - mapping between local and global numbering 8370040bde1SJunchao Zhang . type - IS_GTOLM_MASK - maps global indices with no local value to -1 in the output list (i.e., mask them) 838d4bb536fSBarry Smith IS_GTOLM_DROP - drops the indices with no local value from the output list 839d4bb536fSBarry Smith . n - number of global indices to map 840b9cd556bSLois Curfman McInnes - idx - global indices to map 841d4bb536fSBarry Smith 842d4bb536fSBarry Smith Output Parameters: 843b9cd556bSLois Curfman McInnes + nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n) 844b9cd556bSLois Curfman McInnes - idxout - local index of each global index, one must pass in an array long enough 845e182c471SBarry Smith to hold all the indices. You can call ISGlobalToLocalMappingApply() with 8460298fd71SBarry Smith idxout == NULL to determine the required length (returned in nout) 847e182c471SBarry Smith and then allocate the required space and call ISGlobalToLocalMappingApply() 848e182c471SBarry Smith a second time to set the values. 849d4bb536fSBarry Smith 850b9cd556bSLois Curfman McInnes Notes: 8510298fd71SBarry Smith Either nout or idxout may be NULL. idx and idxout may be identical. 852d4bb536fSBarry Smith 8539a7b7924SJed Brown For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used; 854413f72f0SBarry 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. 855413f72f0SBarry Smith Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used. 8560f5bd95cSBarry Smith 857a997ad1aSLois Curfman McInnes Level: advanced 858a997ad1aSLois Curfman McInnes 85932fd6b96SBarry Smith Developer Note: The manual page states that idx and idxout may be identical but the calling 86032fd6b96SBarry Smith sequence declares idx as const so it cannot be the same as idxout. 86132fd6b96SBarry Smith 862db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingApply()`, `ISGlobalToLocalMappingApplyBlock()`, `ISLocalToGlobalMappingCreate()`, 863db781477SPatrick Sanan `ISLocalToGlobalMappingDestroy()` 864d4bb536fSBarry Smith @*/ 865*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, PetscInt n, const PetscInt idx[], PetscInt *nout, PetscInt idxout[]) 866*d71ae5a4SJacob Faibussowitsch { 8679d90f715SBarry Smith PetscFunctionBegin; 8689d90f715SBarry Smith PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 86948a46eb9SPierre Jolivet if (!mapping->data) PetscCall(ISGlobalToLocalMappingSetUp(mapping)); 870dbbe0bcdSBarry Smith PetscUseTypeMethod(mapping, globaltolocalmappingapply, type, n, idx, nout, idxout); 8719d90f715SBarry Smith PetscFunctionReturn(0); 8729d90f715SBarry Smith } 8739d90f715SBarry Smith 874d4fe737eSStefano Zampini /*@ 875d4fe737eSStefano Zampini ISGlobalToLocalMappingApplyIS - Creates from an IS in the global numbering 876d4fe737eSStefano Zampini a new index set using the local numbering defined in an ISLocalToGlobalMapping 877d4fe737eSStefano Zampini context. 878d4fe737eSStefano Zampini 879d4fe737eSStefano Zampini Not collective 880d4fe737eSStefano Zampini 881d4fe737eSStefano Zampini Input Parameters: 882d4fe737eSStefano Zampini + mapping - mapping between local and global numbering 8830040bde1SJunchao Zhang . type - IS_GTOLM_MASK - maps global indices with no local value to -1 in the output list (i.e., mask them) 8842785b321SStefano Zampini IS_GTOLM_DROP - drops the indices with no local value from the output list 885d4fe737eSStefano Zampini - is - index set in global numbering 886d4fe737eSStefano Zampini 887d4fe737eSStefano Zampini Output Parameters: 888d4fe737eSStefano Zampini . newis - index set in local numbering 889d4fe737eSStefano Zampini 8904cb36875SStefano Zampini Notes: 8914cb36875SStefano Zampini The output IS will be sequential, as it encodes a purely local operation 8924cb36875SStefano Zampini 893d4fe737eSStefano Zampini Level: advanced 894d4fe737eSStefano Zampini 895db781477SPatrick Sanan .seealso: `ISGlobalToLocalMappingApply()`, `ISLocalToGlobalMappingCreate()`, 896db781477SPatrick Sanan `ISLocalToGlobalMappingDestroy()` 897d4fe737eSStefano Zampini @*/ 898*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, IS is, IS *newis) 899*d71ae5a4SJacob Faibussowitsch { 900d4fe737eSStefano Zampini PetscInt n, nout, *idxout; 901d4fe737eSStefano Zampini const PetscInt *idxin; 902d4fe737eSStefano Zampini 903d4fe737eSStefano Zampini PetscFunctionBegin; 904d4fe737eSStefano Zampini PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 905d4fe737eSStefano Zampini PetscValidHeaderSpecific(is, IS_CLASSID, 3); 906d4fe737eSStefano Zampini PetscValidPointer(newis, 4); 907d4fe737eSStefano Zampini 9089566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &n)); 9099566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &idxin)); 910d4fe737eSStefano Zampini if (type == IS_GTOLM_MASK) { 9119566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &idxout)); 912d4fe737eSStefano Zampini } else { 9139566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(mapping, type, n, idxin, &nout, NULL)); 9149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nout, &idxout)); 915d4fe737eSStefano Zampini } 9169566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(mapping, type, n, idxin, &nout, idxout)); 9179566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &idxin)); 9189566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nout, idxout, PETSC_OWN_POINTER, newis)); 919d4fe737eSStefano Zampini PetscFunctionReturn(0); 920d4fe737eSStefano Zampini } 921d4fe737eSStefano Zampini 9229d90f715SBarry Smith /*@ 9239d90f715SBarry Smith ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers 9249d90f715SBarry Smith specified with a block global numbering. 9259d90f715SBarry Smith 9269d90f715SBarry Smith Not collective 9279d90f715SBarry Smith 9289d90f715SBarry Smith Input Parameters: 9299d90f715SBarry Smith + mapping - mapping between local and global numbering 9300040bde1SJunchao Zhang . type - IS_GTOLM_MASK - maps global indices with no local value to -1 in the output list (i.e., mask them) 9319d90f715SBarry Smith IS_GTOLM_DROP - drops the indices with no local value from the output list 9329d90f715SBarry Smith . n - number of global indices to map 9339d90f715SBarry Smith - idx - global indices to map 9349d90f715SBarry Smith 9359d90f715SBarry Smith Output Parameters: 9369d90f715SBarry Smith + nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n) 9379d90f715SBarry Smith - idxout - local index of each global index, one must pass in an array long enough 9389d90f715SBarry Smith to hold all the indices. You can call ISGlobalToLocalMappingApplyBlock() with 9399d90f715SBarry Smith idxout == NULL to determine the required length (returned in nout) 9409d90f715SBarry Smith and then allocate the required space and call ISGlobalToLocalMappingApplyBlock() 9419d90f715SBarry Smith a second time to set the values. 9429d90f715SBarry Smith 9439d90f715SBarry Smith Notes: 9449d90f715SBarry Smith Either nout or idxout may be NULL. idx and idxout may be identical. 9459d90f715SBarry Smith 9469a7b7924SJed Brown For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used; 947413f72f0SBarry 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. 948413f72f0SBarry Smith Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used. 9499d90f715SBarry Smith 9509d90f715SBarry Smith Level: advanced 9519d90f715SBarry Smith 9529d90f715SBarry Smith Developer Note: The manual page states that idx and idxout may be identical but the calling 9539d90f715SBarry Smith sequence declares idx as const so it cannot be the same as idxout. 9549d90f715SBarry Smith 955db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingApply()`, `ISGlobalToLocalMappingApply()`, `ISLocalToGlobalMappingCreate()`, 956db781477SPatrick Sanan `ISLocalToGlobalMappingDestroy()` 9579d90f715SBarry Smith @*/ 958*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, PetscInt n, const PetscInt idx[], PetscInt *nout, PetscInt idxout[]) 959*d71ae5a4SJacob Faibussowitsch { 9603a40ed3dSBarry Smith PetscFunctionBegin; 9610700a824SBarry Smith PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 96248a46eb9SPierre Jolivet if (!mapping->data) PetscCall(ISGlobalToLocalMappingSetUp(mapping)); 963dbbe0bcdSBarry Smith PetscUseTypeMethod(mapping, globaltolocalmappingapplyblock, type, n, idx, nout, idxout); 9643a40ed3dSBarry Smith PetscFunctionReturn(0); 965d4bb536fSBarry Smith } 96690f02eecSBarry Smith 96789d82c54SBarry Smith /*@C 9686a818285SBarry Smith ISLocalToGlobalMappingGetBlockInfo - Gets the neighbor information for each processor and 96989d82c54SBarry Smith each index shared by more than one processor 97089d82c54SBarry Smith 97189d82c54SBarry Smith Collective on ISLocalToGlobalMapping 97289d82c54SBarry Smith 973f899ff85SJose E. Roman Input Parameter: 97489d82c54SBarry Smith . mapping - the mapping from local to global indexing 97589d82c54SBarry Smith 976d8d19677SJose E. Roman Output Parameters: 97789d82c54SBarry Smith + nproc - number of processors that are connected to this one 97889d82c54SBarry Smith . proc - neighboring processors 97907b52d57SBarry Smith . numproc - number of indices for each subdomain (processor) 9803463a7baSJed Brown - indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering) 98189d82c54SBarry Smith 98289d82c54SBarry Smith Level: advanced 98389d82c54SBarry Smith 9842cfcea29SBarry Smith Fortran Usage: 9852cfcea29SBarry Smith $ ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by 9862cfcea29SBarry Smith $ ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc], 9872cfcea29SBarry Smith PetscInt indices[nproc][numprocmax],ierr) 9882cfcea29SBarry Smith There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and 9892cfcea29SBarry Smith indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough. 9902cfcea29SBarry Smith 991db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 992db781477SPatrick Sanan `ISLocalToGlobalMappingRestoreInfo()` 99389d82c54SBarry Smith @*/ 994*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[]) 995*d71ae5a4SJacob Faibussowitsch { 996268a049cSStefano Zampini PetscFunctionBegin; 997268a049cSStefano Zampini PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 998268a049cSStefano Zampini if (mapping->info_cached) { 999268a049cSStefano Zampini *nproc = mapping->info_nproc; 1000268a049cSStefano Zampini *procs = mapping->info_procs; 1001268a049cSStefano Zampini *numprocs = mapping->info_numprocs; 1002268a049cSStefano Zampini *indices = mapping->info_indices; 1003268a049cSStefano Zampini } else { 10049566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockInfo_Private(mapping, nproc, procs, numprocs, indices)); 1005268a049cSStefano Zampini } 1006268a049cSStefano Zampini PetscFunctionReturn(0); 1007268a049cSStefano Zampini } 1008268a049cSStefano Zampini 1009*d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[]) 1010*d71ae5a4SJacob Faibussowitsch { 101197f1f81fSBarry Smith PetscMPIInt size, rank, tag1, tag2, tag3, *len, *source, imdex; 101232dcc486SBarry Smith PetscInt i, n = mapping->n, Ng, ng, max = 0, *lindices = mapping->indices; 101332dcc486SBarry Smith PetscInt *nprocs, *owner, nsends, *sends, j, *starts, nmax, nrecvs, *recvs, proc; 1014c599c493SJunchao Zhang PetscInt cnt, scale, *ownedsenders, *nownedsenders, rstart; 101532dcc486SBarry Smith PetscInt node, nownedm, nt, *sends2, nsends2, *starts2, *lens2, *dest, nrecvs2, *starts3, *recvs2, k, *bprocs, *tmp; 101632dcc486SBarry Smith PetscInt first_procs, first_numprocs, *first_indices; 101789d82c54SBarry Smith MPI_Request *recv_waits, *send_waits; 101830dcb7c9SBarry Smith MPI_Status recv_status, *send_status, *recv_statuses; 1019ce94432eSBarry Smith MPI_Comm comm; 1020ace3abfcSBarry Smith PetscBool debug = PETSC_FALSE; 102189d82c54SBarry Smith 102289d82c54SBarry Smith PetscFunctionBegin; 10239566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)mapping, &comm)); 10249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 10259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 102624cf384cSBarry Smith if (size == 1) { 102724cf384cSBarry Smith *nproc = 0; 10280298fd71SBarry Smith *procs = NULL; 10299566063dSJacob Faibussowitsch PetscCall(PetscNew(numprocs)); 10301e2105dcSBarry Smith (*numprocs)[0] = 0; 10319566063dSJacob Faibussowitsch PetscCall(PetscNew(indices)); 10320298fd71SBarry Smith (*indices)[0] = NULL; 1033268a049cSStefano Zampini /* save info for reuse */ 1034268a049cSStefano Zampini mapping->info_nproc = *nproc; 1035268a049cSStefano Zampini mapping->info_procs = *procs; 1036268a049cSStefano Zampini mapping->info_numprocs = *numprocs; 1037268a049cSStefano Zampini mapping->info_indices = *indices; 1038268a049cSStefano Zampini mapping->info_cached = PETSC_TRUE; 103924cf384cSBarry Smith PetscFunctionReturn(0); 104024cf384cSBarry Smith } 104124cf384cSBarry Smith 10429566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)mapping)->options, NULL, "-islocaltoglobalmappinggetinfo_debug", &debug, NULL)); 104307b52d57SBarry Smith 10443677ff5aSBarry Smith /* 10456a818285SBarry Smith Notes on ISLocalToGlobalMappingGetBlockInfo 10463677ff5aSBarry Smith 10473677ff5aSBarry Smith globally owned node - the nodes that have been assigned to this processor in global 10483677ff5aSBarry Smith numbering, just for this routine. 10493677ff5aSBarry Smith 10503677ff5aSBarry Smith nontrivial globally owned node - node assigned to this processor that is on a subdomain 10513677ff5aSBarry Smith boundary (i.e. is has more than one local owner) 10523677ff5aSBarry Smith 10533677ff5aSBarry Smith locally owned node - node that exists on this processors subdomain 10543677ff5aSBarry Smith 10553677ff5aSBarry Smith nontrivial locally owned node - node that is not in the interior (i.e. has more than one 10563677ff5aSBarry Smith local subdomain 10573677ff5aSBarry Smith */ 10589566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)mapping, &tag1)); 10599566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)mapping, &tag2)); 10609566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)mapping, &tag3)); 106189d82c54SBarry Smith 106289d82c54SBarry Smith for (i = 0; i < n; i++) { 106389d82c54SBarry Smith if (lindices[i] > max) max = lindices[i]; 106489d82c54SBarry Smith } 10651c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&max, &Ng, 1, MPIU_INT, MPI_MAX, comm)); 106678058e43SBarry Smith Ng++; 10679566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 10689566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1069bc8ff85bSBarry Smith scale = Ng / size + 1; 10709371c9d4SSatish Balay ng = scale; 10719371c9d4SSatish Balay if (rank == size - 1) ng = Ng - scale * (size - 1); 10729371c9d4SSatish Balay ng = PetscMax(1, ng); 1073caba0dd0SBarry Smith rstart = scale * rank; 107489d82c54SBarry Smith 107589d82c54SBarry Smith /* determine ownership ranges of global indices */ 10769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * size, &nprocs)); 10779566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(nprocs, 2 * size)); 107889d82c54SBarry Smith 107989d82c54SBarry Smith /* determine owners of each local node */ 10809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &owner)); 108189d82c54SBarry Smith for (i = 0; i < n; i++) { 10823677ff5aSBarry Smith proc = lindices[i] / scale; /* processor that globally owns this index */ 108327c402fcSBarry Smith nprocs[2 * proc + 1] = 1; /* processor globally owns at least one of ours */ 10843677ff5aSBarry Smith owner[i] = proc; 108527c402fcSBarry Smith nprocs[2 * proc]++; /* count of how many that processor globally owns of ours */ 108689d82c54SBarry Smith } 10879371c9d4SSatish Balay nsends = 0; 10889371c9d4SSatish Balay for (i = 0; i < size; i++) nsends += nprocs[2 * i + 1]; 10899566063dSJacob Faibussowitsch PetscCall(PetscInfo(mapping, "Number of global owners for my local data %" PetscInt_FMT "\n", nsends)); 109089d82c54SBarry Smith 109189d82c54SBarry Smith /* inform other processors of number of messages and max length*/ 10929566063dSJacob Faibussowitsch PetscCall(PetscMaxSum(comm, nprocs, &nmax, &nrecvs)); 10939566063dSJacob Faibussowitsch PetscCall(PetscInfo(mapping, "Number of local owners for my global data %" PetscInt_FMT "\n", nrecvs)); 109489d82c54SBarry Smith 109589d82c54SBarry Smith /* post receives for owned rows */ 10969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((2 * nrecvs + 1) * (nmax + 1), &recvs)); 10979566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs + 1, &recv_waits)); 109848a46eb9SPierre Jolivet for (i = 0; i < nrecvs; i++) PetscCallMPI(MPI_Irecv(recvs + 2 * nmax * i, 2 * nmax, MPIU_INT, MPI_ANY_SOURCE, tag1, comm, recv_waits + i)); 109989d82c54SBarry Smith 110089d82c54SBarry Smith /* pack messages containing lists of local nodes to owners */ 11019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * n + 1, &sends)); 11029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size + 1, &starts)); 110389d82c54SBarry Smith starts[0] = 0; 1104f6e5521dSKarl Rupp for (i = 1; i < size; i++) starts[i] = starts[i - 1] + 2 * nprocs[2 * i - 2]; 110589d82c54SBarry Smith for (i = 0; i < n; i++) { 110689d82c54SBarry Smith sends[starts[owner[i]]++] = lindices[i]; 110730dcb7c9SBarry Smith sends[starts[owner[i]]++] = i; 110889d82c54SBarry Smith } 11099566063dSJacob Faibussowitsch PetscCall(PetscFree(owner)); 111089d82c54SBarry Smith starts[0] = 0; 1111f6e5521dSKarl Rupp for (i = 1; i < size; i++) starts[i] = starts[i - 1] + 2 * nprocs[2 * i - 2]; 111289d82c54SBarry Smith 111389d82c54SBarry Smith /* send the messages */ 11149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends + 1, &send_waits)); 11159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends + 1, &dest)); 111689d82c54SBarry Smith cnt = 0; 111789d82c54SBarry Smith for (i = 0; i < size; i++) { 111827c402fcSBarry Smith if (nprocs[2 * i]) { 11199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Isend(sends + starts[i], 2 * nprocs[2 * i], MPIU_INT, i, tag1, comm, send_waits + cnt)); 112030dcb7c9SBarry Smith dest[cnt] = i; 112189d82c54SBarry Smith cnt++; 112289d82c54SBarry Smith } 112389d82c54SBarry Smith } 11249566063dSJacob Faibussowitsch PetscCall(PetscFree(starts)); 112589d82c54SBarry Smith 112689d82c54SBarry Smith /* wait on receives */ 11279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs + 1, &source)); 11289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs + 1, &len)); 112989d82c54SBarry Smith cnt = nrecvs; 11309566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(ng + 1, &nownedsenders)); 113189d82c54SBarry Smith while (cnt) { 11329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitany(nrecvs, recv_waits, &imdex, &recv_status)); 113389d82c54SBarry Smith /* unpack receives into our local space */ 11349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&recv_status, MPIU_INT, &len[imdex])); 113589d82c54SBarry Smith source[imdex] = recv_status.MPI_SOURCE; 113630dcb7c9SBarry Smith len[imdex] = len[imdex] / 2; 1137caba0dd0SBarry Smith /* count how many local owners for each of my global owned indices */ 113830dcb7c9SBarry Smith for (i = 0; i < len[imdex]; i++) nownedsenders[recvs[2 * imdex * nmax + 2 * i] - rstart]++; 113989d82c54SBarry Smith cnt--; 114089d82c54SBarry Smith } 11419566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_waits)); 114289d82c54SBarry Smith 114330dcb7c9SBarry Smith /* count how many globally owned indices are on an edge multiplied by how many processors own them. */ 1144bc8ff85bSBarry Smith nownedm = 0; 1145bc8ff85bSBarry Smith for (i = 0; i < ng; i++) { 1146c599c493SJunchao Zhang if (nownedsenders[i] > 1) nownedm += nownedsenders[i]; 1147bc8ff85bSBarry Smith } 1148bc8ff85bSBarry Smith 1149bc8ff85bSBarry Smith /* create single array to contain rank of all local owners of each globally owned index */ 11509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nownedm + 1, &ownedsenders)); 11519566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ng + 1, &starts)); 1152bc8ff85bSBarry Smith starts[0] = 0; 1153bc8ff85bSBarry Smith for (i = 1; i < ng; i++) { 1154bc8ff85bSBarry Smith if (nownedsenders[i - 1] > 1) starts[i] = starts[i - 1] + nownedsenders[i - 1]; 1155bc8ff85bSBarry Smith else starts[i] = starts[i - 1]; 1156bc8ff85bSBarry Smith } 1157bc8ff85bSBarry Smith 11586aad120cSJose E. Roman /* for each nontrivial globally owned node list all arriving processors */ 1159bc8ff85bSBarry Smith for (i = 0; i < nrecvs; i++) { 1160bc8ff85bSBarry Smith for (j = 0; j < len[i]; j++) { 116130dcb7c9SBarry Smith node = recvs[2 * i * nmax + 2 * j] - rstart; 1162f6e5521dSKarl Rupp if (nownedsenders[node] > 1) ownedsenders[starts[node]++] = source[i]; 1163bc8ff85bSBarry Smith } 1164bc8ff85bSBarry Smith } 1165bc8ff85bSBarry Smith 116607b52d57SBarry Smith if (debug) { /* ----------------------------------- */ 116730dcb7c9SBarry Smith starts[0] = 0; 116830dcb7c9SBarry Smith for (i = 1; i < ng; i++) { 116930dcb7c9SBarry Smith if (nownedsenders[i - 1] > 1) starts[i] = starts[i - 1] + nownedsenders[i - 1]; 117030dcb7c9SBarry Smith else starts[i] = starts[i - 1]; 117130dcb7c9SBarry Smith } 117230dcb7c9SBarry Smith for (i = 0; i < ng; i++) { 117330dcb7c9SBarry Smith if (nownedsenders[i] > 1) { 11749566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm, "[%d] global node %" PetscInt_FMT " local owner processors: ", rank, i + rstart)); 117548a46eb9SPierre Jolivet for (j = 0; j < nownedsenders[i]; j++) PetscCall(PetscSynchronizedPrintf(comm, "%" PetscInt_FMT " ", ownedsenders[starts[i] + j])); 11769566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm, "\n")); 117730dcb7c9SBarry Smith } 117830dcb7c9SBarry Smith } 11799566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 118007b52d57SBarry Smith } /* ----------------------------------- */ 118130dcb7c9SBarry Smith 11823677ff5aSBarry Smith /* wait on original sends */ 11833a96401aSBarry Smith if (nsends) { 11849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends, &send_status)); 11859566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nsends, send_waits, send_status)); 11869566063dSJacob Faibussowitsch PetscCall(PetscFree(send_status)); 11873a96401aSBarry Smith } 11889566063dSJacob Faibussowitsch PetscCall(PetscFree(send_waits)); 11899566063dSJacob Faibussowitsch PetscCall(PetscFree(sends)); 11909566063dSJacob Faibussowitsch PetscCall(PetscFree(nprocs)); 11913677ff5aSBarry Smith 11923677ff5aSBarry Smith /* pack messages to send back to local owners */ 119330dcb7c9SBarry Smith starts[0] = 0; 119430dcb7c9SBarry Smith for (i = 1; i < ng; i++) { 119530dcb7c9SBarry Smith if (nownedsenders[i - 1] > 1) starts[i] = starts[i - 1] + nownedsenders[i - 1]; 119630dcb7c9SBarry Smith else starts[i] = starts[i - 1]; 119730dcb7c9SBarry Smith } 119830dcb7c9SBarry Smith nsends2 = nrecvs; 11999566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends2 + 1, &nprocs)); /* length of each message */ 120030dcb7c9SBarry Smith for (i = 0; i < nrecvs; i++) { 120130dcb7c9SBarry Smith nprocs[i] = 1; 120230dcb7c9SBarry Smith for (j = 0; j < len[i]; j++) { 120330dcb7c9SBarry Smith node = recvs[2 * i * nmax + 2 * j] - rstart; 1204f6e5521dSKarl Rupp if (nownedsenders[node] > 1) nprocs[i] += 2 + nownedsenders[node]; 120530dcb7c9SBarry Smith } 120630dcb7c9SBarry Smith } 1207f6e5521dSKarl Rupp nt = 0; 1208f6e5521dSKarl Rupp for (i = 0; i < nsends2; i++) nt += nprocs[i]; 1209f6e5521dSKarl Rupp 12109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nt + 1, &sends2)); 12119566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends2 + 1, &starts2)); 1212f6e5521dSKarl Rupp 1213f6e5521dSKarl Rupp starts2[0] = 0; 1214f6e5521dSKarl Rupp for (i = 1; i < nsends2; i++) starts2[i] = starts2[i - 1] + nprocs[i - 1]; 121530dcb7c9SBarry Smith /* 121630dcb7c9SBarry Smith Each message is 1 + nprocs[i] long, and consists of 121730dcb7c9SBarry Smith (0) the number of nodes being sent back 121830dcb7c9SBarry Smith (1) the local node number, 121930dcb7c9SBarry Smith (2) the number of processors sharing it, 122030dcb7c9SBarry Smith (3) the processors sharing it 122130dcb7c9SBarry Smith */ 122230dcb7c9SBarry Smith for (i = 0; i < nsends2; i++) { 122330dcb7c9SBarry Smith cnt = 1; 122430dcb7c9SBarry Smith sends2[starts2[i]] = 0; 122530dcb7c9SBarry Smith for (j = 0; j < len[i]; j++) { 122630dcb7c9SBarry Smith node = recvs[2 * i * nmax + 2 * j] - rstart; 122730dcb7c9SBarry Smith if (nownedsenders[node] > 1) { 122830dcb7c9SBarry Smith sends2[starts2[i]]++; 122930dcb7c9SBarry Smith sends2[starts2[i] + cnt++] = recvs[2 * i * nmax + 2 * j + 1]; 123030dcb7c9SBarry Smith sends2[starts2[i] + cnt++] = nownedsenders[node]; 12319566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(&sends2[starts2[i] + cnt], &ownedsenders[starts[node]], nownedsenders[node])); 123230dcb7c9SBarry Smith cnt += nownedsenders[node]; 123330dcb7c9SBarry Smith } 123430dcb7c9SBarry Smith } 123530dcb7c9SBarry Smith } 123630dcb7c9SBarry Smith 123730dcb7c9SBarry Smith /* receive the message lengths */ 123830dcb7c9SBarry Smith nrecvs2 = nsends; 12399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2 + 1, &lens2)); 12409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2 + 1, &starts3)); 12419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2 + 1, &recv_waits)); 124248a46eb9SPierre Jolivet for (i = 0; i < nrecvs2; i++) PetscCallMPI(MPI_Irecv(&lens2[i], 1, MPIU_INT, dest[i], tag2, comm, recv_waits + i)); 1243d44834fbSBarry Smith 12448a8e0b3aSBarry Smith /* send the message lengths */ 124548a46eb9SPierre Jolivet for (i = 0; i < nsends2; i++) PetscCallMPI(MPI_Send(&nprocs[i], 1, MPIU_INT, source[i], tag2, comm)); 12468a8e0b3aSBarry Smith 1247d44834fbSBarry Smith /* wait on receives of lens */ 12480c468ba9SBarry Smith if (nrecvs2) { 12499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2, &recv_statuses)); 12509566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrecvs2, recv_waits, recv_statuses)); 12519566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_statuses)); 12520c468ba9SBarry Smith } 12539566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_waits)); 1254d44834fbSBarry Smith 125530dcb7c9SBarry Smith starts3[0] = 0; 1256d44834fbSBarry Smith nt = 0; 125730dcb7c9SBarry Smith for (i = 0; i < nrecvs2 - 1; i++) { 125830dcb7c9SBarry Smith starts3[i + 1] = starts3[i] + lens2[i]; 1259d44834fbSBarry Smith nt += lens2[i]; 126030dcb7c9SBarry Smith } 126176466f69SStefano Zampini if (nrecvs2) nt += lens2[nrecvs2 - 1]; 1262d44834fbSBarry Smith 12639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nt + 1, &recvs2)); 12649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2 + 1, &recv_waits)); 126548a46eb9SPierre Jolivet for (i = 0; i < nrecvs2; i++) PetscCallMPI(MPI_Irecv(recvs2 + starts3[i], lens2[i], MPIU_INT, dest[i], tag3, comm, recv_waits + i)); 126630dcb7c9SBarry Smith 126730dcb7c9SBarry Smith /* send the messages */ 12689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends2 + 1, &send_waits)); 126948a46eb9SPierre Jolivet for (i = 0; i < nsends2; i++) PetscCallMPI(MPI_Isend(sends2 + starts2[i], nprocs[i], MPIU_INT, source[i], tag3, comm, send_waits + i)); 127030dcb7c9SBarry Smith 127130dcb7c9SBarry Smith /* wait on receives */ 12720c468ba9SBarry Smith if (nrecvs2) { 12739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2, &recv_statuses)); 12749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrecvs2, recv_waits, recv_statuses)); 12759566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_statuses)); 12760c468ba9SBarry Smith } 12779566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_waits)); 12789566063dSJacob Faibussowitsch PetscCall(PetscFree(nprocs)); 127930dcb7c9SBarry Smith 128007b52d57SBarry Smith if (debug) { /* ----------------------------------- */ 128130dcb7c9SBarry Smith cnt = 0; 128230dcb7c9SBarry Smith for (i = 0; i < nrecvs2; i++) { 128330dcb7c9SBarry Smith nt = recvs2[cnt++]; 128430dcb7c9SBarry Smith for (j = 0; j < nt; j++) { 12859566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm, "[%d] local node %" PetscInt_FMT " number of subdomains %" PetscInt_FMT ": ", rank, recvs2[cnt], recvs2[cnt + 1])); 128648a46eb9SPierre Jolivet for (k = 0; k < recvs2[cnt + 1]; k++) PetscCall(PetscSynchronizedPrintf(comm, "%" PetscInt_FMT " ", recvs2[cnt + 2 + k])); 128730dcb7c9SBarry Smith cnt += 2 + recvs2[cnt + 1]; 12889566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm, "\n")); 128930dcb7c9SBarry Smith } 129030dcb7c9SBarry Smith } 12919566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 129207b52d57SBarry Smith } /* ----------------------------------- */ 129330dcb7c9SBarry Smith 129430dcb7c9SBarry Smith /* count number subdomains for each local node */ 12959566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(size, &nprocs)); 129630dcb7c9SBarry Smith cnt = 0; 129730dcb7c9SBarry Smith for (i = 0; i < nrecvs2; i++) { 129830dcb7c9SBarry Smith nt = recvs2[cnt++]; 129930dcb7c9SBarry Smith for (j = 0; j < nt; j++) { 1300f6e5521dSKarl Rupp for (k = 0; k < recvs2[cnt + 1]; k++) nprocs[recvs2[cnt + 2 + k]]++; 130130dcb7c9SBarry Smith cnt += 2 + recvs2[cnt + 1]; 130230dcb7c9SBarry Smith } 130330dcb7c9SBarry Smith } 13049371c9d4SSatish Balay nt = 0; 13059371c9d4SSatish Balay for (i = 0; i < size; i++) nt += (nprocs[i] > 0); 130630dcb7c9SBarry Smith *nproc = nt; 13079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nt + 1, procs)); 13089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nt + 1, numprocs)); 13099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nt + 1, indices)); 13100298fd71SBarry Smith for (i = 0; i < nt + 1; i++) (*indices)[i] = NULL; 13119566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &bprocs)); 131230dcb7c9SBarry Smith cnt = 0; 131330dcb7c9SBarry Smith for (i = 0; i < size; i++) { 131430dcb7c9SBarry Smith if (nprocs[i] > 0) { 131530dcb7c9SBarry Smith bprocs[i] = cnt; 131630dcb7c9SBarry Smith (*procs)[cnt] = i; 131730dcb7c9SBarry Smith (*numprocs)[cnt] = nprocs[i]; 13189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nprocs[i], &(*indices)[cnt])); 131930dcb7c9SBarry Smith cnt++; 132030dcb7c9SBarry Smith } 132130dcb7c9SBarry Smith } 132230dcb7c9SBarry Smith 132330dcb7c9SBarry Smith /* make the list of subdomains for each nontrivial local node */ 13249566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(*numprocs, nt)); 132530dcb7c9SBarry Smith cnt = 0; 132630dcb7c9SBarry Smith for (i = 0; i < nrecvs2; i++) { 132730dcb7c9SBarry Smith nt = recvs2[cnt++]; 132830dcb7c9SBarry Smith for (j = 0; j < nt; j++) { 1329f6e5521dSKarl Rupp for (k = 0; k < recvs2[cnt + 1]; k++) (*indices)[bprocs[recvs2[cnt + 2 + k]]][(*numprocs)[bprocs[recvs2[cnt + 2 + k]]]++] = recvs2[cnt]; 133030dcb7c9SBarry Smith cnt += 2 + recvs2[cnt + 1]; 133130dcb7c9SBarry Smith } 133230dcb7c9SBarry Smith } 13339566063dSJacob Faibussowitsch PetscCall(PetscFree(bprocs)); 13349566063dSJacob Faibussowitsch PetscCall(PetscFree(recvs2)); 133530dcb7c9SBarry Smith 133607b52d57SBarry Smith /* sort the node indexing by their global numbers */ 133707b52d57SBarry Smith nt = *nproc; 133807b52d57SBarry Smith for (i = 0; i < nt; i++) { 13399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((*numprocs)[i], &tmp)); 1340f6e5521dSKarl Rupp for (j = 0; j < (*numprocs)[i]; j++) tmp[j] = lindices[(*indices)[i][j]]; 13419566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithArray((*numprocs)[i], tmp, (*indices)[i])); 13429566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp)); 134307b52d57SBarry Smith } 134407b52d57SBarry Smith 134507b52d57SBarry Smith if (debug) { /* ----------------------------------- */ 134630dcb7c9SBarry Smith nt = *nproc; 134730dcb7c9SBarry Smith for (i = 0; i < nt; i++) { 13489566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm, "[%d] subdomain %" PetscInt_FMT " number of indices %" PetscInt_FMT ": ", rank, (*procs)[i], (*numprocs)[i])); 134948a46eb9SPierre Jolivet for (j = 0; j < (*numprocs)[i]; j++) PetscCall(PetscSynchronizedPrintf(comm, "%" PetscInt_FMT " ", (*indices)[i][j])); 13509566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm, "\n")); 135130dcb7c9SBarry Smith } 13529566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 135307b52d57SBarry Smith } /* ----------------------------------- */ 135430dcb7c9SBarry Smith 135530dcb7c9SBarry Smith /* wait on sends */ 135630dcb7c9SBarry Smith if (nsends2) { 13579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends2, &send_status)); 13589566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nsends2, send_waits, send_status)); 13599566063dSJacob Faibussowitsch PetscCall(PetscFree(send_status)); 136030dcb7c9SBarry Smith } 136130dcb7c9SBarry Smith 13629566063dSJacob Faibussowitsch PetscCall(PetscFree(starts3)); 13639566063dSJacob Faibussowitsch PetscCall(PetscFree(dest)); 13649566063dSJacob Faibussowitsch PetscCall(PetscFree(send_waits)); 13653677ff5aSBarry Smith 13669566063dSJacob Faibussowitsch PetscCall(PetscFree(nownedsenders)); 13679566063dSJacob Faibussowitsch PetscCall(PetscFree(ownedsenders)); 13689566063dSJacob Faibussowitsch PetscCall(PetscFree(starts)); 13699566063dSJacob Faibussowitsch PetscCall(PetscFree(starts2)); 13709566063dSJacob Faibussowitsch PetscCall(PetscFree(lens2)); 137189d82c54SBarry Smith 13729566063dSJacob Faibussowitsch PetscCall(PetscFree(source)); 13739566063dSJacob Faibussowitsch PetscCall(PetscFree(len)); 13749566063dSJacob Faibussowitsch PetscCall(PetscFree(recvs)); 13759566063dSJacob Faibussowitsch PetscCall(PetscFree(nprocs)); 13769566063dSJacob Faibussowitsch PetscCall(PetscFree(sends2)); 137724cf384cSBarry Smith 137824cf384cSBarry Smith /* put the information about myself as the first entry in the list */ 137924cf384cSBarry Smith first_procs = (*procs)[0]; 138024cf384cSBarry Smith first_numprocs = (*numprocs)[0]; 138124cf384cSBarry Smith first_indices = (*indices)[0]; 138224cf384cSBarry Smith for (i = 0; i < *nproc; i++) { 138324cf384cSBarry Smith if ((*procs)[i] == rank) { 138424cf384cSBarry Smith (*procs)[0] = (*procs)[i]; 138524cf384cSBarry Smith (*numprocs)[0] = (*numprocs)[i]; 138624cf384cSBarry Smith (*indices)[0] = (*indices)[i]; 138724cf384cSBarry Smith (*procs)[i] = first_procs; 138824cf384cSBarry Smith (*numprocs)[i] = first_numprocs; 138924cf384cSBarry Smith (*indices)[i] = first_indices; 139024cf384cSBarry Smith break; 139124cf384cSBarry Smith } 139224cf384cSBarry Smith } 1393268a049cSStefano Zampini 1394268a049cSStefano Zampini /* save info for reuse */ 1395268a049cSStefano Zampini mapping->info_nproc = *nproc; 1396268a049cSStefano Zampini mapping->info_procs = *procs; 1397268a049cSStefano Zampini mapping->info_numprocs = *numprocs; 1398268a049cSStefano Zampini mapping->info_indices = *indices; 1399268a049cSStefano Zampini mapping->info_cached = PETSC_TRUE; 140089d82c54SBarry Smith PetscFunctionReturn(0); 140189d82c54SBarry Smith } 140289d82c54SBarry Smith 14036a818285SBarry Smith /*@C 14046a818285SBarry Smith ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by ISLocalToGlobalMappingGetBlockInfo() 14056a818285SBarry Smith 14066a818285SBarry Smith Collective on ISLocalToGlobalMapping 14076a818285SBarry Smith 1408f899ff85SJose E. Roman Input Parameter: 14096a818285SBarry Smith . mapping - the mapping from local to global indexing 14106a818285SBarry Smith 1411d8d19677SJose E. Roman Output Parameters: 14126a818285SBarry Smith + nproc - number of processors that are connected to this one 14136a818285SBarry Smith . proc - neighboring processors 14146a818285SBarry Smith . numproc - number of indices for each processor 14156a818285SBarry Smith - indices - indices of local nodes shared with neighbor (sorted by global numbering) 14166a818285SBarry Smith 14176a818285SBarry Smith Level: advanced 14186a818285SBarry Smith 1419db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 1420db781477SPatrick Sanan `ISLocalToGlobalMappingGetInfo()` 14216a818285SBarry Smith @*/ 1422*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[]) 1423*d71ae5a4SJacob Faibussowitsch { 14246a818285SBarry Smith PetscFunctionBegin; 1425cbc1caf0SMatthew G. Knepley PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 1426268a049cSStefano Zampini if (mapping->info_free) { 14279566063dSJacob Faibussowitsch PetscCall(PetscFree(*numprocs)); 14286a818285SBarry Smith if (*indices) { 1429268a049cSStefano Zampini PetscInt i; 1430268a049cSStefano Zampini 14319566063dSJacob Faibussowitsch PetscCall(PetscFree((*indices)[0])); 143248a46eb9SPierre Jolivet for (i = 1; i < *nproc; i++) PetscCall(PetscFree((*indices)[i])); 14339566063dSJacob Faibussowitsch PetscCall(PetscFree(*indices)); 14346a818285SBarry Smith } 1435268a049cSStefano Zampini } 1436268a049cSStefano Zampini *nproc = 0; 1437268a049cSStefano Zampini *procs = NULL; 1438268a049cSStefano Zampini *numprocs = NULL; 1439268a049cSStefano Zampini *indices = NULL; 14406a818285SBarry Smith PetscFunctionReturn(0); 14416a818285SBarry Smith } 14426a818285SBarry Smith 14436a818285SBarry Smith /*@C 14446a818285SBarry Smith ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each processor and 14456a818285SBarry Smith each index shared by more than one processor 14466a818285SBarry Smith 14476a818285SBarry Smith Collective on ISLocalToGlobalMapping 14486a818285SBarry Smith 1449f899ff85SJose E. Roman Input Parameter: 14506a818285SBarry Smith . mapping - the mapping from local to global indexing 14516a818285SBarry Smith 1452d8d19677SJose E. Roman Output Parameters: 14536a818285SBarry Smith + nproc - number of processors that are connected to this one 14546a818285SBarry Smith . proc - neighboring processors 14556a818285SBarry Smith . numproc - number of indices for each subdomain (processor) 14566a818285SBarry Smith - indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering) 14576a818285SBarry Smith 14586a818285SBarry Smith Level: advanced 14596a818285SBarry Smith 14601bd0b88eSStefano Zampini Notes: The user needs to call ISLocalToGlobalMappingRestoreInfo when the data is no longer needed. 14611bd0b88eSStefano Zampini 14626a818285SBarry Smith Fortran Usage: 14636a818285SBarry Smith $ ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by 14646a818285SBarry Smith $ ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc], 14656a818285SBarry Smith PetscInt indices[nproc][numprocmax],ierr) 14666a818285SBarry Smith There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and 14676a818285SBarry Smith indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough. 14686a818285SBarry Smith 1469db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 1470db781477SPatrick Sanan `ISLocalToGlobalMappingRestoreInfo()` 14716a818285SBarry Smith @*/ 1472*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[]) 1473*d71ae5a4SJacob Faibussowitsch { 14748a1f772fSStefano Zampini PetscInt **bindices = NULL, *bnumprocs = NULL, bs, i, j, k; 14756a818285SBarry Smith 14766a818285SBarry Smith PetscFunctionBegin; 14776a818285SBarry Smith PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 14788a1f772fSStefano Zampini bs = mapping->bs; 14799566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockInfo(mapping, nproc, procs, &bnumprocs, &bindices)); 1480268a049cSStefano Zampini if (bs > 1) { /* we need to expand the cached info */ 14819566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(*nproc, &*indices)); 14829566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(*nproc, &*numprocs)); 14836a818285SBarry Smith for (i = 0; i < *nproc; i++) { 14849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * bnumprocs[i], &(*indices)[i])); 1485268a049cSStefano Zampini for (j = 0; j < bnumprocs[i]; j++) { 1486ad540459SPierre Jolivet for (k = 0; k < bs; k++) (*indices)[i][j * bs + k] = bs * bindices[i][j] + k; 14876a818285SBarry Smith } 1488268a049cSStefano Zampini (*numprocs)[i] = bnumprocs[i] * bs; 14896a818285SBarry Smith } 1490268a049cSStefano Zampini mapping->info_free = PETSC_TRUE; 1491268a049cSStefano Zampini } else { 1492268a049cSStefano Zampini *numprocs = bnumprocs; 1493268a049cSStefano Zampini *indices = bindices; 14946a818285SBarry Smith } 14956a818285SBarry Smith PetscFunctionReturn(0); 14966a818285SBarry Smith } 14976a818285SBarry Smith 149807b52d57SBarry Smith /*@C 149907b52d57SBarry Smith ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by ISLocalToGlobalMappingGetInfo() 150089d82c54SBarry Smith 150107b52d57SBarry Smith Collective on ISLocalToGlobalMapping 150207b52d57SBarry Smith 1503f899ff85SJose E. Roman Input Parameter: 150407b52d57SBarry Smith . mapping - the mapping from local to global indexing 150507b52d57SBarry Smith 1506d8d19677SJose E. Roman Output Parameters: 150707b52d57SBarry Smith + nproc - number of processors that are connected to this one 150807b52d57SBarry Smith . proc - neighboring processors 150907b52d57SBarry Smith . numproc - number of indices for each processor 151007b52d57SBarry Smith - indices - indices of local nodes shared with neighbor (sorted by global numbering) 151107b52d57SBarry Smith 151207b52d57SBarry Smith Level: advanced 151307b52d57SBarry Smith 1514db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 1515db781477SPatrick Sanan `ISLocalToGlobalMappingGetInfo()` 151607b52d57SBarry Smith @*/ 1517*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[]) 1518*d71ae5a4SJacob Faibussowitsch { 151907b52d57SBarry Smith PetscFunctionBegin; 15209566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockInfo(mapping, nproc, procs, numprocs, indices)); 152107b52d57SBarry Smith PetscFunctionReturn(0); 152207b52d57SBarry Smith } 152386994e45SJed Brown 152486994e45SJed Brown /*@C 15251bd0b88eSStefano Zampini ISLocalToGlobalMappingGetNodeInfo - Gets the neighbor information for each node 15261bd0b88eSStefano Zampini 15271bd0b88eSStefano Zampini Collective on ISLocalToGlobalMapping 15281bd0b88eSStefano Zampini 1529f899ff85SJose E. Roman Input Parameter: 15301bd0b88eSStefano Zampini . mapping - the mapping from local to global indexing 15311bd0b88eSStefano Zampini 1532d8d19677SJose E. Roman Output Parameters: 15331bd0b88eSStefano Zampini + nnodes - number of local nodes (same ISLocalToGlobalMappingGetSize()) 15341bd0b88eSStefano Zampini . count - number of neighboring processors per node 15351bd0b88eSStefano Zampini - indices - indices of processes sharing the node (sorted) 15361bd0b88eSStefano Zampini 15371bd0b88eSStefano Zampini Level: advanced 15381bd0b88eSStefano Zampini 15391bd0b88eSStefano Zampini Notes: The user needs to call ISLocalToGlobalMappingRestoreInfo when the data is no longer needed. 15401bd0b88eSStefano Zampini 1541db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 1542db781477SPatrick Sanan `ISLocalToGlobalMappingGetInfo()`, `ISLocalToGlobalMappingRestoreNodeInfo()` 15431bd0b88eSStefano Zampini @*/ 1544*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetNodeInfo(ISLocalToGlobalMapping mapping, PetscInt *nnodes, PetscInt *count[], PetscInt **indices[]) 1545*d71ae5a4SJacob Faibussowitsch { 15461bd0b88eSStefano Zampini PetscInt n; 15471bd0b88eSStefano Zampini 15481bd0b88eSStefano Zampini PetscFunctionBegin; 15491bd0b88eSStefano Zampini PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 15509566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(mapping, &n)); 15511bd0b88eSStefano Zampini if (!mapping->info_nodec) { 15521bd0b88eSStefano Zampini PetscInt i, m, n_neigh, *neigh, *n_shared, **shared; 15531bd0b88eSStefano Zampini 15549566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n + 1, &mapping->info_nodec, n, &mapping->info_nodei)); 15559566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetInfo(mapping, &n_neigh, &neigh, &n_shared, &shared)); 1556ad540459SPierre Jolivet for (i = 0; i < n; i++) mapping->info_nodec[i] = 1; 1557071fcb05SBarry Smith m = n; 1558071fcb05SBarry Smith mapping->info_nodec[n] = 0; 15591bd0b88eSStefano Zampini for (i = 1; i < n_neigh; i++) { 15601bd0b88eSStefano Zampini PetscInt j; 15611bd0b88eSStefano Zampini 15621bd0b88eSStefano Zampini m += n_shared[i]; 15631bd0b88eSStefano Zampini for (j = 0; j < n_shared[i]; j++) mapping->info_nodec[shared[i][j]] += 1; 15641bd0b88eSStefano Zampini } 15659566063dSJacob Faibussowitsch if (n) PetscCall(PetscMalloc1(m, &mapping->info_nodei[0])); 15661bd0b88eSStefano Zampini for (i = 1; i < n; i++) mapping->info_nodei[i] = mapping->info_nodei[i - 1] + mapping->info_nodec[i - 1]; 15679566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(mapping->info_nodec, n)); 15689371c9d4SSatish Balay for (i = 0; i < n; i++) { 15699371c9d4SSatish Balay mapping->info_nodec[i] = 1; 15709371c9d4SSatish Balay mapping->info_nodei[i][0] = neigh[0]; 15719371c9d4SSatish Balay } 15721bd0b88eSStefano Zampini for (i = 1; i < n_neigh; i++) { 15731bd0b88eSStefano Zampini PetscInt j; 15741bd0b88eSStefano Zampini 15751bd0b88eSStefano Zampini for (j = 0; j < n_shared[i]; j++) { 15761bd0b88eSStefano Zampini PetscInt k = shared[i][j]; 15771bd0b88eSStefano Zampini 15781bd0b88eSStefano Zampini mapping->info_nodei[k][mapping->info_nodec[k]] = neigh[i]; 15791bd0b88eSStefano Zampini mapping->info_nodec[k] += 1; 15801bd0b88eSStefano Zampini } 15811bd0b88eSStefano Zampini } 15829566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(PetscSortRemoveDupsInt(&mapping->info_nodec[i], mapping->info_nodei[i])); 15839566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreInfo(mapping, &n_neigh, &neigh, &n_shared, &shared)); 15841bd0b88eSStefano Zampini } 15851bd0b88eSStefano Zampini if (nnodes) *nnodes = n; 15861bd0b88eSStefano Zampini if (count) *count = mapping->info_nodec; 15871bd0b88eSStefano Zampini if (indices) *indices = mapping->info_nodei; 15881bd0b88eSStefano Zampini PetscFunctionReturn(0); 15891bd0b88eSStefano Zampini } 15901bd0b88eSStefano Zampini 15911bd0b88eSStefano Zampini /*@C 15921bd0b88eSStefano Zampini ISLocalToGlobalMappingRestoreNodeInfo - Frees the memory allocated by ISLocalToGlobalMappingGetNodeInfo() 15931bd0b88eSStefano Zampini 15941bd0b88eSStefano Zampini Collective on ISLocalToGlobalMapping 15951bd0b88eSStefano Zampini 1596f899ff85SJose E. Roman Input Parameter: 15971bd0b88eSStefano Zampini . mapping - the mapping from local to global indexing 15981bd0b88eSStefano Zampini 1599d8d19677SJose E. Roman Output Parameters: 16001bd0b88eSStefano Zampini + nnodes - number of local nodes 16011bd0b88eSStefano Zampini . count - number of neighboring processors per node 16021bd0b88eSStefano Zampini - indices - indices of processes sharing the node (sorted) 16031bd0b88eSStefano Zampini 16041bd0b88eSStefano Zampini Level: advanced 16051bd0b88eSStefano Zampini 1606db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 1607db781477SPatrick Sanan `ISLocalToGlobalMappingGetInfo()` 16081bd0b88eSStefano Zampini @*/ 1609*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreNodeInfo(ISLocalToGlobalMapping mapping, PetscInt *nnodes, PetscInt *count[], PetscInt **indices[]) 1610*d71ae5a4SJacob Faibussowitsch { 16111bd0b88eSStefano Zampini PetscFunctionBegin; 16121bd0b88eSStefano Zampini PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 16131bd0b88eSStefano Zampini if (nnodes) *nnodes = 0; 16141bd0b88eSStefano Zampini if (count) *count = NULL; 16151bd0b88eSStefano Zampini if (indices) *indices = NULL; 16161bd0b88eSStefano Zampini PetscFunctionReturn(0); 16171bd0b88eSStefano Zampini } 16181bd0b88eSStefano Zampini 16191bd0b88eSStefano Zampini /*@C 1620107e9a97SBarry Smith ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped 162186994e45SJed Brown 162286994e45SJed Brown Not Collective 162386994e45SJed Brown 16244165533cSJose E. Roman Input Parameter: 162586994e45SJed Brown . ltog - local to global mapping 162686994e45SJed Brown 16274165533cSJose E. Roman Output Parameter: 1628565245c5SBarry Smith . array - array of indices, the length of this array may be obtained with ISLocalToGlobalMappingGetSize() 162986994e45SJed Brown 163086994e45SJed Brown Level: advanced 163186994e45SJed Brown 163295452b02SPatrick Sanan Notes: 163395452b02SPatrick Sanan ISLocalToGlobalMappingGetSize() returns the length the this array 1634107e9a97SBarry Smith 1635db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingRestoreIndices()`, `ISLocalToGlobalMappingGetBlockIndices()`, `ISLocalToGlobalMappingRestoreBlockIndices()` 163686994e45SJed Brown @*/ 1637*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog, const PetscInt **array) 1638*d71ae5a4SJacob Faibussowitsch { 163986994e45SJed Brown PetscFunctionBegin; 164086994e45SJed Brown PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1); 164186994e45SJed Brown PetscValidPointer(array, 2); 164245b6f7e9SBarry Smith if (ltog->bs == 1) { 164386994e45SJed Brown *array = ltog->indices; 164445b6f7e9SBarry Smith } else { 164545b6f7e9SBarry Smith PetscInt *jj, k, i, j, n = ltog->n, bs = ltog->bs; 164645b6f7e9SBarry Smith const PetscInt *ii; 164745b6f7e9SBarry Smith 16489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * n, &jj)); 164945b6f7e9SBarry Smith *array = jj; 165045b6f7e9SBarry Smith k = 0; 165145b6f7e9SBarry Smith ii = ltog->indices; 165245b6f7e9SBarry Smith for (i = 0; i < n; i++) 16539371c9d4SSatish Balay for (j = 0; j < bs; j++) jj[k++] = bs * ii[i] + j; 165445b6f7e9SBarry Smith } 165586994e45SJed Brown PetscFunctionReturn(0); 165686994e45SJed Brown } 165786994e45SJed Brown 165886994e45SJed Brown /*@C 1659193a2b41SJulian Andrej ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with ISLocalToGlobalMappingGetIndices() 166086994e45SJed Brown 166186994e45SJed Brown Not Collective 166286994e45SJed Brown 16634165533cSJose E. Roman Input Parameters: 166486994e45SJed Brown + ltog - local to global mapping 166586994e45SJed Brown - array - array of indices 166686994e45SJed Brown 166786994e45SJed Brown Level: advanced 166886994e45SJed Brown 1669db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingGetIndices()` 167086994e45SJed Brown @*/ 1671*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog, const PetscInt **array) 1672*d71ae5a4SJacob Faibussowitsch { 167386994e45SJed Brown PetscFunctionBegin; 167486994e45SJed Brown PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1); 167586994e45SJed Brown PetscValidPointer(array, 2); 1676c9cc58a2SBarry Smith PetscCheck(ltog->bs != 1 || *array == ltog->indices, PETSC_COMM_SELF, PETSC_ERR_ARG_BADPTR, "Trying to return mismatched pointer"); 167745b6f7e9SBarry Smith 167848a46eb9SPierre Jolivet if (ltog->bs > 1) PetscCall(PetscFree(*(void **)array)); 167945b6f7e9SBarry Smith PetscFunctionReturn(0); 168045b6f7e9SBarry Smith } 168145b6f7e9SBarry Smith 168245b6f7e9SBarry Smith /*@C 168345b6f7e9SBarry Smith ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block 168445b6f7e9SBarry Smith 168545b6f7e9SBarry Smith Not Collective 168645b6f7e9SBarry Smith 16874165533cSJose E. Roman Input Parameter: 168845b6f7e9SBarry Smith . ltog - local to global mapping 168945b6f7e9SBarry Smith 16904165533cSJose E. Roman Output Parameter: 169145b6f7e9SBarry Smith . array - array of indices 169245b6f7e9SBarry Smith 169345b6f7e9SBarry Smith Level: advanced 169445b6f7e9SBarry Smith 1695db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingRestoreBlockIndices()` 169645b6f7e9SBarry Smith @*/ 1697*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog, const PetscInt **array) 1698*d71ae5a4SJacob Faibussowitsch { 169945b6f7e9SBarry Smith PetscFunctionBegin; 170045b6f7e9SBarry Smith PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1); 170145b6f7e9SBarry Smith PetscValidPointer(array, 2); 170245b6f7e9SBarry Smith *array = ltog->indices; 170345b6f7e9SBarry Smith PetscFunctionReturn(0); 170445b6f7e9SBarry Smith } 170545b6f7e9SBarry Smith 170645b6f7e9SBarry Smith /*@C 170745b6f7e9SBarry Smith ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with ISLocalToGlobalMappingGetBlockIndices() 170845b6f7e9SBarry Smith 170945b6f7e9SBarry Smith Not Collective 171045b6f7e9SBarry Smith 17114165533cSJose E. Roman Input Parameters: 171245b6f7e9SBarry Smith + ltog - local to global mapping 171345b6f7e9SBarry Smith - array - array of indices 171445b6f7e9SBarry Smith 171545b6f7e9SBarry Smith Level: advanced 171645b6f7e9SBarry Smith 1717db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingGetIndices()` 171845b6f7e9SBarry Smith @*/ 1719*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog, const PetscInt **array) 1720*d71ae5a4SJacob Faibussowitsch { 172145b6f7e9SBarry Smith PetscFunctionBegin; 172245b6f7e9SBarry Smith PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1); 172345b6f7e9SBarry Smith PetscValidPointer(array, 2); 172408401ef6SPierre Jolivet PetscCheck(*array == ltog->indices, PETSC_COMM_SELF, PETSC_ERR_ARG_BADPTR, "Trying to return mismatched pointer"); 17250298fd71SBarry Smith *array = NULL; 172686994e45SJed Brown PetscFunctionReturn(0); 172786994e45SJed Brown } 1728f7efa3c7SJed Brown 1729f7efa3c7SJed Brown /*@C 1730f7efa3c7SJed Brown ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings 1731f7efa3c7SJed Brown 1732f7efa3c7SJed Brown Not Collective 1733f7efa3c7SJed Brown 17344165533cSJose E. Roman Input Parameters: 1735f7efa3c7SJed Brown + comm - communicator for the new mapping, must contain the communicator of every mapping to concatenate 1736f7efa3c7SJed Brown . n - number of mappings to concatenate 1737f7efa3c7SJed Brown - ltogs - local to global mappings 1738f7efa3c7SJed Brown 17394165533cSJose E. Roman Output Parameter: 1740f7efa3c7SJed Brown . ltogcat - new mapping 1741f7efa3c7SJed Brown 17429d90f715SBarry Smith Note: this currently always returns a mapping with block size of 1 17439d90f715SBarry Smith 17449d90f715SBarry Smith Developer Note: If all the input mapping have the same block size we could easily handle that as a special case 17459d90f715SBarry Smith 1746f7efa3c7SJed Brown Level: advanced 1747f7efa3c7SJed Brown 1748db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()` 1749f7efa3c7SJed Brown @*/ 1750*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm, PetscInt n, const ISLocalToGlobalMapping ltogs[], ISLocalToGlobalMapping *ltogcat) 1751*d71ae5a4SJacob Faibussowitsch { 1752f7efa3c7SJed Brown PetscInt i, cnt, m, *idx; 1753f7efa3c7SJed Brown 1754f7efa3c7SJed Brown PetscFunctionBegin; 175508401ef6SPierre Jolivet PetscCheck(n >= 0, comm, PETSC_ERR_ARG_OUTOFRANGE, "Must have a non-negative number of mappings, given %" PetscInt_FMT, n); 1756f7efa3c7SJed Brown if (n > 0) PetscValidPointer(ltogs, 3); 1757f7efa3c7SJed Brown for (i = 0; i < n; i++) PetscValidHeaderSpecific(ltogs[i], IS_LTOGM_CLASSID, 3); 1758f7efa3c7SJed Brown PetscValidPointer(ltogcat, 4); 1759f7efa3c7SJed Brown for (cnt = 0, i = 0; i < n; i++) { 17609566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(ltogs[i], &m)); 1761f7efa3c7SJed Brown cnt += m; 1762f7efa3c7SJed Brown } 17639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(cnt, &idx)); 1764f7efa3c7SJed Brown for (cnt = 0, i = 0; i < n; i++) { 1765f7efa3c7SJed Brown const PetscInt *subidx; 17669566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(ltogs[i], &m)); 17679566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(ltogs[i], &subidx)); 17689566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(&idx[cnt], subidx, m)); 17699566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(ltogs[i], &subidx)); 1770f7efa3c7SJed Brown cnt += m; 1771f7efa3c7SJed Brown } 17729566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm, 1, cnt, idx, PETSC_OWN_POINTER, ltogcat)); 1773f7efa3c7SJed Brown PetscFunctionReturn(0); 1774f7efa3c7SJed Brown } 177504a59952SBarry Smith 1776413f72f0SBarry Smith /*MC 1777413f72f0SBarry Smith ISLOCALTOGLOBALMAPPINGBASIC - basic implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is 1778413f72f0SBarry Smith used this is good for only small and moderate size problems. 1779413f72f0SBarry Smith 1780413f72f0SBarry Smith Options Database Keys: 1781a2b725a8SWilliam Gropp . -islocaltoglobalmapping_type basic - select this method 1782413f72f0SBarry Smith 1783413f72f0SBarry Smith Level: beginner 1784413f72f0SBarry Smith 1785db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`, `ISLOCALTOGLOBALMAPPINGHASH` 1786413f72f0SBarry Smith M*/ 1787*d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Basic(ISLocalToGlobalMapping ltog) 1788*d71ae5a4SJacob Faibussowitsch { 1789413f72f0SBarry Smith PetscFunctionBegin; 1790413f72f0SBarry Smith ltog->ops->globaltolocalmappingapply = ISGlobalToLocalMappingApply_Basic; 1791413f72f0SBarry Smith ltog->ops->globaltolocalmappingsetup = ISGlobalToLocalMappingSetUp_Basic; 1792413f72f0SBarry Smith ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Basic; 1793413f72f0SBarry Smith ltog->ops->destroy = ISLocalToGlobalMappingDestroy_Basic; 1794413f72f0SBarry Smith PetscFunctionReturn(0); 1795413f72f0SBarry Smith } 1796413f72f0SBarry Smith 1797413f72f0SBarry Smith /*MC 1798413f72f0SBarry Smith ISLOCALTOGLOBALMAPPINGHASH - hash implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is 1799ed56e8eaSBarry Smith used this is good for large memory problems. 1800413f72f0SBarry Smith 1801413f72f0SBarry Smith Options Database Keys: 1802a2b725a8SWilliam Gropp . -islocaltoglobalmapping_type hash - select this method 1803413f72f0SBarry Smith 180495452b02SPatrick Sanan Notes: 180595452b02SPatrick Sanan This is selected automatically for large problems if the user does not set the type. 1806ed56e8eaSBarry Smith 1807413f72f0SBarry Smith Level: beginner 1808413f72f0SBarry Smith 1809db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`, `ISLOCALTOGLOBALMAPPINGHASH` 1810413f72f0SBarry Smith M*/ 1811*d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Hash(ISLocalToGlobalMapping ltog) 1812*d71ae5a4SJacob Faibussowitsch { 1813413f72f0SBarry Smith PetscFunctionBegin; 1814413f72f0SBarry Smith ltog->ops->globaltolocalmappingapply = ISGlobalToLocalMappingApply_Hash; 1815413f72f0SBarry Smith ltog->ops->globaltolocalmappingsetup = ISGlobalToLocalMappingSetUp_Hash; 1816413f72f0SBarry Smith ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Hash; 1817413f72f0SBarry Smith ltog->ops->destroy = ISLocalToGlobalMappingDestroy_Hash; 1818413f72f0SBarry Smith PetscFunctionReturn(0); 1819413f72f0SBarry Smith } 1820413f72f0SBarry Smith 1821413f72f0SBarry Smith /*@C 1822413f72f0SBarry Smith ISLocalToGlobalMappingRegister - Adds a method for applying a global to local mapping with an ISLocalToGlobalMapping 1823413f72f0SBarry Smith 1824413f72f0SBarry Smith Not Collective 1825413f72f0SBarry Smith 1826413f72f0SBarry Smith Input Parameters: 1827413f72f0SBarry Smith + sname - name of a new method 1828413f72f0SBarry Smith - routine_create - routine to create method context 1829413f72f0SBarry Smith 1830413f72f0SBarry Smith Notes: 1831ed56e8eaSBarry Smith ISLocalToGlobalMappingRegister() may be called multiple times to add several user-defined mappings. 1832413f72f0SBarry Smith 1833413f72f0SBarry Smith Sample usage: 1834413f72f0SBarry Smith .vb 1835413f72f0SBarry Smith ISLocalToGlobalMappingRegister("my_mapper",MyCreate); 1836413f72f0SBarry Smith .ve 1837413f72f0SBarry Smith 1838ed56e8eaSBarry Smith Then, your mapping can be chosen with the procedural interface via 1839413f72f0SBarry Smith $ ISLocalToGlobalMappingSetType(ltog,"my_mapper") 1840413f72f0SBarry Smith or at runtime via the option 1841ed56e8eaSBarry Smith $ -islocaltoglobalmapping_type my_mapper 1842413f72f0SBarry Smith 1843413f72f0SBarry Smith Level: advanced 1844413f72f0SBarry Smith 1845db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingRegisterAll()`, `ISLocalToGlobalMappingRegisterDestroy()`, `ISLOCALTOGLOBALMAPPINGBASIC`, `ISLOCALTOGLOBALMAPPINGHASH` 1846413f72f0SBarry Smith 1847413f72f0SBarry Smith @*/ 1848*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRegister(const char sname[], PetscErrorCode (*function)(ISLocalToGlobalMapping)) 1849*d71ae5a4SJacob Faibussowitsch { 1850413f72f0SBarry Smith PetscFunctionBegin; 18519566063dSJacob Faibussowitsch PetscCall(ISInitializePackage()); 18529566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&ISLocalToGlobalMappingList, sname, function)); 1853413f72f0SBarry Smith PetscFunctionReturn(0); 1854413f72f0SBarry Smith } 1855413f72f0SBarry Smith 1856413f72f0SBarry Smith /*@C 1857ed56e8eaSBarry Smith ISLocalToGlobalMappingSetType - Builds ISLocalToGlobalMapping for a particular global to local mapping approach. 1858413f72f0SBarry Smith 1859413f72f0SBarry Smith Logically Collective on ISLocalToGlobalMapping 1860413f72f0SBarry Smith 1861413f72f0SBarry Smith Input Parameters: 1862413f72f0SBarry Smith + ltog - the ISLocalToGlobalMapping object 1863413f72f0SBarry Smith - type - a known method 1864413f72f0SBarry Smith 1865413f72f0SBarry Smith Options Database Key: 1866ed56e8eaSBarry Smith . -islocaltoglobalmapping_type <method> - Sets the method; use -help for a list 1867413f72f0SBarry Smith of available methods (for instance, basic or hash) 1868413f72f0SBarry Smith 1869413f72f0SBarry Smith Notes: 1870413f72f0SBarry Smith See "petsc/include/petscis.h" for available methods 1871413f72f0SBarry Smith 1872413f72f0SBarry Smith Normally, it is best to use the ISLocalToGlobalMappingSetFromOptions() command and 1873413f72f0SBarry Smith then set the ISLocalToGlobalMapping type from the options database rather than by using 1874413f72f0SBarry Smith this routine. 1875413f72f0SBarry Smith 1876413f72f0SBarry Smith Level: intermediate 1877413f72f0SBarry Smith 1878413f72f0SBarry Smith Developer Note: ISLocalToGlobalMappingRegister() is used to add new types to ISLocalToGlobalMappingList from which they 1879413f72f0SBarry Smith are accessed by ISLocalToGlobalMappingSetType(). 1880413f72f0SBarry Smith 1881db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingRegister()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingGetType()` 1882413f72f0SBarry Smith @*/ 1883*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType type) 1884*d71ae5a4SJacob Faibussowitsch { 1885413f72f0SBarry Smith PetscBool match; 18865f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(ISLocalToGlobalMapping) = NULL; 1887413f72f0SBarry Smith 1888413f72f0SBarry Smith PetscFunctionBegin; 1889413f72f0SBarry Smith PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1); 1890a0d79125SStefano Zampini if (type) PetscValidCharPointer(type, 2); 1891413f72f0SBarry Smith 18929566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)ltog, type, &match)); 1893413f72f0SBarry Smith if (match) PetscFunctionReturn(0); 1894413f72f0SBarry Smith 1895a0d79125SStefano Zampini /* L2G maps defer type setup at globaltolocal calls, allow passing NULL here */ 1896a0d79125SStefano Zampini if (type) { 18979566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(ISLocalToGlobalMappingList, type, &r)); 1898a0d79125SStefano Zampini PetscCheck(r, PetscObjectComm((PetscObject)ltog), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested ISLocalToGlobalMapping type %s", type); 1899a0d79125SStefano Zampini } 1900413f72f0SBarry Smith /* Destroy the previous private LTOG context */ 1901dbbe0bcdSBarry Smith PetscTryTypeMethod(ltog, destroy); 1902413f72f0SBarry Smith ltog->ops->destroy = NULL; 1903dbbe0bcdSBarry Smith 19049566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)ltog, type)); 19059566063dSJacob Faibussowitsch if (r) PetscCall((*r)(ltog)); 1906a0d79125SStefano Zampini PetscFunctionReturn(0); 1907a0d79125SStefano Zampini } 1908a0d79125SStefano Zampini 1909a0d79125SStefano Zampini /*@C 1910a0d79125SStefano Zampini ISLocalToGlobalMappingGetType - Get the type of the l2g map 1911a0d79125SStefano Zampini 1912a0d79125SStefano Zampini Not Collective 1913a0d79125SStefano Zampini 1914a0d79125SStefano Zampini Input Parameter: 1915a0d79125SStefano Zampini . ltog - the ISLocalToGlobalMapping object 1916a0d79125SStefano Zampini 1917a0d79125SStefano Zampini Output Parameter: 1918a0d79125SStefano Zampini . type - the type 1919a0d79125SStefano Zampini 192049762cbcSSatish Balay Level: intermediate 192149762cbcSSatish Balay 1922db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingRegister()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()` 1923a0d79125SStefano Zampini @*/ 1924*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType *type) 1925*d71ae5a4SJacob Faibussowitsch { 1926a0d79125SStefano Zampini PetscFunctionBegin; 1927a0d79125SStefano Zampini PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1); 1928a0d79125SStefano Zampini PetscValidPointer(type, 2); 1929a0d79125SStefano Zampini *type = ((PetscObject)ltog)->type_name; 1930413f72f0SBarry Smith PetscFunctionReturn(0); 1931413f72f0SBarry Smith } 1932413f72f0SBarry Smith 1933413f72f0SBarry Smith PetscBool ISLocalToGlobalMappingRegisterAllCalled = PETSC_FALSE; 1934413f72f0SBarry Smith 1935413f72f0SBarry Smith /*@C 1936413f72f0SBarry Smith ISLocalToGlobalMappingRegisterAll - Registers all of the local to global mapping components in the IS package. 1937413f72f0SBarry Smith 1938413f72f0SBarry Smith Not Collective 1939413f72f0SBarry Smith 1940413f72f0SBarry Smith Level: advanced 1941413f72f0SBarry Smith 1942db781477SPatrick Sanan .seealso: `ISRegister()`, `ISLocalToGlobalRegister()` 1943413f72f0SBarry Smith @*/ 1944*d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRegisterAll(void) 1945*d71ae5a4SJacob Faibussowitsch { 1946413f72f0SBarry Smith PetscFunctionBegin; 1947413f72f0SBarry Smith if (ISLocalToGlobalMappingRegisterAllCalled) PetscFunctionReturn(0); 1948413f72f0SBarry Smith ISLocalToGlobalMappingRegisterAllCalled = PETSC_TRUE; 19499566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGBASIC, ISLocalToGlobalMappingCreate_Basic)); 19509566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGHASH, ISLocalToGlobalMappingCreate_Hash)); 1951413f72f0SBarry Smith PetscFunctionReturn(0); 1952413f72f0SBarry Smith } 1953