12362add9SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/isimpl.h> /*I "petscis.h" I*/ 3e8f14785SLisandro Dalcin #include <petsc/private/hashmapi.h> 40c312b8eSJed Brown #include <petscsf.h> 5665c2dedSJed Brown #include <petscviewer.h> 62362add9SBarry Smith 77087cfbeSBarry Smith PetscClassId IS_LTOGM_CLASSID; 8268a049cSStefano Zampini static PetscErrorCode ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping, PetscInt *, PetscInt **, PetscInt **, PetscInt ***); 98e58c17dSMatthew Knepley 10413f72f0SBarry Smith typedef struct { 11413f72f0SBarry Smith PetscInt *globals; 12413f72f0SBarry Smith } ISLocalToGlobalMapping_Basic; 13413f72f0SBarry Smith 14413f72f0SBarry Smith typedef struct { 15e8f14785SLisandro Dalcin PetscHMapI globalht; 16413f72f0SBarry Smith } ISLocalToGlobalMapping_Hash; 17413f72f0SBarry Smith 186528b96dSMatthew G. Knepley /*@C 19cab54364SBarry Smith ISGetPointRange - Returns a description of the points in an `IS` suitable for traversal 20413f72f0SBarry Smith 216528b96dSMatthew G. Knepley Not collective 226528b96dSMatthew G. Knepley 236528b96dSMatthew G. Knepley Input Parameter: 24cab54364SBarry Smith . pointIS - The `IS` object 256528b96dSMatthew G. Knepley 266528b96dSMatthew G. Knepley Output Parameters: 276528b96dSMatthew G. Knepley + pStart - The first index, see notes 286528b96dSMatthew G. Knepley . pEnd - One past the last index, see notes 296528b96dSMatthew G. Knepley - points - The indices, see notes 306528b96dSMatthew G. Knepley 316528b96dSMatthew G. Knepley Level: intermediate 326528b96dSMatthew G. Knepley 33cab54364SBarry Smith Notes: 34cab54364SBarry Smith If the `IS` contains contiguous indices in an `ISSTRIDE`, then the indices are contained in [pStart, pEnd) and points = NULL. Otherwise, pStart = 0, pEnd = numIndices, and points is an array of the indices. This supports the following pattern 35cab54364SBarry Smith .vb 36cab54364SBarry Smith ISGetPointRange(is, &pStart, &pEnd, &points); 37cab54364SBarry Smith for (p = pStart; p < pEnd; ++p) { 38cab54364SBarry Smith const PetscInt point = points ? points[p] : p; 39cab54364SBarry Smith } 40cab54364SBarry Smith ISRestorePointRange(is, &pstart, &pEnd, &points); 41cab54364SBarry Smith .ve 42cab54364SBarry Smith 43cab54364SBarry Smith .seealso: [](sec_scatter), `IS`, `ISRestorePointRange()`, `ISGetPointSubrange()`, `ISGetIndices()`, `ISCreateStride()` 446528b96dSMatthew G. Knepley @*/ 45d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGetPointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points) 46d71ae5a4SJacob Faibussowitsch { 479305a4c7SMatthew G. Knepley PetscInt numCells, step = 1; 489305a4c7SMatthew G. Knepley PetscBool isStride; 499305a4c7SMatthew G. Knepley 509305a4c7SMatthew G. Knepley PetscFunctionBeginHot; 519305a4c7SMatthew G. Knepley *pStart = 0; 529305a4c7SMatthew G. Knepley *points = NULL; 539566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(pointIS, &numCells)); 549566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pointIS, ISSTRIDE, &isStride)); 559566063dSJacob Faibussowitsch if (isStride) PetscCall(ISStrideGetInfo(pointIS, pStart, &step)); 569305a4c7SMatthew G. Knepley *pEnd = *pStart + numCells; 579566063dSJacob Faibussowitsch if (!isStride || step != 1) PetscCall(ISGetIndices(pointIS, points)); 589305a4c7SMatthew G. Knepley PetscFunctionReturn(0); 599305a4c7SMatthew G. Knepley } 609305a4c7SMatthew G. Knepley 616528b96dSMatthew G. Knepley /*@C 62cab54364SBarry Smith ISRestorePointRange - Destroys the traversal description created with `ISGetPointRange()` 636528b96dSMatthew G. Knepley 646528b96dSMatthew G. Knepley Not collective 656528b96dSMatthew G. Knepley 666528b96dSMatthew G. Knepley Input Parameters: 67cab54364SBarry Smith + pointIS - The `IS` object 68cab54364SBarry Smith . pStart - The first index, from `ISGetPointRange()` 69cab54364SBarry Smith . pEnd - One past the last index, from `ISGetPointRange()` 70cab54364SBarry Smith - points - The indices, from `ISGetPointRange()` 716528b96dSMatthew G. Knepley 726528b96dSMatthew G. Knepley Level: intermediate 736528b96dSMatthew G. Knepley 74cab54364SBarry Smith .seealso: [](sec_scatter), `IS`, `ISGetPointRange()`, `ISGetPointSubrange()`, `ISGetIndices()`, `ISCreateStride()` 756528b96dSMatthew G. Knepley @*/ 76d71ae5a4SJacob Faibussowitsch PetscErrorCode ISRestorePointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points) 77d71ae5a4SJacob Faibussowitsch { 789305a4c7SMatthew G. Knepley PetscInt step = 1; 799305a4c7SMatthew G. Knepley PetscBool isStride; 809305a4c7SMatthew G. Knepley 819305a4c7SMatthew G. Knepley PetscFunctionBeginHot; 829566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pointIS, ISSTRIDE, &isStride)); 839566063dSJacob Faibussowitsch if (isStride) PetscCall(ISStrideGetInfo(pointIS, pStart, &step)); 849566063dSJacob Faibussowitsch if (!isStride || step != 1) PetscCall(ISGetIndices(pointIS, points)); 859305a4c7SMatthew G. Knepley PetscFunctionReturn(0); 869305a4c7SMatthew G. Knepley } 879305a4c7SMatthew G. Knepley 886528b96dSMatthew G. Knepley /*@C 89cab54364SBarry Smith ISGetPointSubrange - Configures the input `IS` to be a subrange for the traversal information given 906528b96dSMatthew G. Knepley 916528b96dSMatthew G. Knepley Not collective 926528b96dSMatthew G. Knepley 936528b96dSMatthew G. Knepley Input Parameters: 94cab54364SBarry Smith + subpointIS - The `IS` object to be configured 956528b96dSMatthew G. Knepley . pStar t - The first index of the subrange 966528b96dSMatthew G. Knepley . pEnd - One past the last index for the subrange 97cab54364SBarry Smith - points - The indices for the entire range, from `ISGetPointRange()` 986528b96dSMatthew G. Knepley 996528b96dSMatthew G. Knepley Output Parameters: 100cab54364SBarry Smith . subpointIS - The `IS` object now configured to be a subrange 1016528b96dSMatthew G. Knepley 1026528b96dSMatthew G. Knepley Level: intermediate 1036528b96dSMatthew G. Knepley 104cab54364SBarry Smith Note: 105cab54364SBarry Smith The input `IS` will now respond properly to calls to `ISGetPointRange()` and return the subrange. 106cab54364SBarry Smith 107cab54364SBarry Smith .seealso: [](sec_scatter), `IS`, `ISGetPointRange()`, `ISRestorePointRange()`, `ISGetIndices()`, `ISCreateStride()` 1086528b96dSMatthew G. Knepley @*/ 109d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGetPointSubrange(IS subpointIS, PetscInt pStart, PetscInt pEnd, const PetscInt *points) 110d71ae5a4SJacob Faibussowitsch { 1119305a4c7SMatthew G. Knepley PetscFunctionBeginHot; 1129305a4c7SMatthew G. Knepley if (points) { 1139566063dSJacob Faibussowitsch PetscCall(ISSetType(subpointIS, ISGENERAL)); 1149566063dSJacob Faibussowitsch PetscCall(ISGeneralSetIndices(subpointIS, pEnd - pStart, &points[pStart], PETSC_USE_POINTER)); 1159305a4c7SMatthew G. Knepley } else { 1169566063dSJacob Faibussowitsch PetscCall(ISSetType(subpointIS, ISSTRIDE)); 1179566063dSJacob Faibussowitsch PetscCall(ISStrideSetStride(subpointIS, pEnd - pStart, pStart, 1)); 1189305a4c7SMatthew G. Knepley } 1199305a4c7SMatthew G. Knepley PetscFunctionReturn(0); 1209305a4c7SMatthew G. Knepley } 1219305a4c7SMatthew G. Knepley 122413f72f0SBarry Smith /* -----------------------------------------------------------------------------------------*/ 123413f72f0SBarry Smith 124413f72f0SBarry Smith /* 125413f72f0SBarry Smith Creates the global mapping information in the ISLocalToGlobalMapping structure 126413f72f0SBarry Smith 127413f72f0SBarry Smith If the user has not selected how to handle the global to local mapping then use HASH for "large" problems 128413f72f0SBarry Smith */ 129d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISGlobalToLocalMappingSetUp(ISLocalToGlobalMapping mapping) 130d71ae5a4SJacob Faibussowitsch { 131413f72f0SBarry Smith PetscInt i, *idx = mapping->indices, n = mapping->n, end, start; 132413f72f0SBarry Smith 133413f72f0SBarry Smith PetscFunctionBegin; 134413f72f0SBarry Smith if (mapping->data) PetscFunctionReturn(0); 135413f72f0SBarry Smith end = 0; 136413f72f0SBarry Smith start = PETSC_MAX_INT; 137413f72f0SBarry Smith 138413f72f0SBarry Smith for (i = 0; i < n; i++) { 139413f72f0SBarry Smith if (idx[i] < 0) continue; 140413f72f0SBarry Smith if (idx[i] < start) start = idx[i]; 141413f72f0SBarry Smith if (idx[i] > end) end = idx[i]; 142413f72f0SBarry Smith } 1439371c9d4SSatish Balay if (start > end) { 1449371c9d4SSatish Balay start = 0; 1459371c9d4SSatish Balay end = -1; 1469371c9d4SSatish Balay } 147413f72f0SBarry Smith mapping->globalstart = start; 148413f72f0SBarry Smith mapping->globalend = end; 149413f72f0SBarry Smith if (!((PetscObject)mapping)->type_name) { 150413f72f0SBarry Smith if ((end - start) > PetscMax(4 * n, 1000000)) { 1519566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingSetType(mapping, ISLOCALTOGLOBALMAPPINGHASH)); 152413f72f0SBarry Smith } else { 1539566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingSetType(mapping, ISLOCALTOGLOBALMAPPINGBASIC)); 154413f72f0SBarry Smith } 155413f72f0SBarry Smith } 156dbbe0bcdSBarry Smith PetscTryTypeMethod(mapping, globaltolocalmappingsetup); 157413f72f0SBarry Smith PetscFunctionReturn(0); 158413f72f0SBarry Smith } 159413f72f0SBarry Smith 160d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISGlobalToLocalMappingSetUp_Basic(ISLocalToGlobalMapping mapping) 161d71ae5a4SJacob Faibussowitsch { 162413f72f0SBarry Smith PetscInt i, *idx = mapping->indices, n = mapping->n, end, start, *globals; 163413f72f0SBarry Smith ISLocalToGlobalMapping_Basic *map; 164413f72f0SBarry Smith 165413f72f0SBarry Smith PetscFunctionBegin; 166413f72f0SBarry Smith start = mapping->globalstart; 167413f72f0SBarry Smith end = mapping->globalend; 1689566063dSJacob Faibussowitsch PetscCall(PetscNew(&map)); 1699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(end - start + 2, &globals)); 170413f72f0SBarry Smith map->globals = globals; 171413f72f0SBarry Smith for (i = 0; i < end - start + 1; i++) globals[i] = -1; 172413f72f0SBarry Smith for (i = 0; i < n; i++) { 173413f72f0SBarry Smith if (idx[i] < 0) continue; 174413f72f0SBarry Smith globals[idx[i] - start] = i; 175413f72f0SBarry Smith } 176413f72f0SBarry Smith mapping->data = (void *)map; 177413f72f0SBarry Smith PetscFunctionReturn(0); 178413f72f0SBarry Smith } 179413f72f0SBarry Smith 180d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISGlobalToLocalMappingSetUp_Hash(ISLocalToGlobalMapping mapping) 181d71ae5a4SJacob Faibussowitsch { 182413f72f0SBarry Smith PetscInt i, *idx = mapping->indices, n = mapping->n; 183413f72f0SBarry Smith ISLocalToGlobalMapping_Hash *map; 184413f72f0SBarry Smith 185413f72f0SBarry Smith PetscFunctionBegin; 1869566063dSJacob Faibussowitsch PetscCall(PetscNew(&map)); 1879566063dSJacob Faibussowitsch PetscCall(PetscHMapICreate(&map->globalht)); 188413f72f0SBarry Smith for (i = 0; i < n; i++) { 189413f72f0SBarry Smith if (idx[i] < 0) continue; 1909566063dSJacob Faibussowitsch PetscCall(PetscHMapISet(map->globalht, idx[i], i)); 191413f72f0SBarry Smith } 192413f72f0SBarry Smith mapping->data = (void *)map; 193413f72f0SBarry Smith PetscFunctionReturn(0); 194413f72f0SBarry Smith } 195413f72f0SBarry Smith 196d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISLocalToGlobalMappingDestroy_Basic(ISLocalToGlobalMapping mapping) 197d71ae5a4SJacob Faibussowitsch { 198413f72f0SBarry Smith ISLocalToGlobalMapping_Basic *map = (ISLocalToGlobalMapping_Basic *)mapping->data; 199413f72f0SBarry Smith 200413f72f0SBarry Smith PetscFunctionBegin; 201413f72f0SBarry Smith if (!map) PetscFunctionReturn(0); 2029566063dSJacob Faibussowitsch PetscCall(PetscFree(map->globals)); 2039566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->data)); 204413f72f0SBarry Smith PetscFunctionReturn(0); 205413f72f0SBarry Smith } 206413f72f0SBarry Smith 207d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISLocalToGlobalMappingDestroy_Hash(ISLocalToGlobalMapping mapping) 208d71ae5a4SJacob Faibussowitsch { 209413f72f0SBarry Smith ISLocalToGlobalMapping_Hash *map = (ISLocalToGlobalMapping_Hash *)mapping->data; 210413f72f0SBarry Smith 211413f72f0SBarry Smith PetscFunctionBegin; 212413f72f0SBarry Smith if (!map) PetscFunctionReturn(0); 2139566063dSJacob Faibussowitsch PetscCall(PetscHMapIDestroy(&map->globalht)); 2149566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->data)); 215413f72f0SBarry Smith PetscFunctionReturn(0); 216413f72f0SBarry Smith } 217413f72f0SBarry Smith 218413f72f0SBarry Smith #define GTOLTYPE _Basic 219413f72f0SBarry Smith #define GTOLNAME _Basic 220541bf97eSAdrian Croucher #define GTOLBS mapping->bs 2219371c9d4SSatish Balay #define GTOL(g, local) \ 2229371c9d4SSatish Balay do { \ 223413f72f0SBarry Smith local = map->globals[g / bs - start]; \ 2240040bde1SJunchao Zhang if (local >= 0) local = bs * local + (g % bs); \ 225413f72f0SBarry Smith } while (0) 226541bf97eSAdrian Croucher 227413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h> 228413f72f0SBarry Smith 229413f72f0SBarry Smith #define GTOLTYPE _Basic 230413f72f0SBarry Smith #define GTOLNAME Block_Basic 231541bf97eSAdrian Croucher #define GTOLBS 1 2329371c9d4SSatish Balay #define GTOL(g, local) \ 233d71ae5a4SJacob Faibussowitsch do { \ 234d71ae5a4SJacob Faibussowitsch local = map->globals[g - start]; \ 235d71ae5a4SJacob Faibussowitsch } while (0) 236413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h> 237413f72f0SBarry Smith 238413f72f0SBarry Smith #define GTOLTYPE _Hash 239413f72f0SBarry Smith #define GTOLNAME _Hash 240541bf97eSAdrian Croucher #define GTOLBS mapping->bs 2419371c9d4SSatish Balay #define GTOL(g, local) \ 2429371c9d4SSatish Balay do { \ 243e8f14785SLisandro Dalcin (void)PetscHMapIGet(map->globalht, g / bs, &local); \ 2440040bde1SJunchao Zhang if (local >= 0) local = bs * local + (g % bs); \ 245413f72f0SBarry Smith } while (0) 246413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h> 247413f72f0SBarry Smith 248413f72f0SBarry Smith #define GTOLTYPE _Hash 249413f72f0SBarry Smith #define GTOLNAME Block_Hash 250541bf97eSAdrian Croucher #define GTOLBS 1 2519371c9d4SSatish Balay #define GTOL(g, local) \ 252d71ae5a4SJacob Faibussowitsch do { \ 253d71ae5a4SJacob Faibussowitsch (void)PetscHMapIGet(map->globalht, g, &local); \ 254d71ae5a4SJacob Faibussowitsch } while (0) 255413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h> 256413f72f0SBarry Smith 2576658fb44Sstefano_zampini /*@ 2586658fb44Sstefano_zampini ISLocalToGlobalMappingDuplicate - Duplicates the local to global mapping object 2596658fb44Sstefano_zampini 2606658fb44Sstefano_zampini Not Collective 2616658fb44Sstefano_zampini 2626658fb44Sstefano_zampini Input Parameter: 2636658fb44Sstefano_zampini . ltog - local to global mapping 2646658fb44Sstefano_zampini 2656658fb44Sstefano_zampini Output Parameter: 2666658fb44Sstefano_zampini . nltog - the duplicated local to global mapping 2676658fb44Sstefano_zampini 2686658fb44Sstefano_zampini Level: advanced 2696658fb44Sstefano_zampini 270cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()` 2716658fb44Sstefano_zampini @*/ 272d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingDuplicate(ISLocalToGlobalMapping ltog, ISLocalToGlobalMapping *nltog) 273d71ae5a4SJacob Faibussowitsch { 274a0d79125SStefano Zampini ISLocalToGlobalMappingType l2gtype; 2756658fb44Sstefano_zampini 2766658fb44Sstefano_zampini PetscFunctionBegin; 2776658fb44Sstefano_zampini PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1); 2789566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)ltog), ltog->bs, ltog->n, ltog->indices, PETSC_COPY_VALUES, nltog)); 2799566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetType(ltog, &l2gtype)); 2809566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingSetType(*nltog, l2gtype)); 2816658fb44Sstefano_zampini PetscFunctionReturn(0); 2826658fb44Sstefano_zampini } 2836658fb44Sstefano_zampini 284565245c5SBarry Smith /*@ 285107e9a97SBarry Smith ISLocalToGlobalMappingGetSize - Gets the local size of a local to global mapping 2863b9aefa3SBarry Smith 2873b9aefa3SBarry Smith Not Collective 2883b9aefa3SBarry Smith 2893b9aefa3SBarry Smith Input Parameter: 2903b9aefa3SBarry Smith . ltog - local to global mapping 2913b9aefa3SBarry Smith 2923b9aefa3SBarry Smith Output Parameter: 293cab54364SBarry Smith . n - the number of entries in the local mapping, `ISLocalToGlobalMappingGetIndices()` returns an array of this length 2943b9aefa3SBarry Smith 2953b9aefa3SBarry Smith Level: advanced 2963b9aefa3SBarry Smith 297cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()` 2983b9aefa3SBarry Smith @*/ 299d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping mapping, PetscInt *n) 300d71ae5a4SJacob Faibussowitsch { 3013b9aefa3SBarry Smith PetscFunctionBegin; 3020700a824SBarry Smith PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 3034482741eSBarry Smith PetscValidIntPointer(n, 2); 304107e9a97SBarry Smith *n = mapping->bs * mapping->n; 3053b9aefa3SBarry Smith PetscFunctionReturn(0); 3063b9aefa3SBarry Smith } 3073b9aefa3SBarry Smith 3085a5d4f66SBarry Smith /*@C 309cab54364SBarry Smith ISLocalToGlobalMappingViewFromOptions - View an `ISLocalToGlobalMapping` based on values in the options database 310fe2efc57SMark 311c3339decSBarry Smith Collective 312fe2efc57SMark 313fe2efc57SMark Input Parameters: 314fe2efc57SMark + A - the local to global mapping object 315736c3998SJose E. Roman . obj - Optional object 316736c3998SJose E. Roman - name - command line option 317fe2efc57SMark 318fe2efc57SMark Level: intermediate 319cab54364SBarry Smith 320cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingView`, `PetscObjectViewFromOptions()`, `ISLocalToGlobalMappingCreate()` 321fe2efc57SMark @*/ 322d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingViewFromOptions(ISLocalToGlobalMapping A, PetscObject obj, const char name[]) 323d71ae5a4SJacob Faibussowitsch { 324fe2efc57SMark PetscFunctionBegin; 325fe2efc57SMark PetscValidHeaderSpecific(A, IS_LTOGM_CLASSID, 1); 3269566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name)); 327fe2efc57SMark PetscFunctionReturn(0); 328fe2efc57SMark } 329fe2efc57SMark 330fe2efc57SMark /*@C 3315a5d4f66SBarry Smith ISLocalToGlobalMappingView - View a local to global mapping 3325a5d4f66SBarry Smith 333b9cd556bSLois Curfman McInnes Not Collective 334b9cd556bSLois Curfman McInnes 3355a5d4f66SBarry Smith Input Parameters: 3363b9aefa3SBarry Smith + ltog - local to global mapping 3373b9aefa3SBarry Smith - viewer - viewer 3385a5d4f66SBarry Smith 339a997ad1aSLois Curfman McInnes Level: advanced 340a997ad1aSLois Curfman McInnes 341cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()` 3425a5d4f66SBarry Smith @*/ 343d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingView(ISLocalToGlobalMapping mapping, PetscViewer viewer) 344d71ae5a4SJacob Faibussowitsch { 34532dcc486SBarry Smith PetscInt i; 34632dcc486SBarry Smith PetscMPIInt rank; 347ace3abfcSBarry Smith PetscBool iascii; 3485a5d4f66SBarry Smith 3495a5d4f66SBarry Smith PetscFunctionBegin; 3500700a824SBarry Smith PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 35148a46eb9SPierre Jolivet if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mapping), &viewer)); 3520700a824SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 3535a5d4f66SBarry Smith 3549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mapping), &rank)); 3559566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 35632077d6dSBarry Smith if (iascii) { 3579566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)mapping, viewer)); 3589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 35948a46eb9SPierre Jolivet for (i = 0; i < mapping->n; i++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %" PetscInt_FMT " %" PetscInt_FMT "\n", rank, i, mapping->indices[i])); 3609566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 3619566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 3621575c14dSBarry Smith } 3635a5d4f66SBarry Smith PetscFunctionReturn(0); 3645a5d4f66SBarry Smith } 3655a5d4f66SBarry Smith 3661f428162SBarry Smith /*@ 3672bdab257SBarry Smith ISLocalToGlobalMappingCreateIS - Creates a mapping between a local (0 to n) 3682bdab257SBarry Smith ordering and a global parallel ordering. 3692bdab257SBarry Smith 3700f5bd95cSBarry Smith Not collective 371b9cd556bSLois Curfman McInnes 372a997ad1aSLois Curfman McInnes Input Parameter: 3738c03b21aSDmitry Karpeev . is - index set containing the global numbers for each local number 3742bdab257SBarry Smith 375a997ad1aSLois Curfman McInnes Output Parameter: 3762bdab257SBarry Smith . mapping - new mapping data structure 3772bdab257SBarry Smith 378a997ad1aSLois Curfman McInnes Level: advanced 379a997ad1aSLois Curfman McInnes 380cab54364SBarry Smith Note: 381cab54364SBarry Smith the block size of the `IS` determines the block size of the mapping 382cab54364SBarry Smith 383cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetFromOptions()` 3842bdab257SBarry Smith @*/ 385d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingCreateIS(IS is, ISLocalToGlobalMapping *mapping) 386d71ae5a4SJacob Faibussowitsch { 3873bbf0e92SBarry Smith PetscInt n, bs; 3885d0c19d7SBarry Smith const PetscInt *indices; 3892bdab257SBarry Smith MPI_Comm comm; 3903bbf0e92SBarry Smith PetscBool isblock; 3913a40ed3dSBarry Smith 3923a40ed3dSBarry Smith PetscFunctionBegin; 3930700a824SBarry Smith PetscValidHeaderSpecific(is, IS_CLASSID, 1); 3944482741eSBarry Smith PetscValidPointer(mapping, 2); 3952bdab257SBarry Smith 3969566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)is, &comm)); 3979566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &n)); 3989566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)is, ISBLOCK, &isblock)); 3996006e8d2SBarry Smith if (!isblock) { 4009566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &indices)); 4019566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm, 1, n, indices, PETSC_COPY_VALUES, mapping)); 4029566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &indices)); 4036006e8d2SBarry Smith } else { 4049566063dSJacob Faibussowitsch PetscCall(ISGetBlockSize(is, &bs)); 4059566063dSJacob Faibussowitsch PetscCall(ISBlockGetIndices(is, &indices)); 4069566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm, bs, n / bs, indices, PETSC_COPY_VALUES, mapping)); 4079566063dSJacob Faibussowitsch PetscCall(ISBlockRestoreIndices(is, &indices)); 4086006e8d2SBarry Smith } 4093a40ed3dSBarry Smith PetscFunctionReturn(0); 4102bdab257SBarry Smith } 4115a5d4f66SBarry Smith 412a4d96a55SJed Brown /*@C 413a4d96a55SJed Brown ISLocalToGlobalMappingCreateSF - Creates a mapping between a local (0 to n) 414a4d96a55SJed Brown ordering and a global parallel ordering. 415a4d96a55SJed Brown 416a4d96a55SJed Brown Collective 417a4d96a55SJed Brown 418d8d19677SJose E. Roman Input Parameters: 419a4d96a55SJed Brown + sf - star forest mapping contiguous local indices to (rank, offset) 420cab54364SBarry Smith - start - first global index on this process, or `PETSC_DECIDE` to compute contiguous global numbering automatically 421a4d96a55SJed Brown 422a4d96a55SJed Brown Output Parameter: 423a4d96a55SJed Brown . mapping - new mapping data structure 424a4d96a55SJed Brown 425a4d96a55SJed Brown Level: advanced 426a4d96a55SJed Brown 4279a535bafSVaclav Hapla Notes: 428cab54364SBarry Smith If any processor calls this with start = `PETSC_DECIDE` then all processors must, otherwise the program will hang. 4299a535bafSVaclav Hapla 430cab54364SBarry Smith .seealso: [](sec_scatter), `PetscSF`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingSetFromOptions()` 431a4d96a55SJed Brown @*/ 432d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF sf, PetscInt start, ISLocalToGlobalMapping *mapping) 433d71ae5a4SJacob Faibussowitsch { 434a4d96a55SJed Brown PetscInt i, maxlocal, nroots, nleaves, *globals, *ltog; 435a4d96a55SJed Brown MPI_Comm comm; 436a4d96a55SJed Brown 437a4d96a55SJed Brown PetscFunctionBegin; 438a4d96a55SJed Brown PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 439a4d96a55SJed Brown PetscValidPointer(mapping, 3); 4409566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 44141f4c31fSVaclav Hapla PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL)); 4429a535bafSVaclav Hapla if (start == PETSC_DECIDE) { 4439a535bafSVaclav Hapla start = 0; 4449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Exscan(&nroots, &start, 1, MPIU_INT, MPI_SUM, comm)); 44541f4c31fSVaclav Hapla } else PetscCheck(start >= 0, comm, PETSC_ERR_ARG_OUTOFRANGE, "start must be nonnegative or PETSC_DECIDE"); 44641f4c31fSVaclav Hapla PetscCall(PetscSFGetLeafRange(sf, NULL, &maxlocal)); 44741f4c31fSVaclav Hapla ++maxlocal; 4489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nroots, &globals)); 4499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxlocal, <og)); 450a4d96a55SJed Brown for (i = 0; i < nroots; i++) globals[i] = start + i; 451a4d96a55SJed Brown for (i = 0; i < maxlocal; i++) ltog[i] = -1; 4529566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf, MPIU_INT, globals, ltog, MPI_REPLACE)); 4539566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf, MPIU_INT, globals, ltog, MPI_REPLACE)); 4549566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm, 1, maxlocal, ltog, PETSC_OWN_POINTER, mapping)); 4559566063dSJacob Faibussowitsch PetscCall(PetscFree(globals)); 456a4d96a55SJed Brown PetscFunctionReturn(0); 457a4d96a55SJed Brown } 458b46b645bSBarry Smith 45963fa5c83Sstefano_zampini /*@ 46063fa5c83Sstefano_zampini ISLocalToGlobalMappingSetBlockSize - Sets the blocksize of the mapping 46163fa5c83Sstefano_zampini 46263fa5c83Sstefano_zampini Not collective 46363fa5c83Sstefano_zampini 46463fa5c83Sstefano_zampini Input Parameters: 465a2b725a8SWilliam Gropp + mapping - mapping data structure 466a2b725a8SWilliam Gropp - bs - the blocksize 46763fa5c83Sstefano_zampini 46863fa5c83Sstefano_zampini Level: advanced 46963fa5c83Sstefano_zampini 470cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()` 47163fa5c83Sstefano_zampini @*/ 472d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingSetBlockSize(ISLocalToGlobalMapping mapping, PetscInt bs) 473d71ae5a4SJacob Faibussowitsch { 474a59f3c4dSstefano_zampini PetscInt *nid; 475a59f3c4dSstefano_zampini const PetscInt *oid; 476a59f3c4dSstefano_zampini PetscInt i, cn, on, obs, nn; 47763fa5c83Sstefano_zampini 47863fa5c83Sstefano_zampini PetscFunctionBegin; 47963fa5c83Sstefano_zampini PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 48008401ef6SPierre Jolivet PetscCheck(bs >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid block size %" PetscInt_FMT, bs); 48163fa5c83Sstefano_zampini if (bs == mapping->bs) PetscFunctionReturn(0); 48263fa5c83Sstefano_zampini on = mapping->n; 48363fa5c83Sstefano_zampini obs = mapping->bs; 48463fa5c83Sstefano_zampini oid = mapping->indices; 48563fa5c83Sstefano_zampini nn = (on * obs) / bs; 48608401ef6SPierre 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); 487a59f3c4dSstefano_zampini 4889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nn, &nid)); 4899566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(mapping, &oid)); 490a59f3c4dSstefano_zampini for (i = 0; i < nn; i++) { 491a59f3c4dSstefano_zampini PetscInt j; 492a59f3c4dSstefano_zampini for (j = 0, cn = 0; j < bs - 1; j++) { 4939371c9d4SSatish Balay if (oid[i * bs + j] < 0) { 4949371c9d4SSatish Balay cn++; 4959371c9d4SSatish Balay continue; 4969371c9d4SSatish Balay } 49708401ef6SPierre Jolivet PetscCheck(oid[i * bs + j] == oid[i * bs + j + 1] - 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Block sizes %" PetscInt_FMT " and %" PetscInt_FMT " are incompatible with the block indices: non consecutive indices %" PetscInt_FMT " %" PetscInt_FMT, bs, obs, oid[i * bs + j], oid[i * bs + j + 1]); 498a59f3c4dSstefano_zampini } 499a59f3c4dSstefano_zampini if (oid[i * bs + j] < 0) cn++; 5008b7cb0e6Sstefano_zampini if (cn) { 50108401ef6SPierre Jolivet PetscCheck(cn == bs, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Block sizes %" PetscInt_FMT " and %" PetscInt_FMT " are incompatible with the block indices: invalid number of negative entries in block %" PetscInt_FMT, bs, obs, cn); 502a59f3c4dSstefano_zampini nid[i] = -1; 5038b7cb0e6Sstefano_zampini } else { 504a59f3c4dSstefano_zampini nid[i] = oid[i * bs] / bs; 50563fa5c83Sstefano_zampini } 50663fa5c83Sstefano_zampini } 5079566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(mapping, &oid)); 508a59f3c4dSstefano_zampini 50963fa5c83Sstefano_zampini mapping->n = nn; 51063fa5c83Sstefano_zampini mapping->bs = bs; 5119566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->indices)); 51263fa5c83Sstefano_zampini mapping->indices = nid; 513c9345713Sstefano_zampini mapping->globalstart = 0; 514c9345713Sstefano_zampini mapping->globalend = 0; 5151bd0b88eSStefano Zampini 5161bd0b88eSStefano Zampini /* reset the cached information */ 5179566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->info_procs)); 5189566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->info_numprocs)); 5191bd0b88eSStefano Zampini if (mapping->info_indices) { 5201bd0b88eSStefano Zampini PetscInt i; 5211bd0b88eSStefano Zampini 5229566063dSJacob Faibussowitsch PetscCall(PetscFree((mapping->info_indices)[0])); 52348a46eb9SPierre Jolivet for (i = 1; i < mapping->info_nproc; i++) PetscCall(PetscFree(mapping->info_indices[i])); 5249566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->info_indices)); 5251bd0b88eSStefano Zampini } 5261bd0b88eSStefano Zampini mapping->info_cached = PETSC_FALSE; 5271bd0b88eSStefano Zampini 528dbbe0bcdSBarry Smith PetscTryTypeMethod(mapping, destroy); 52963fa5c83Sstefano_zampini PetscFunctionReturn(0); 53063fa5c83Sstefano_zampini } 53163fa5c83Sstefano_zampini 53245b6f7e9SBarry Smith /*@ 53345b6f7e9SBarry Smith ISLocalToGlobalMappingGetBlockSize - Gets the blocksize of the mapping 53445b6f7e9SBarry Smith ordering and a global parallel ordering. 53545b6f7e9SBarry Smith 53645b6f7e9SBarry Smith Not Collective 53745b6f7e9SBarry Smith 53845b6f7e9SBarry Smith Input Parameters: 53945b6f7e9SBarry Smith . mapping - mapping data structure 54045b6f7e9SBarry Smith 54145b6f7e9SBarry Smith Output Parameter: 54245b6f7e9SBarry Smith . bs - the blocksize 54345b6f7e9SBarry Smith 54445b6f7e9SBarry Smith Level: advanced 54545b6f7e9SBarry Smith 546cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()` 54745b6f7e9SBarry Smith @*/ 548d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping mapping, PetscInt *bs) 549d71ae5a4SJacob Faibussowitsch { 55045b6f7e9SBarry Smith PetscFunctionBegin; 551cbc1caf0SMatthew G. Knepley PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 55245b6f7e9SBarry Smith *bs = mapping->bs; 55345b6f7e9SBarry Smith PetscFunctionReturn(0); 55445b6f7e9SBarry Smith } 55545b6f7e9SBarry Smith 556ba5bb76aSSatish Balay /*@ 55790f02eecSBarry Smith ISLocalToGlobalMappingCreate - Creates a mapping between a local (0 to n) 55890f02eecSBarry Smith ordering and a global parallel ordering. 5592362add9SBarry Smith 56089d82c54SBarry Smith Not Collective, but communicator may have more than one process 561b9cd556bSLois Curfman McInnes 5622362add9SBarry Smith Input Parameters: 56389d82c54SBarry Smith + comm - MPI communicator 564f0413b6fSBarry Smith . bs - the block size 56528bc9809SBarry Smith . n - the number of local elements divided by the block size, or equivalently the number of block indices 56628bc9809SBarry 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 567d5ad8652SBarry Smith - mode - see PetscCopyMode 5682362add9SBarry Smith 569a997ad1aSLois Curfman McInnes Output Parameter: 57090f02eecSBarry Smith . mapping - new mapping data structure 5712362add9SBarry Smith 572cab54364SBarry Smith Level: advanced 573cab54364SBarry Smith 57495452b02SPatrick Sanan Notes: 57595452b02SPatrick Sanan There is one integer value in indices per block and it represents the actual indices bs*idx + j, where j=0,..,bs-1 576413f72f0SBarry Smith 577cab54364SBarry Smith For "small" problems when using `ISGlobalToLocalMappingApply()` and `ISGlobalToLocalMappingApplyBlock()`, the `ISLocalToGlobalMappingType` 578cab54364SBarry Smith of `ISLOCALTOGLOBALMAPPINGBASIC` will be used; this uses more memory but is faster; this approach is not scalable for extremely large mappings. 579413f72f0SBarry Smith 580cab54364SBarry Smith For large problems `ISLOCALTOGLOBALMAPPINGHASH` is used, this is scalable. 581cab54364SBarry Smith Use `ISLocalToGlobalMappingSetType()` or call `ISLocalToGlobalMappingSetFromOptions()` with the option -islocaltoglobalmapping_type <basic,hash> to control which is used. 582a997ad1aSLois Curfman McInnes 583cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingSetFromOptions()`, `ISLOCALTOGLOBALMAPPINGBASIC`, `ISLOCALTOGLOBALMAPPINGHASH` 584db781477SPatrick Sanan `ISLocalToGlobalMappingSetType()`, `ISLocalToGlobalMappingType` 5852362add9SBarry Smith @*/ 586d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingCreate(MPI_Comm comm, PetscInt bs, PetscInt n, const PetscInt indices[], PetscCopyMode mode, ISLocalToGlobalMapping *mapping) 587d71ae5a4SJacob Faibussowitsch { 58832dcc486SBarry Smith PetscInt *in; 589b46b645bSBarry Smith 590b46b645bSBarry Smith PetscFunctionBegin; 591064a246eSJacob Faibussowitsch if (n) PetscValidIntPointer(indices, 4); 592064a246eSJacob Faibussowitsch PetscValidPointer(mapping, 6); 593b46b645bSBarry Smith 5940298fd71SBarry Smith *mapping = NULL; 5959566063dSJacob Faibussowitsch PetscCall(ISInitializePackage()); 5962362add9SBarry Smith 5979566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(*mapping, IS_LTOGM_CLASSID, "ISLocalToGlobalMapping", "Local to global mapping", "IS", comm, ISLocalToGlobalMappingDestroy, ISLocalToGlobalMappingView)); 598d4bb536fSBarry Smith (*mapping)->n = n; 599f0413b6fSBarry Smith (*mapping)->bs = bs; 600d5ad8652SBarry Smith if (mode == PETSC_COPY_VALUES) { 6019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &in)); 6029566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(in, indices, n)); 603d5ad8652SBarry Smith (*mapping)->indices = in; 60471910c26SVaclav Hapla (*mapping)->dealloc_indices = PETSC_TRUE; 6056389a1a1SBarry Smith } else if (mode == PETSC_OWN_POINTER) { 6066389a1a1SBarry Smith (*mapping)->indices = (PetscInt *)indices; 60771910c26SVaclav Hapla (*mapping)->dealloc_indices = PETSC_TRUE; 60871910c26SVaclav Hapla } else if (mode == PETSC_USE_POINTER) { 60971910c26SVaclav Hapla (*mapping)->indices = (PetscInt *)indices; 6109371c9d4SSatish Balay } else SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid mode %d", mode); 6113a40ed3dSBarry Smith PetscFunctionReturn(0); 6122362add9SBarry Smith } 6132362add9SBarry Smith 614413f72f0SBarry Smith PetscFunctionList ISLocalToGlobalMappingList = NULL; 615413f72f0SBarry Smith 61690f02eecSBarry Smith /*@ 6177e99dc12SLawrence Mitchell ISLocalToGlobalMappingSetFromOptions - Set mapping options from the options database. 6187e99dc12SLawrence Mitchell 6197e99dc12SLawrence Mitchell Not collective 6207e99dc12SLawrence Mitchell 6217e99dc12SLawrence Mitchell Input Parameters: 6227e99dc12SLawrence Mitchell . mapping - mapping data structure 6237e99dc12SLawrence Mitchell 6247e99dc12SLawrence Mitchell Level: advanced 6257e99dc12SLawrence Mitchell 626cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingSetFromOptions()`, `ISLOCALTOGLOBALMAPPINGBASIC`, `ISLOCALTOGLOBALMAPPINGHASH` 627cab54364SBarry Smith `ISLocalToGlobalMappingSetType()`, `ISLocalToGlobalMappingType` 6287e99dc12SLawrence Mitchell @*/ 629d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingSetFromOptions(ISLocalToGlobalMapping mapping) 630d71ae5a4SJacob 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 656cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingCreate()` 65790f02eecSBarry Smith @*/ 658d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping *mapping) 659d71ae5a4SJacob 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 /*@ 686cab54364SBarry Smith ISLocalToGlobalMappingApplyIS - Creates from an `IS` in the local numbering 687cab54364SBarry Smith a new index set using the global numbering defined in an `ISLocalToGlobalMapping` 6883acfe500SLois Curfman McInnes context. 68990f02eecSBarry Smith 690c3339decSBarry Smith Collective 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 696cab54364SBarry Smith Output Parameter: 69790f02eecSBarry Smith . newis - index set in global numbering 69890f02eecSBarry Smith 699a997ad1aSLois Curfman McInnes Level: advanced 700a997ad1aSLois Curfman McInnes 701cab54364SBarry Smith Note: 702cab54364SBarry Smith The output `IS` will have the same communicator of the input `IS`. 703cab54364SBarry Smith 704cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingCreate()`, 705db781477SPatrick Sanan `ISLocalToGlobalMappingDestroy()`, `ISGlobalToLocalMappingApply()` 70690f02eecSBarry Smith @*/ 707d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping, IS is, IS *newis) 708d71ae5a4SJacob 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 740a997ad1aSLois Curfman McInnes Level: advanced 741a997ad1aSLois Curfman McInnes 742cab54364SBarry Smith Note: 743cab54364SBarry Smith The in and out array parameters may be identical. 744cab54364SBarry Smith 745cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApplyBlock()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingDestroy()`, 746c2e3fba1SPatrick Sanan `ISLocalToGlobalMappingApplyIS()`, `AOCreateBasic()`, `AOApplicationToPetsc()`, 747db781477SPatrick Sanan `AOPetscToApplication()`, `ISGlobalToLocalMappingApply()` 748afcb2eb5SJed Brown @*/ 749d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping, PetscInt N, const PetscInt in[], PetscInt out[]) 750d71ae5a4SJacob Faibussowitsch { 751cbc1caf0SMatthew G. Knepley PetscInt i, bs, Nmax; 75245b6f7e9SBarry Smith 75345b6f7e9SBarry Smith PetscFunctionBegin; 754cbc1caf0SMatthew G. Knepley PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 755cbc1caf0SMatthew G. Knepley bs = mapping->bs; 756cbc1caf0SMatthew G. Knepley Nmax = bs * mapping->n; 75745b6f7e9SBarry Smith if (bs == 1) { 758cbc1caf0SMatthew G. Knepley const PetscInt *idx = mapping->indices; 75945b6f7e9SBarry Smith for (i = 0; i < N; i++) { 76045b6f7e9SBarry Smith if (in[i] < 0) { 76145b6f7e9SBarry Smith out[i] = in[i]; 76245b6f7e9SBarry Smith continue; 76345b6f7e9SBarry Smith } 76408401ef6SPierre 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); 76545b6f7e9SBarry Smith out[i] = idx[in[i]]; 76645b6f7e9SBarry Smith } 76745b6f7e9SBarry Smith } else { 768cbc1caf0SMatthew G. Knepley const PetscInt *idx = mapping->indices; 76945b6f7e9SBarry Smith for (i = 0; i < N; i++) { 77045b6f7e9SBarry Smith if (in[i] < 0) { 77145b6f7e9SBarry Smith out[i] = in[i]; 77245b6f7e9SBarry Smith continue; 77345b6f7e9SBarry Smith } 77408401ef6SPierre 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); 77545b6f7e9SBarry Smith out[i] = idx[in[i] / bs] * bs + (in[i] % bs); 77645b6f7e9SBarry Smith } 77745b6f7e9SBarry Smith } 77845b6f7e9SBarry Smith PetscFunctionReturn(0); 77945b6f7e9SBarry Smith } 78045b6f7e9SBarry Smith 78145b6f7e9SBarry Smith /*@ 7826006e8d2SBarry Smith ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local block numbering and converts them to the global block numbering 78345b6f7e9SBarry Smith 78445b6f7e9SBarry Smith Not collective 78545b6f7e9SBarry Smith 78645b6f7e9SBarry Smith Input Parameters: 78745b6f7e9SBarry Smith + mapping - the local to global mapping context 78845b6f7e9SBarry Smith . N - number of integers 7896006e8d2SBarry Smith - in - input indices in local block numbering 79045b6f7e9SBarry Smith 79145b6f7e9SBarry Smith Output Parameter: 7926006e8d2SBarry Smith . out - indices in global block numbering 79345b6f7e9SBarry Smith 794cab54364SBarry Smith Level: advanced 795cab54364SBarry Smith 796cab54364SBarry Smith Note: 79745b6f7e9SBarry Smith The in and out array parameters may be identical. 79845b6f7e9SBarry Smith 7996006e8d2SBarry Smith Example: 800cab54364SBarry 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 8016006e8d2SBarry Smith (the first block) would produce 0 and the mapping applied to 1 (the second block) would produce 3. 8026006e8d2SBarry Smith 803cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingDestroy()`, 804c2e3fba1SPatrick Sanan `ISLocalToGlobalMappingApplyIS()`, `AOCreateBasic()`, `AOApplicationToPetsc()`, 805db781477SPatrick Sanan `AOPetscToApplication()`, `ISGlobalToLocalMappingApply()` 80645b6f7e9SBarry Smith @*/ 807d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping, PetscInt N, const PetscInt in[], PetscInt out[]) 808d71ae5a4SJacob Faibussowitsch { 8098a1f772fSStefano Zampini PetscInt i, Nmax; 8108a1f772fSStefano Zampini const PetscInt *idx; 811d4bb536fSBarry Smith 812a0d79125SStefano Zampini PetscFunctionBegin; 813a0d79125SStefano Zampini PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 8148a1f772fSStefano Zampini Nmax = mapping->n; 8158a1f772fSStefano Zampini idx = mapping->indices; 816afcb2eb5SJed Brown for (i = 0; i < N; i++) { 817afcb2eb5SJed Brown if (in[i] < 0) { 818afcb2eb5SJed Brown out[i] = in[i]; 819afcb2eb5SJed Brown continue; 820afcb2eb5SJed Brown } 82108401ef6SPierre 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); 822afcb2eb5SJed Brown out[i] = idx[in[i]]; 823afcb2eb5SJed Brown } 824afcb2eb5SJed Brown PetscFunctionReturn(0); 825afcb2eb5SJed Brown } 826d4bb536fSBarry Smith 8277e99dc12SLawrence Mitchell /*@ 828a997ad1aSLois Curfman McInnes ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers 829a997ad1aSLois Curfman McInnes specified with a global numbering. 830d4bb536fSBarry Smith 831b9cd556bSLois Curfman McInnes Not collective 832b9cd556bSLois Curfman McInnes 833d4bb536fSBarry Smith Input Parameters: 834b9cd556bSLois Curfman McInnes + mapping - mapping between local and global numbering 835cab54364SBarry Smith . type - `IS_GTOLM_MASK` - maps global indices with no local value to -1 in the output list (i.e., mask them) 836cab54364SBarry Smith `IS_GTOLM_DROP` - drops the indices with no local value from the output list 837d4bb536fSBarry Smith . n - number of global indices to map 838b9cd556bSLois Curfman McInnes - idx - global indices to map 839d4bb536fSBarry Smith 840d4bb536fSBarry Smith Output Parameters: 841cab54364SBarry Smith + nout - number of indices in output array (if type == `IS_GTOLM_MASK` then nout = n) 842b9cd556bSLois Curfman McInnes - idxout - local index of each global index, one must pass in an array long enough 843cab54364SBarry Smith to hold all the indices. You can call `ISGlobalToLocalMappingApply()` with 8440298fd71SBarry Smith idxout == NULL to determine the required length (returned in nout) 845cab54364SBarry Smith and then allocate the required space and call `ISGlobalToLocalMappingApply()` 846e182c471SBarry Smith a second time to set the values. 847d4bb536fSBarry Smith 848cab54364SBarry Smith Level: advanced 849cab54364SBarry Smith 850b9cd556bSLois Curfman McInnes Notes: 8510298fd71SBarry Smith Either nout or idxout may be NULL. idx and idxout may be identical. 852d4bb536fSBarry Smith 853cab54364SBarry Smith For "small" problems when using `ISGlobalToLocalMappingApply()` and `ISGlobalToLocalMappingApplyBlock()`, the `ISLocalToGlobalMappingType` of `ISLOCALTOGLOBALMAPPINGBASIC` will be used; 854cab54364SBarry 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. 855cab54364SBarry Smith Use `ISLocalToGlobalMappingSetType()` or call `ISLocalToGlobalMappingSetFromOptions()` with the option -islocaltoglobalmapping_type <basic,hash> to control which is used. 8560f5bd95cSBarry Smith 857cab54364SBarry Smith Developer Note: 858cab54364SBarry Smith The manual page states that idx and idxout may be identical but the calling 85932fd6b96SBarry Smith sequence declares idx as const so it cannot be the same as idxout. 86032fd6b96SBarry Smith 861cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApply()`, `ISGlobalToLocalMappingApplyBlock()`, `ISLocalToGlobalMappingCreate()`, 862db781477SPatrick Sanan `ISLocalToGlobalMappingDestroy()` 863d4bb536fSBarry Smith @*/ 864d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, PetscInt n, const PetscInt idx[], PetscInt *nout, PetscInt idxout[]) 865d71ae5a4SJacob Faibussowitsch { 8669d90f715SBarry Smith PetscFunctionBegin; 8679d90f715SBarry Smith PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 86848a46eb9SPierre Jolivet if (!mapping->data) PetscCall(ISGlobalToLocalMappingSetUp(mapping)); 869dbbe0bcdSBarry Smith PetscUseTypeMethod(mapping, globaltolocalmappingapply, type, n, idx, nout, idxout); 8709d90f715SBarry Smith PetscFunctionReturn(0); 8719d90f715SBarry Smith } 8729d90f715SBarry Smith 873d4fe737eSStefano Zampini /*@ 874cab54364SBarry Smith ISGlobalToLocalMappingApplyIS - Creates from an `IS` in the global numbering 875cab54364SBarry Smith a new index set using the local numbering defined in an `ISLocalToGlobalMapping` 876d4fe737eSStefano Zampini context. 877d4fe737eSStefano Zampini 878d4fe737eSStefano Zampini Not collective 879d4fe737eSStefano Zampini 880d4fe737eSStefano Zampini Input Parameters: 881d4fe737eSStefano Zampini + mapping - mapping between local and global numbering 882cab54364SBarry Smith . type - `IS_GTOLM_MASK` - maps global indices with no local value to -1 in the output list (i.e., mask them) 883cab54364SBarry Smith `IS_GTOLM_DROP` - drops the indices with no local value from the output list 884d4fe737eSStefano Zampini - is - index set in global numbering 885d4fe737eSStefano Zampini 886d4fe737eSStefano Zampini Output Parameters: 887d4fe737eSStefano Zampini . newis - index set in local numbering 888d4fe737eSStefano Zampini 889d4fe737eSStefano Zampini Level: advanced 890d4fe737eSStefano Zampini 891cab54364SBarry Smith Note: 892cab54364SBarry Smith The output `IS` will be sequential, as it encodes a purely local operation 893cab54364SBarry Smith 894cab54364SBarry Smith .seealso: [](sec_scatter), `ISGlobalToLocalMapping`, `ISGlobalToLocalMappingApply()`, `ISLocalToGlobalMappingCreate()`, 895db781477SPatrick Sanan `ISLocalToGlobalMappingDestroy()` 896d4fe737eSStefano Zampini @*/ 897d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, IS is, IS *newis) 898d71ae5a4SJacob Faibussowitsch { 899d4fe737eSStefano Zampini PetscInt n, nout, *idxout; 900d4fe737eSStefano Zampini const PetscInt *idxin; 901d4fe737eSStefano Zampini 902d4fe737eSStefano Zampini PetscFunctionBegin; 903d4fe737eSStefano Zampini PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 904d4fe737eSStefano Zampini PetscValidHeaderSpecific(is, IS_CLASSID, 3); 905d4fe737eSStefano Zampini PetscValidPointer(newis, 4); 906d4fe737eSStefano Zampini 9079566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &n)); 9089566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &idxin)); 909d4fe737eSStefano Zampini if (type == IS_GTOLM_MASK) { 9109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &idxout)); 911d4fe737eSStefano Zampini } else { 9129566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(mapping, type, n, idxin, &nout, NULL)); 9139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nout, &idxout)); 914d4fe737eSStefano Zampini } 9159566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(mapping, type, n, idxin, &nout, idxout)); 9169566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &idxin)); 9179566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nout, idxout, PETSC_OWN_POINTER, newis)); 918d4fe737eSStefano Zampini PetscFunctionReturn(0); 919d4fe737eSStefano Zampini } 920d4fe737eSStefano Zampini 9219d90f715SBarry Smith /*@ 9229d90f715SBarry Smith ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers 9239d90f715SBarry Smith specified with a block global numbering. 9249d90f715SBarry Smith 9259d90f715SBarry Smith Not collective 9269d90f715SBarry Smith 9279d90f715SBarry Smith Input Parameters: 9289d90f715SBarry Smith + mapping - mapping between local and global numbering 929cab54364SBarry Smith . type - `IS_GTOLM_MASK` - maps global indices with no local value to -1 in the output list (i.e., mask them) 930cab54364SBarry Smith `IS_GTOLM_DROP` - drops the indices with no local value from the output list 9319d90f715SBarry Smith . n - number of global indices to map 9329d90f715SBarry Smith - idx - global indices to map 9339d90f715SBarry Smith 9349d90f715SBarry Smith Output Parameters: 935cab54364SBarry Smith + nout - number of indices in output array (if type == `IS_GTOLM_MASK` then nout = n) 9369d90f715SBarry Smith - idxout - local index of each global index, one must pass in an array long enough 937cab54364SBarry Smith to hold all the indices. You can call `ISGlobalToLocalMappingApplyBlock()` with 9389d90f715SBarry Smith idxout == NULL to determine the required length (returned in nout) 939cab54364SBarry Smith and then allocate the required space and call `ISGlobalToLocalMappingApplyBlock()` 9409d90f715SBarry Smith a second time to set the values. 9419d90f715SBarry Smith 942cab54364SBarry Smith Level: advanced 943cab54364SBarry Smith 9449d90f715SBarry Smith Notes: 9459d90f715SBarry Smith Either nout or idxout may be NULL. idx and idxout may be identical. 9469d90f715SBarry Smith 947cab54364SBarry Smith For "small" problems when using `ISGlobalToLocalMappingApply()` and `ISGlobalToLocalMappingApplyBlock()`, the `ISLocalToGlobalMappingType` of `ISLOCALTOGLOBALMAPPINGBASIC` will be used; 948cab54364SBarry 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. 949cab54364SBarry Smith Use `ISLocalToGlobalMappingSetType()` or call `ISLocalToGlobalMappingSetFromOptions()` with the option -islocaltoglobalmapping_type <basic,hash> to control which is used. 9509d90f715SBarry Smith 951cab54364SBarry Smith Developer Note: 952cab54364SBarry Smith 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 955cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApply()`, `ISGlobalToLocalMappingApply()`, `ISLocalToGlobalMappingCreate()`, 956db781477SPatrick Sanan `ISLocalToGlobalMappingDestroy()` 9579d90f715SBarry Smith @*/ 958d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, PetscInt n, const PetscInt idx[], PetscInt *nout, PetscInt idxout[]) 959d71ae5a4SJacob 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 971c3339decSBarry Smith Collective 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: 985cab54364SBarry Smith .vb 9862cfcea29SBarry Smith PetscInt indices[nproc][numprocmax],ierr) 987cab54364SBarry Smith ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by 988cab54364SBarry Smith ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc], 989cab54364SBarry Smith .ve 990cab54364SBarry Smith 991cab54364SBarry Smith Fortran Note: 992cab54364SBarry Smith There is no `ISLocalToGlobalMappingRestoreInfo()` in Fortran. You must make sure that procs[], numprocs[] and 9932cfcea29SBarry Smith indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough. 9942cfcea29SBarry Smith 995cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 996db781477SPatrick Sanan `ISLocalToGlobalMappingRestoreInfo()` 99789d82c54SBarry Smith @*/ 998d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[]) 999d71ae5a4SJacob Faibussowitsch { 1000268a049cSStefano Zampini PetscFunctionBegin; 1001268a049cSStefano Zampini PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 1002268a049cSStefano Zampini if (mapping->info_cached) { 1003268a049cSStefano Zampini *nproc = mapping->info_nproc; 1004268a049cSStefano Zampini *procs = mapping->info_procs; 1005268a049cSStefano Zampini *numprocs = mapping->info_numprocs; 1006268a049cSStefano Zampini *indices = mapping->info_indices; 1007268a049cSStefano Zampini } else { 10089566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockInfo_Private(mapping, nproc, procs, numprocs, indices)); 1009268a049cSStefano Zampini } 1010268a049cSStefano Zampini PetscFunctionReturn(0); 1011268a049cSStefano Zampini } 1012268a049cSStefano Zampini 1013d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[]) 1014d71ae5a4SJacob Faibussowitsch { 101597f1f81fSBarry Smith PetscMPIInt size, rank, tag1, tag2, tag3, *len, *source, imdex; 101632dcc486SBarry Smith PetscInt i, n = mapping->n, Ng, ng, max = 0, *lindices = mapping->indices; 101732dcc486SBarry Smith PetscInt *nprocs, *owner, nsends, *sends, j, *starts, nmax, nrecvs, *recvs, proc; 1018c599c493SJunchao Zhang PetscInt cnt, scale, *ownedsenders, *nownedsenders, rstart; 101932dcc486SBarry Smith PetscInt node, nownedm, nt, *sends2, nsends2, *starts2, *lens2, *dest, nrecvs2, *starts3, *recvs2, k, *bprocs, *tmp; 102032dcc486SBarry Smith PetscInt first_procs, first_numprocs, *first_indices; 102189d82c54SBarry Smith MPI_Request *recv_waits, *send_waits; 102230dcb7c9SBarry Smith MPI_Status recv_status, *send_status, *recv_statuses; 1023ce94432eSBarry Smith MPI_Comm comm; 1024ace3abfcSBarry Smith PetscBool debug = PETSC_FALSE; 102589d82c54SBarry Smith 102689d82c54SBarry Smith PetscFunctionBegin; 10279566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)mapping, &comm)); 10289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 10299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 103024cf384cSBarry Smith if (size == 1) { 103124cf384cSBarry Smith *nproc = 0; 10320298fd71SBarry Smith *procs = NULL; 10339566063dSJacob Faibussowitsch PetscCall(PetscNew(numprocs)); 10341e2105dcSBarry Smith (*numprocs)[0] = 0; 10359566063dSJacob Faibussowitsch PetscCall(PetscNew(indices)); 10360298fd71SBarry Smith (*indices)[0] = NULL; 1037268a049cSStefano Zampini /* save info for reuse */ 1038268a049cSStefano Zampini mapping->info_nproc = *nproc; 1039268a049cSStefano Zampini mapping->info_procs = *procs; 1040268a049cSStefano Zampini mapping->info_numprocs = *numprocs; 1041268a049cSStefano Zampini mapping->info_indices = *indices; 1042268a049cSStefano Zampini mapping->info_cached = PETSC_TRUE; 104324cf384cSBarry Smith PetscFunctionReturn(0); 104424cf384cSBarry Smith } 104524cf384cSBarry Smith 10469566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)mapping)->options, NULL, "-islocaltoglobalmappinggetinfo_debug", &debug, NULL)); 104707b52d57SBarry Smith 10483677ff5aSBarry Smith /* 10496a818285SBarry Smith Notes on ISLocalToGlobalMappingGetBlockInfo 10503677ff5aSBarry Smith 10513677ff5aSBarry Smith globally owned node - the nodes that have been assigned to this processor in global 10523677ff5aSBarry Smith numbering, just for this routine. 10533677ff5aSBarry Smith 10543677ff5aSBarry Smith nontrivial globally owned node - node assigned to this processor that is on a subdomain 10553677ff5aSBarry Smith boundary (i.e. is has more than one local owner) 10563677ff5aSBarry Smith 10573677ff5aSBarry Smith locally owned node - node that exists on this processors subdomain 10583677ff5aSBarry Smith 10593677ff5aSBarry Smith nontrivial locally owned node - node that is not in the interior (i.e. has more than one 10603677ff5aSBarry Smith local subdomain 10613677ff5aSBarry Smith */ 10629566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)mapping, &tag1)); 10639566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)mapping, &tag2)); 10649566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)mapping, &tag3)); 106589d82c54SBarry Smith 106689d82c54SBarry Smith for (i = 0; i < n; i++) { 106789d82c54SBarry Smith if (lindices[i] > max) max = lindices[i]; 106889d82c54SBarry Smith } 10691c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&max, &Ng, 1, MPIU_INT, MPI_MAX, comm)); 107078058e43SBarry Smith Ng++; 10719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 10729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1073bc8ff85bSBarry Smith scale = Ng / size + 1; 10749371c9d4SSatish Balay ng = scale; 10759371c9d4SSatish Balay if (rank == size - 1) ng = Ng - scale * (size - 1); 10769371c9d4SSatish Balay ng = PetscMax(1, ng); 1077caba0dd0SBarry Smith rstart = scale * rank; 107889d82c54SBarry Smith 107989d82c54SBarry Smith /* determine ownership ranges of global indices */ 10809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * size, &nprocs)); 10819566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(nprocs, 2 * size)); 108289d82c54SBarry Smith 108389d82c54SBarry Smith /* determine owners of each local node */ 10849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &owner)); 108589d82c54SBarry Smith for (i = 0; i < n; i++) { 10863677ff5aSBarry Smith proc = lindices[i] / scale; /* processor that globally owns this index */ 108727c402fcSBarry Smith nprocs[2 * proc + 1] = 1; /* processor globally owns at least one of ours */ 10883677ff5aSBarry Smith owner[i] = proc; 108927c402fcSBarry Smith nprocs[2 * proc]++; /* count of how many that processor globally owns of ours */ 109089d82c54SBarry Smith } 10919371c9d4SSatish Balay nsends = 0; 10929371c9d4SSatish Balay for (i = 0; i < size; i++) nsends += nprocs[2 * i + 1]; 10939566063dSJacob Faibussowitsch PetscCall(PetscInfo(mapping, "Number of global owners for my local data %" PetscInt_FMT "\n", nsends)); 109489d82c54SBarry Smith 109589d82c54SBarry Smith /* inform other processors of number of messages and max length*/ 10969566063dSJacob Faibussowitsch PetscCall(PetscMaxSum(comm, nprocs, &nmax, &nrecvs)); 10979566063dSJacob Faibussowitsch PetscCall(PetscInfo(mapping, "Number of local owners for my global data %" PetscInt_FMT "\n", nrecvs)); 109889d82c54SBarry Smith 109989d82c54SBarry Smith /* post receives for owned rows */ 11009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((2 * nrecvs + 1) * (nmax + 1), &recvs)); 11019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs + 1, &recv_waits)); 110248a46eb9SPierre 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)); 110389d82c54SBarry Smith 110489d82c54SBarry Smith /* pack messages containing lists of local nodes to owners */ 11059566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * n + 1, &sends)); 11069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size + 1, &starts)); 110789d82c54SBarry Smith starts[0] = 0; 1108f6e5521dSKarl Rupp for (i = 1; i < size; i++) starts[i] = starts[i - 1] + 2 * nprocs[2 * i - 2]; 110989d82c54SBarry Smith for (i = 0; i < n; i++) { 111089d82c54SBarry Smith sends[starts[owner[i]]++] = lindices[i]; 111130dcb7c9SBarry Smith sends[starts[owner[i]]++] = i; 111289d82c54SBarry Smith } 11139566063dSJacob Faibussowitsch PetscCall(PetscFree(owner)); 111489d82c54SBarry Smith starts[0] = 0; 1115f6e5521dSKarl Rupp for (i = 1; i < size; i++) starts[i] = starts[i - 1] + 2 * nprocs[2 * i - 2]; 111689d82c54SBarry Smith 111789d82c54SBarry Smith /* send the messages */ 11189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends + 1, &send_waits)); 11199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends + 1, &dest)); 112089d82c54SBarry Smith cnt = 0; 112189d82c54SBarry Smith for (i = 0; i < size; i++) { 112227c402fcSBarry Smith if (nprocs[2 * i]) { 11239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Isend(sends + starts[i], 2 * nprocs[2 * i], MPIU_INT, i, tag1, comm, send_waits + cnt)); 112430dcb7c9SBarry Smith dest[cnt] = i; 112589d82c54SBarry Smith cnt++; 112689d82c54SBarry Smith } 112789d82c54SBarry Smith } 11289566063dSJacob Faibussowitsch PetscCall(PetscFree(starts)); 112989d82c54SBarry Smith 113089d82c54SBarry Smith /* wait on receives */ 11319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs + 1, &source)); 11329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs + 1, &len)); 113389d82c54SBarry Smith cnt = nrecvs; 11349566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(ng + 1, &nownedsenders)); 113589d82c54SBarry Smith while (cnt) { 11369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitany(nrecvs, recv_waits, &imdex, &recv_status)); 113789d82c54SBarry Smith /* unpack receives into our local space */ 11389566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&recv_status, MPIU_INT, &len[imdex])); 113989d82c54SBarry Smith source[imdex] = recv_status.MPI_SOURCE; 114030dcb7c9SBarry Smith len[imdex] = len[imdex] / 2; 1141caba0dd0SBarry Smith /* count how many local owners for each of my global owned indices */ 114230dcb7c9SBarry Smith for (i = 0; i < len[imdex]; i++) nownedsenders[recvs[2 * imdex * nmax + 2 * i] - rstart]++; 114389d82c54SBarry Smith cnt--; 114489d82c54SBarry Smith } 11459566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_waits)); 114689d82c54SBarry Smith 114730dcb7c9SBarry Smith /* count how many globally owned indices are on an edge multiplied by how many processors own them. */ 1148bc8ff85bSBarry Smith nownedm = 0; 1149bc8ff85bSBarry Smith for (i = 0; i < ng; i++) { 1150c599c493SJunchao Zhang if (nownedsenders[i] > 1) nownedm += nownedsenders[i]; 1151bc8ff85bSBarry Smith } 1152bc8ff85bSBarry Smith 1153bc8ff85bSBarry Smith /* create single array to contain rank of all local owners of each globally owned index */ 11549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nownedm + 1, &ownedsenders)); 11559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ng + 1, &starts)); 1156bc8ff85bSBarry Smith starts[0] = 0; 1157bc8ff85bSBarry Smith for (i = 1; i < ng; i++) { 1158bc8ff85bSBarry Smith if (nownedsenders[i - 1] > 1) starts[i] = starts[i - 1] + nownedsenders[i - 1]; 1159bc8ff85bSBarry Smith else starts[i] = starts[i - 1]; 1160bc8ff85bSBarry Smith } 1161bc8ff85bSBarry Smith 11626aad120cSJose E. Roman /* for each nontrivial globally owned node list all arriving processors */ 1163bc8ff85bSBarry Smith for (i = 0; i < nrecvs; i++) { 1164bc8ff85bSBarry Smith for (j = 0; j < len[i]; j++) { 116530dcb7c9SBarry Smith node = recvs[2 * i * nmax + 2 * j] - rstart; 1166f6e5521dSKarl Rupp if (nownedsenders[node] > 1) ownedsenders[starts[node]++] = source[i]; 1167bc8ff85bSBarry Smith } 1168bc8ff85bSBarry Smith } 1169bc8ff85bSBarry Smith 117007b52d57SBarry Smith if (debug) { /* ----------------------------------- */ 117130dcb7c9SBarry Smith starts[0] = 0; 117230dcb7c9SBarry Smith for (i = 1; i < ng; i++) { 117330dcb7c9SBarry Smith if (nownedsenders[i - 1] > 1) starts[i] = starts[i - 1] + nownedsenders[i - 1]; 117430dcb7c9SBarry Smith else starts[i] = starts[i - 1]; 117530dcb7c9SBarry Smith } 117630dcb7c9SBarry Smith for (i = 0; i < ng; i++) { 117730dcb7c9SBarry Smith if (nownedsenders[i] > 1) { 11789566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm, "[%d] global node %" PetscInt_FMT " local owner processors: ", rank, i + rstart)); 117948a46eb9SPierre Jolivet for (j = 0; j < nownedsenders[i]; j++) PetscCall(PetscSynchronizedPrintf(comm, "%" PetscInt_FMT " ", ownedsenders[starts[i] + j])); 11809566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm, "\n")); 118130dcb7c9SBarry Smith } 118230dcb7c9SBarry Smith } 11839566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 118407b52d57SBarry Smith } /* ----------------------------------- */ 118530dcb7c9SBarry Smith 11863677ff5aSBarry Smith /* wait on original sends */ 11873a96401aSBarry Smith if (nsends) { 11889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends, &send_status)); 11899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nsends, send_waits, send_status)); 11909566063dSJacob Faibussowitsch PetscCall(PetscFree(send_status)); 11913a96401aSBarry Smith } 11929566063dSJacob Faibussowitsch PetscCall(PetscFree(send_waits)); 11939566063dSJacob Faibussowitsch PetscCall(PetscFree(sends)); 11949566063dSJacob Faibussowitsch PetscCall(PetscFree(nprocs)); 11953677ff5aSBarry Smith 11963677ff5aSBarry Smith /* pack messages to send back to local owners */ 119730dcb7c9SBarry Smith starts[0] = 0; 119830dcb7c9SBarry Smith for (i = 1; i < ng; i++) { 119930dcb7c9SBarry Smith if (nownedsenders[i - 1] > 1) starts[i] = starts[i - 1] + nownedsenders[i - 1]; 120030dcb7c9SBarry Smith else starts[i] = starts[i - 1]; 120130dcb7c9SBarry Smith } 120230dcb7c9SBarry Smith nsends2 = nrecvs; 12039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends2 + 1, &nprocs)); /* length of each message */ 120430dcb7c9SBarry Smith for (i = 0; i < nrecvs; i++) { 120530dcb7c9SBarry Smith nprocs[i] = 1; 120630dcb7c9SBarry Smith for (j = 0; j < len[i]; j++) { 120730dcb7c9SBarry Smith node = recvs[2 * i * nmax + 2 * j] - rstart; 1208f6e5521dSKarl Rupp if (nownedsenders[node] > 1) nprocs[i] += 2 + nownedsenders[node]; 120930dcb7c9SBarry Smith } 121030dcb7c9SBarry Smith } 1211f6e5521dSKarl Rupp nt = 0; 1212f6e5521dSKarl Rupp for (i = 0; i < nsends2; i++) nt += nprocs[i]; 1213f6e5521dSKarl Rupp 12149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nt + 1, &sends2)); 12159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends2 + 1, &starts2)); 1216f6e5521dSKarl Rupp 1217f6e5521dSKarl Rupp starts2[0] = 0; 1218f6e5521dSKarl Rupp for (i = 1; i < nsends2; i++) starts2[i] = starts2[i - 1] + nprocs[i - 1]; 121930dcb7c9SBarry Smith /* 122030dcb7c9SBarry Smith Each message is 1 + nprocs[i] long, and consists of 122130dcb7c9SBarry Smith (0) the number of nodes being sent back 122230dcb7c9SBarry Smith (1) the local node number, 122330dcb7c9SBarry Smith (2) the number of processors sharing it, 122430dcb7c9SBarry Smith (3) the processors sharing it 122530dcb7c9SBarry Smith */ 122630dcb7c9SBarry Smith for (i = 0; i < nsends2; i++) { 122730dcb7c9SBarry Smith cnt = 1; 122830dcb7c9SBarry Smith sends2[starts2[i]] = 0; 122930dcb7c9SBarry Smith for (j = 0; j < len[i]; j++) { 123030dcb7c9SBarry Smith node = recvs[2 * i * nmax + 2 * j] - rstart; 123130dcb7c9SBarry Smith if (nownedsenders[node] > 1) { 123230dcb7c9SBarry Smith sends2[starts2[i]]++; 123330dcb7c9SBarry Smith sends2[starts2[i] + cnt++] = recvs[2 * i * nmax + 2 * j + 1]; 123430dcb7c9SBarry Smith sends2[starts2[i] + cnt++] = nownedsenders[node]; 12359566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(&sends2[starts2[i] + cnt], &ownedsenders[starts[node]], nownedsenders[node])); 123630dcb7c9SBarry Smith cnt += nownedsenders[node]; 123730dcb7c9SBarry Smith } 123830dcb7c9SBarry Smith } 123930dcb7c9SBarry Smith } 124030dcb7c9SBarry Smith 124130dcb7c9SBarry Smith /* receive the message lengths */ 124230dcb7c9SBarry Smith nrecvs2 = nsends; 12439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2 + 1, &lens2)); 12449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2 + 1, &starts3)); 12459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2 + 1, &recv_waits)); 124648a46eb9SPierre Jolivet for (i = 0; i < nrecvs2; i++) PetscCallMPI(MPI_Irecv(&lens2[i], 1, MPIU_INT, dest[i], tag2, comm, recv_waits + i)); 1247d44834fbSBarry Smith 12488a8e0b3aSBarry Smith /* send the message lengths */ 124948a46eb9SPierre Jolivet for (i = 0; i < nsends2; i++) PetscCallMPI(MPI_Send(&nprocs[i], 1, MPIU_INT, source[i], tag2, comm)); 12508a8e0b3aSBarry Smith 1251d44834fbSBarry Smith /* wait on receives of lens */ 12520c468ba9SBarry Smith if (nrecvs2) { 12539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2, &recv_statuses)); 12549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrecvs2, recv_waits, recv_statuses)); 12559566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_statuses)); 12560c468ba9SBarry Smith } 12579566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_waits)); 1258d44834fbSBarry Smith 125930dcb7c9SBarry Smith starts3[0] = 0; 1260d44834fbSBarry Smith nt = 0; 126130dcb7c9SBarry Smith for (i = 0; i < nrecvs2 - 1; i++) { 126230dcb7c9SBarry Smith starts3[i + 1] = starts3[i] + lens2[i]; 1263d44834fbSBarry Smith nt += lens2[i]; 126430dcb7c9SBarry Smith } 126576466f69SStefano Zampini if (nrecvs2) nt += lens2[nrecvs2 - 1]; 1266d44834fbSBarry Smith 12679566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nt + 1, &recvs2)); 12689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2 + 1, &recv_waits)); 126948a46eb9SPierre Jolivet for (i = 0; i < nrecvs2; i++) PetscCallMPI(MPI_Irecv(recvs2 + starts3[i], lens2[i], MPIU_INT, dest[i], tag3, comm, recv_waits + i)); 127030dcb7c9SBarry Smith 127130dcb7c9SBarry Smith /* send the messages */ 12729566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends2 + 1, &send_waits)); 127348a46eb9SPierre Jolivet for (i = 0; i < nsends2; i++) PetscCallMPI(MPI_Isend(sends2 + starts2[i], nprocs[i], MPIU_INT, source[i], tag3, comm, send_waits + i)); 127430dcb7c9SBarry Smith 127530dcb7c9SBarry Smith /* wait on receives */ 12760c468ba9SBarry Smith if (nrecvs2) { 12779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2, &recv_statuses)); 12789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrecvs2, recv_waits, recv_statuses)); 12799566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_statuses)); 12800c468ba9SBarry Smith } 12819566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_waits)); 12829566063dSJacob Faibussowitsch PetscCall(PetscFree(nprocs)); 128330dcb7c9SBarry Smith 128407b52d57SBarry Smith if (debug) { /* ----------------------------------- */ 128530dcb7c9SBarry Smith cnt = 0; 128630dcb7c9SBarry Smith for (i = 0; i < nrecvs2; i++) { 128730dcb7c9SBarry Smith nt = recvs2[cnt++]; 128830dcb7c9SBarry Smith for (j = 0; j < nt; j++) { 12899566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm, "[%d] local node %" PetscInt_FMT " number of subdomains %" PetscInt_FMT ": ", rank, recvs2[cnt], recvs2[cnt + 1])); 129048a46eb9SPierre Jolivet for (k = 0; k < recvs2[cnt + 1]; k++) PetscCall(PetscSynchronizedPrintf(comm, "%" PetscInt_FMT " ", recvs2[cnt + 2 + k])); 129130dcb7c9SBarry Smith cnt += 2 + recvs2[cnt + 1]; 12929566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm, "\n")); 129330dcb7c9SBarry Smith } 129430dcb7c9SBarry Smith } 12959566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 129607b52d57SBarry Smith } /* ----------------------------------- */ 129730dcb7c9SBarry Smith 129830dcb7c9SBarry Smith /* count number subdomains for each local node */ 12999566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(size, &nprocs)); 130030dcb7c9SBarry Smith cnt = 0; 130130dcb7c9SBarry Smith for (i = 0; i < nrecvs2; i++) { 130230dcb7c9SBarry Smith nt = recvs2[cnt++]; 130330dcb7c9SBarry Smith for (j = 0; j < nt; j++) { 1304f6e5521dSKarl Rupp for (k = 0; k < recvs2[cnt + 1]; k++) nprocs[recvs2[cnt + 2 + k]]++; 130530dcb7c9SBarry Smith cnt += 2 + recvs2[cnt + 1]; 130630dcb7c9SBarry Smith } 130730dcb7c9SBarry Smith } 13089371c9d4SSatish Balay nt = 0; 13099371c9d4SSatish Balay for (i = 0; i < size; i++) nt += (nprocs[i] > 0); 131030dcb7c9SBarry Smith *nproc = nt; 13119566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nt + 1, procs)); 13129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nt + 1, numprocs)); 13139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nt + 1, indices)); 13140298fd71SBarry Smith for (i = 0; i < nt + 1; i++) (*indices)[i] = NULL; 13159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &bprocs)); 131630dcb7c9SBarry Smith cnt = 0; 131730dcb7c9SBarry Smith for (i = 0; i < size; i++) { 131830dcb7c9SBarry Smith if (nprocs[i] > 0) { 131930dcb7c9SBarry Smith bprocs[i] = cnt; 132030dcb7c9SBarry Smith (*procs)[cnt] = i; 132130dcb7c9SBarry Smith (*numprocs)[cnt] = nprocs[i]; 13229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nprocs[i], &(*indices)[cnt])); 132330dcb7c9SBarry Smith cnt++; 132430dcb7c9SBarry Smith } 132530dcb7c9SBarry Smith } 132630dcb7c9SBarry Smith 132730dcb7c9SBarry Smith /* make the list of subdomains for each nontrivial local node */ 13289566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(*numprocs, nt)); 132930dcb7c9SBarry Smith cnt = 0; 133030dcb7c9SBarry Smith for (i = 0; i < nrecvs2; i++) { 133130dcb7c9SBarry Smith nt = recvs2[cnt++]; 133230dcb7c9SBarry Smith for (j = 0; j < nt; j++) { 1333f6e5521dSKarl Rupp for (k = 0; k < recvs2[cnt + 1]; k++) (*indices)[bprocs[recvs2[cnt + 2 + k]]][(*numprocs)[bprocs[recvs2[cnt + 2 + k]]]++] = recvs2[cnt]; 133430dcb7c9SBarry Smith cnt += 2 + recvs2[cnt + 1]; 133530dcb7c9SBarry Smith } 133630dcb7c9SBarry Smith } 13379566063dSJacob Faibussowitsch PetscCall(PetscFree(bprocs)); 13389566063dSJacob Faibussowitsch PetscCall(PetscFree(recvs2)); 133930dcb7c9SBarry Smith 134007b52d57SBarry Smith /* sort the node indexing by their global numbers */ 134107b52d57SBarry Smith nt = *nproc; 134207b52d57SBarry Smith for (i = 0; i < nt; i++) { 13439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((*numprocs)[i], &tmp)); 1344f6e5521dSKarl Rupp for (j = 0; j < (*numprocs)[i]; j++) tmp[j] = lindices[(*indices)[i][j]]; 13459566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithArray((*numprocs)[i], tmp, (*indices)[i])); 13469566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp)); 134707b52d57SBarry Smith } 134807b52d57SBarry Smith 134907b52d57SBarry Smith if (debug) { /* ----------------------------------- */ 135030dcb7c9SBarry Smith nt = *nproc; 135130dcb7c9SBarry Smith for (i = 0; i < nt; i++) { 13529566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm, "[%d] subdomain %" PetscInt_FMT " number of indices %" PetscInt_FMT ": ", rank, (*procs)[i], (*numprocs)[i])); 135348a46eb9SPierre Jolivet for (j = 0; j < (*numprocs)[i]; j++) PetscCall(PetscSynchronizedPrintf(comm, "%" PetscInt_FMT " ", (*indices)[i][j])); 13549566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm, "\n")); 135530dcb7c9SBarry Smith } 13569566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 135707b52d57SBarry Smith } /* ----------------------------------- */ 135830dcb7c9SBarry Smith 135930dcb7c9SBarry Smith /* wait on sends */ 136030dcb7c9SBarry Smith if (nsends2) { 13619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends2, &send_status)); 13629566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nsends2, send_waits, send_status)); 13639566063dSJacob Faibussowitsch PetscCall(PetscFree(send_status)); 136430dcb7c9SBarry Smith } 136530dcb7c9SBarry Smith 13669566063dSJacob Faibussowitsch PetscCall(PetscFree(starts3)); 13679566063dSJacob Faibussowitsch PetscCall(PetscFree(dest)); 13689566063dSJacob Faibussowitsch PetscCall(PetscFree(send_waits)); 13693677ff5aSBarry Smith 13709566063dSJacob Faibussowitsch PetscCall(PetscFree(nownedsenders)); 13719566063dSJacob Faibussowitsch PetscCall(PetscFree(ownedsenders)); 13729566063dSJacob Faibussowitsch PetscCall(PetscFree(starts)); 13739566063dSJacob Faibussowitsch PetscCall(PetscFree(starts2)); 13749566063dSJacob Faibussowitsch PetscCall(PetscFree(lens2)); 137589d82c54SBarry Smith 13769566063dSJacob Faibussowitsch PetscCall(PetscFree(source)); 13779566063dSJacob Faibussowitsch PetscCall(PetscFree(len)); 13789566063dSJacob Faibussowitsch PetscCall(PetscFree(recvs)); 13799566063dSJacob Faibussowitsch PetscCall(PetscFree(nprocs)); 13809566063dSJacob Faibussowitsch PetscCall(PetscFree(sends2)); 138124cf384cSBarry Smith 138224cf384cSBarry Smith /* put the information about myself as the first entry in the list */ 138324cf384cSBarry Smith first_procs = (*procs)[0]; 138424cf384cSBarry Smith first_numprocs = (*numprocs)[0]; 138524cf384cSBarry Smith first_indices = (*indices)[0]; 138624cf384cSBarry Smith for (i = 0; i < *nproc; i++) { 138724cf384cSBarry Smith if ((*procs)[i] == rank) { 138824cf384cSBarry Smith (*procs)[0] = (*procs)[i]; 138924cf384cSBarry Smith (*numprocs)[0] = (*numprocs)[i]; 139024cf384cSBarry Smith (*indices)[0] = (*indices)[i]; 139124cf384cSBarry Smith (*procs)[i] = first_procs; 139224cf384cSBarry Smith (*numprocs)[i] = first_numprocs; 139324cf384cSBarry Smith (*indices)[i] = first_indices; 139424cf384cSBarry Smith break; 139524cf384cSBarry Smith } 139624cf384cSBarry Smith } 1397268a049cSStefano Zampini 1398268a049cSStefano Zampini /* save info for reuse */ 1399268a049cSStefano Zampini mapping->info_nproc = *nproc; 1400268a049cSStefano Zampini mapping->info_procs = *procs; 1401268a049cSStefano Zampini mapping->info_numprocs = *numprocs; 1402268a049cSStefano Zampini mapping->info_indices = *indices; 1403268a049cSStefano Zampini mapping->info_cached = PETSC_TRUE; 140489d82c54SBarry Smith PetscFunctionReturn(0); 140589d82c54SBarry Smith } 140689d82c54SBarry Smith 14076a818285SBarry Smith /*@C 1408cab54364SBarry Smith ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by `ISLocalToGlobalMappingGetBlockInfo()` 14096a818285SBarry Smith 1410c3339decSBarry Smith Collective 14116a818285SBarry Smith 1412f899ff85SJose E. Roman Input Parameter: 14136a818285SBarry Smith . mapping - the mapping from local to global indexing 14146a818285SBarry Smith 1415d8d19677SJose E. Roman Output Parameters: 14166a818285SBarry Smith + nproc - number of processors that are connected to this one 14176a818285SBarry Smith . proc - neighboring processors 14186a818285SBarry Smith . numproc - number of indices for each processor 14196a818285SBarry Smith - indices - indices of local nodes shared with neighbor (sorted by global numbering) 14206a818285SBarry Smith 14216a818285SBarry Smith Level: advanced 14226a818285SBarry Smith 1423cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 1424db781477SPatrick Sanan `ISLocalToGlobalMappingGetInfo()` 14256a818285SBarry Smith @*/ 1426d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[]) 1427d71ae5a4SJacob Faibussowitsch { 14286a818285SBarry Smith PetscFunctionBegin; 1429cbc1caf0SMatthew G. Knepley PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 1430268a049cSStefano Zampini if (mapping->info_free) { 14319566063dSJacob Faibussowitsch PetscCall(PetscFree(*numprocs)); 14326a818285SBarry Smith if (*indices) { 1433268a049cSStefano Zampini PetscInt i; 1434268a049cSStefano Zampini 14359566063dSJacob Faibussowitsch PetscCall(PetscFree((*indices)[0])); 143648a46eb9SPierre Jolivet for (i = 1; i < *nproc; i++) PetscCall(PetscFree((*indices)[i])); 14379566063dSJacob Faibussowitsch PetscCall(PetscFree(*indices)); 14386a818285SBarry Smith } 1439268a049cSStefano Zampini } 1440268a049cSStefano Zampini *nproc = 0; 1441268a049cSStefano Zampini *procs = NULL; 1442268a049cSStefano Zampini *numprocs = NULL; 1443268a049cSStefano Zampini *indices = NULL; 14446a818285SBarry Smith PetscFunctionReturn(0); 14456a818285SBarry Smith } 14466a818285SBarry Smith 14476a818285SBarry Smith /*@C 14486a818285SBarry Smith ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each processor and 14496a818285SBarry Smith each index shared by more than one processor 14506a818285SBarry Smith 1451c3339decSBarry Smith Collective 14526a818285SBarry Smith 1453f899ff85SJose E. Roman Input Parameter: 14546a818285SBarry Smith . mapping - the mapping from local to global indexing 14556a818285SBarry Smith 1456d8d19677SJose E. Roman Output Parameters: 14576a818285SBarry Smith + nproc - number of processors that are connected to this one 14586a818285SBarry Smith . proc - neighboring processors 14596a818285SBarry Smith . numproc - number of indices for each subdomain (processor) 14606a818285SBarry Smith - indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering) 14616a818285SBarry Smith 14626a818285SBarry Smith Level: advanced 14636a818285SBarry Smith 1464cab54364SBarry Smith Note: 1465cab54364SBarry Smith The user needs to call `ISLocalToGlobalMappingRestoreInfo()` when the data is no longer needed. 14661bd0b88eSStefano Zampini 14676a818285SBarry Smith Fortran Usage: 1468cab54364SBarry Smith .vb 14696a818285SBarry Smith PetscInt indices[nproc][numprocmax],ierr) 1470cab54364SBarry Smith ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by 1471cab54364SBarry Smith ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc], 1472cab54364SBarry Smith .ve 1473cab54364SBarry Smith 1474cab54364SBarry Smith Fortran Note: 1475cab54364SBarry Smith There is no `ISLocalToGlobalMappingRestoreInfo()` in Fortran. You must make sure that procs[], numprocs[] and 14766a818285SBarry Smith indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough. 14776a818285SBarry Smith 1478cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 1479db781477SPatrick Sanan `ISLocalToGlobalMappingRestoreInfo()` 14806a818285SBarry Smith @*/ 1481d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[]) 1482d71ae5a4SJacob Faibussowitsch { 14838a1f772fSStefano Zampini PetscInt **bindices = NULL, *bnumprocs = NULL, bs, i, j, k; 14846a818285SBarry Smith 14856a818285SBarry Smith PetscFunctionBegin; 14866a818285SBarry Smith PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 14878a1f772fSStefano Zampini bs = mapping->bs; 14889566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockInfo(mapping, nproc, procs, &bnumprocs, &bindices)); 1489268a049cSStefano Zampini if (bs > 1) { /* we need to expand the cached info */ 14909566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(*nproc, &*indices)); 14919566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(*nproc, &*numprocs)); 14926a818285SBarry Smith for (i = 0; i < *nproc; i++) { 14939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * bnumprocs[i], &(*indices)[i])); 1494268a049cSStefano Zampini for (j = 0; j < bnumprocs[i]; j++) { 1495ad540459SPierre Jolivet for (k = 0; k < bs; k++) (*indices)[i][j * bs + k] = bs * bindices[i][j] + k; 14966a818285SBarry Smith } 1497268a049cSStefano Zampini (*numprocs)[i] = bnumprocs[i] * bs; 14986a818285SBarry Smith } 1499268a049cSStefano Zampini mapping->info_free = PETSC_TRUE; 1500268a049cSStefano Zampini } else { 1501268a049cSStefano Zampini *numprocs = bnumprocs; 1502268a049cSStefano Zampini *indices = bindices; 15036a818285SBarry Smith } 15046a818285SBarry Smith PetscFunctionReturn(0); 15056a818285SBarry Smith } 15066a818285SBarry Smith 150707b52d57SBarry Smith /*@C 1508cab54364SBarry Smith ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by `ISLocalToGlobalMappingGetInfo()` 150989d82c54SBarry Smith 1510c3339decSBarry Smith Collective 151107b52d57SBarry Smith 1512f899ff85SJose E. Roman Input Parameter: 151307b52d57SBarry Smith . mapping - the mapping from local to global indexing 151407b52d57SBarry Smith 1515d8d19677SJose E. Roman Output Parameters: 151607b52d57SBarry Smith + nproc - number of processors that are connected to this one 151707b52d57SBarry Smith . proc - neighboring processors 151807b52d57SBarry Smith . numproc - number of indices for each processor 151907b52d57SBarry Smith - indices - indices of local nodes shared with neighbor (sorted by global numbering) 152007b52d57SBarry Smith 152107b52d57SBarry Smith Level: advanced 152207b52d57SBarry Smith 1523cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 1524db781477SPatrick Sanan `ISLocalToGlobalMappingGetInfo()` 152507b52d57SBarry Smith @*/ 1526d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[]) 1527d71ae5a4SJacob Faibussowitsch { 152807b52d57SBarry Smith PetscFunctionBegin; 15299566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockInfo(mapping, nproc, procs, numprocs, indices)); 153007b52d57SBarry Smith PetscFunctionReturn(0); 153107b52d57SBarry Smith } 153286994e45SJed Brown 153386994e45SJed Brown /*@C 1534cab54364SBarry Smith ISLocalToGlobalMappingGetNodeInfo - Gets the neighbor information for each MPI rank 15351bd0b88eSStefano Zampini 1536c3339decSBarry Smith Collective 15371bd0b88eSStefano Zampini 1538f899ff85SJose E. Roman Input Parameter: 15391bd0b88eSStefano Zampini . mapping - the mapping from local to global indexing 15401bd0b88eSStefano Zampini 1541d8d19677SJose E. Roman Output Parameters: 1542cab54364SBarry Smith + nnodes - number of local nodes (same `ISLocalToGlobalMappingGetSize()`) 15431bd0b88eSStefano Zampini . count - number of neighboring processors per node 15441bd0b88eSStefano Zampini - indices - indices of processes sharing the node (sorted) 15451bd0b88eSStefano Zampini 15461bd0b88eSStefano Zampini Level: advanced 15471bd0b88eSStefano Zampini 1548cab54364SBarry Smith Note: 1549cab54364SBarry Smith The user needs to call `ISLocalToGlobalMappingRestoreInfo()` when the data is no longer needed. 15501bd0b88eSStefano Zampini 1551cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 1552db781477SPatrick Sanan `ISLocalToGlobalMappingGetInfo()`, `ISLocalToGlobalMappingRestoreNodeInfo()` 15531bd0b88eSStefano Zampini @*/ 1554d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetNodeInfo(ISLocalToGlobalMapping mapping, PetscInt *nnodes, PetscInt *count[], PetscInt **indices[]) 1555d71ae5a4SJacob Faibussowitsch { 15561bd0b88eSStefano Zampini PetscInt n; 15571bd0b88eSStefano Zampini 15581bd0b88eSStefano Zampini PetscFunctionBegin; 15591bd0b88eSStefano Zampini PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 15609566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(mapping, &n)); 15611bd0b88eSStefano Zampini if (!mapping->info_nodec) { 15621bd0b88eSStefano Zampini PetscInt i, m, n_neigh, *neigh, *n_shared, **shared; 15631bd0b88eSStefano Zampini 15649566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n + 1, &mapping->info_nodec, n, &mapping->info_nodei)); 15659566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetInfo(mapping, &n_neigh, &neigh, &n_shared, &shared)); 1566ad540459SPierre Jolivet for (i = 0; i < n; i++) mapping->info_nodec[i] = 1; 1567071fcb05SBarry Smith m = n; 1568071fcb05SBarry Smith mapping->info_nodec[n] = 0; 15691bd0b88eSStefano Zampini for (i = 1; i < n_neigh; i++) { 15701bd0b88eSStefano Zampini PetscInt j; 15711bd0b88eSStefano Zampini 15721bd0b88eSStefano Zampini m += n_shared[i]; 15731bd0b88eSStefano Zampini for (j = 0; j < n_shared[i]; j++) mapping->info_nodec[shared[i][j]] += 1; 15741bd0b88eSStefano Zampini } 15759566063dSJacob Faibussowitsch if (n) PetscCall(PetscMalloc1(m, &mapping->info_nodei[0])); 15761bd0b88eSStefano Zampini for (i = 1; i < n; i++) mapping->info_nodei[i] = mapping->info_nodei[i - 1] + mapping->info_nodec[i - 1]; 15779566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(mapping->info_nodec, n)); 15789371c9d4SSatish Balay for (i = 0; i < n; i++) { 15799371c9d4SSatish Balay mapping->info_nodec[i] = 1; 15809371c9d4SSatish Balay mapping->info_nodei[i][0] = neigh[0]; 15819371c9d4SSatish Balay } 15821bd0b88eSStefano Zampini for (i = 1; i < n_neigh; i++) { 15831bd0b88eSStefano Zampini PetscInt j; 15841bd0b88eSStefano Zampini 15851bd0b88eSStefano Zampini for (j = 0; j < n_shared[i]; j++) { 15861bd0b88eSStefano Zampini PetscInt k = shared[i][j]; 15871bd0b88eSStefano Zampini 15881bd0b88eSStefano Zampini mapping->info_nodei[k][mapping->info_nodec[k]] = neigh[i]; 15891bd0b88eSStefano Zampini mapping->info_nodec[k] += 1; 15901bd0b88eSStefano Zampini } 15911bd0b88eSStefano Zampini } 15929566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(PetscSortRemoveDupsInt(&mapping->info_nodec[i], mapping->info_nodei[i])); 15939566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreInfo(mapping, &n_neigh, &neigh, &n_shared, &shared)); 15941bd0b88eSStefano Zampini } 15951bd0b88eSStefano Zampini if (nnodes) *nnodes = n; 15961bd0b88eSStefano Zampini if (count) *count = mapping->info_nodec; 15971bd0b88eSStefano Zampini if (indices) *indices = mapping->info_nodei; 15981bd0b88eSStefano Zampini PetscFunctionReturn(0); 15991bd0b88eSStefano Zampini } 16001bd0b88eSStefano Zampini 16011bd0b88eSStefano Zampini /*@C 1602cab54364SBarry Smith ISLocalToGlobalMappingRestoreNodeInfo - Frees the memory allocated by `ISLocalToGlobalMappingGetNodeInfo()` 16031bd0b88eSStefano Zampini 1604c3339decSBarry Smith Collective 16051bd0b88eSStefano Zampini 1606f899ff85SJose E. Roman Input Parameter: 16071bd0b88eSStefano Zampini . mapping - the mapping from local to global indexing 16081bd0b88eSStefano Zampini 1609d8d19677SJose E. Roman Output Parameters: 16101bd0b88eSStefano Zampini + nnodes - number of local nodes 16111bd0b88eSStefano Zampini . count - number of neighboring processors per node 16121bd0b88eSStefano Zampini - indices - indices of processes sharing the node (sorted) 16131bd0b88eSStefano Zampini 16141bd0b88eSStefano Zampini Level: advanced 16151bd0b88eSStefano Zampini 1616cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 1617db781477SPatrick Sanan `ISLocalToGlobalMappingGetInfo()` 16181bd0b88eSStefano Zampini @*/ 1619d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreNodeInfo(ISLocalToGlobalMapping mapping, PetscInt *nnodes, PetscInt *count[], PetscInt **indices[]) 1620d71ae5a4SJacob Faibussowitsch { 16211bd0b88eSStefano Zampini PetscFunctionBegin; 16221bd0b88eSStefano Zampini PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 16231bd0b88eSStefano Zampini if (nnodes) *nnodes = 0; 16241bd0b88eSStefano Zampini if (count) *count = NULL; 16251bd0b88eSStefano Zampini if (indices) *indices = NULL; 16261bd0b88eSStefano Zampini PetscFunctionReturn(0); 16271bd0b88eSStefano Zampini } 16281bd0b88eSStefano Zampini 1629*6ce40d10SBarry Smith /*MC 1630*6ce40d10SBarry Smith ISLocalToGlobalMappingGetIndicesF90 - Get global indices for every local point that is mapped 1631*6ce40d10SBarry Smith 1632*6ce40d10SBarry Smith Synopsis: 1633*6ce40d10SBarry Smith ISLocalToGlobalMappingGetIndicesF90(ISLocalToGlobalMapping ltog, PetscInt, pointer :: array(:)}, integer ierr) 1634*6ce40d10SBarry Smith 1635*6ce40d10SBarry Smith Not Collective 1636*6ce40d10SBarry Smith 1637*6ce40d10SBarry Smith Input Parameter: 1638*6ce40d10SBarry Smith . A - the matrix 1639*6ce40d10SBarry Smith 1640*6ce40d10SBarry Smith Output Parameters: 1641*6ce40d10SBarry Smith . array - array of indices, the length of this array may be obtained with `ISLocalToGlobalMappingGetSize()` 1642*6ce40d10SBarry Smith 1643*6ce40d10SBarry Smith Level: advanced 1644*6ce40d10SBarry Smith 1645*6ce40d10SBarry Smith Note: 1646*6ce40d10SBarry Smith Use `ISLocalToGlobalMappingGetIndicesF90()` when you no longer need access to the data 1647*6ce40d10SBarry Smith 1648*6ce40d10SBarry Smith .seealso: [](sec_fortranarrays), [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingGetIndices()`, `ISLocalToGlobalMappingRestoreIndices()`, `ISLocalToGlobalMappingRestoreIndicesF90()` 1649*6ce40d10SBarry Smith M*/ 1650*6ce40d10SBarry Smith 1651*6ce40d10SBarry Smith /*MC 1652*6ce40d10SBarry Smith ISLocalToGlobalMappingRestoreIndicesF90 - restores the global indices for every local point that is mapped obtained with `ISLocalToGlobalMappingGetIndicesF90()` 1653*6ce40d10SBarry Smith 1654*6ce40d10SBarry Smith Synopsis: 1655*6ce40d10SBarry Smith ISLocalToGlobalMappingRestoreIndicesF90(ISLocalToGlobalMapping ltog, PetscInt, pointer :: array(:)}, integer ierr) 1656*6ce40d10SBarry Smith 1657*6ce40d10SBarry Smith Not Collective 1658*6ce40d10SBarry Smith 1659*6ce40d10SBarry Smith Input Parameters: 1660*6ce40d10SBarry Smith + A - the matrix 1661*6ce40d10SBarry Smith - array - array of indices, the length of this array may be obtained with `ISLocalToGlobalMappingGetSize()` 1662*6ce40d10SBarry Smith 1663*6ce40d10SBarry Smith Level: advanced 1664*6ce40d10SBarry Smith 1665*6ce40d10SBarry Smith .seealso: [](sec_fortranarrays), [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingGetIndices()`, `ISLocalToGlobalMappingRestoreIndices()`, `ISLocalToGlobalMappingGetIndicesF90()` 1666*6ce40d10SBarry Smith M*/ 1667*6ce40d10SBarry Smith 16681bd0b88eSStefano Zampini /*@C 1669107e9a97SBarry Smith ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped 167086994e45SJed Brown 167186994e45SJed Brown Not Collective 167286994e45SJed Brown 16734165533cSJose E. Roman Input Parameter: 167486994e45SJed Brown . ltog - local to global mapping 167586994e45SJed Brown 16764165533cSJose E. Roman Output Parameter: 1677cab54364SBarry Smith . array - array of indices, the length of this array may be obtained with `ISLocalToGlobalMappingGetSize()` 167886994e45SJed Brown 167986994e45SJed Brown Level: advanced 168086994e45SJed Brown 1681cab54364SBarry Smith Note: 1682cab54364SBarry Smith `ISLocalToGlobalMappingGetSize()` returns the length the this array 1683107e9a97SBarry Smith 1684cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingRestoreIndices()`, `ISLocalToGlobalMappingGetBlockIndices()`, `ISLocalToGlobalMappingRestoreBlockIndices()` 168586994e45SJed Brown @*/ 1686d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog, const PetscInt **array) 1687d71ae5a4SJacob Faibussowitsch { 168886994e45SJed Brown PetscFunctionBegin; 168986994e45SJed Brown PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1); 169086994e45SJed Brown PetscValidPointer(array, 2); 169145b6f7e9SBarry Smith if (ltog->bs == 1) { 169286994e45SJed Brown *array = ltog->indices; 169345b6f7e9SBarry Smith } else { 169445b6f7e9SBarry Smith PetscInt *jj, k, i, j, n = ltog->n, bs = ltog->bs; 169545b6f7e9SBarry Smith const PetscInt *ii; 169645b6f7e9SBarry Smith 16979566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * n, &jj)); 169845b6f7e9SBarry Smith *array = jj; 169945b6f7e9SBarry Smith k = 0; 170045b6f7e9SBarry Smith ii = ltog->indices; 170145b6f7e9SBarry Smith for (i = 0; i < n; i++) 17029371c9d4SSatish Balay for (j = 0; j < bs; j++) jj[k++] = bs * ii[i] + j; 170345b6f7e9SBarry Smith } 170486994e45SJed Brown PetscFunctionReturn(0); 170586994e45SJed Brown } 170686994e45SJed Brown 170786994e45SJed Brown /*@C 1708cab54364SBarry Smith ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with `ISLocalToGlobalMappingGetIndices()` 170986994e45SJed Brown 171086994e45SJed Brown Not Collective 171186994e45SJed Brown 17124165533cSJose E. Roman Input Parameters: 171386994e45SJed Brown + ltog - local to global mapping 171486994e45SJed Brown - array - array of indices 171586994e45SJed Brown 171686994e45SJed Brown Level: advanced 171786994e45SJed Brown 1718cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingGetIndices()` 171986994e45SJed Brown @*/ 1720d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog, const PetscInt **array) 1721d71ae5a4SJacob Faibussowitsch { 172286994e45SJed Brown PetscFunctionBegin; 172386994e45SJed Brown PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1); 172486994e45SJed Brown PetscValidPointer(array, 2); 1725c9cc58a2SBarry Smith PetscCheck(ltog->bs != 1 || *array == ltog->indices, PETSC_COMM_SELF, PETSC_ERR_ARG_BADPTR, "Trying to return mismatched pointer"); 172645b6f7e9SBarry Smith 172748a46eb9SPierre Jolivet if (ltog->bs > 1) PetscCall(PetscFree(*(void **)array)); 172845b6f7e9SBarry Smith PetscFunctionReturn(0); 172945b6f7e9SBarry Smith } 173045b6f7e9SBarry Smith 173145b6f7e9SBarry Smith /*@C 173245b6f7e9SBarry Smith ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block 173345b6f7e9SBarry Smith 173445b6f7e9SBarry Smith Not Collective 173545b6f7e9SBarry Smith 17364165533cSJose E. Roman Input Parameter: 173745b6f7e9SBarry Smith . ltog - local to global mapping 173845b6f7e9SBarry Smith 17394165533cSJose E. Roman Output Parameter: 174045b6f7e9SBarry Smith . array - array of indices 174145b6f7e9SBarry Smith 174245b6f7e9SBarry Smith Level: advanced 174345b6f7e9SBarry Smith 1744cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingRestoreBlockIndices()` 174545b6f7e9SBarry Smith @*/ 1746d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog, const PetscInt **array) 1747d71ae5a4SJacob Faibussowitsch { 174845b6f7e9SBarry Smith PetscFunctionBegin; 174945b6f7e9SBarry Smith PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1); 175045b6f7e9SBarry Smith PetscValidPointer(array, 2); 175145b6f7e9SBarry Smith *array = ltog->indices; 175245b6f7e9SBarry Smith PetscFunctionReturn(0); 175345b6f7e9SBarry Smith } 175445b6f7e9SBarry Smith 175545b6f7e9SBarry Smith /*@C 1756cab54364SBarry Smith ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with `ISLocalToGlobalMappingGetBlockIndices()` 175745b6f7e9SBarry Smith 175845b6f7e9SBarry Smith Not Collective 175945b6f7e9SBarry Smith 17604165533cSJose E. Roman Input Parameters: 176145b6f7e9SBarry Smith + ltog - local to global mapping 176245b6f7e9SBarry Smith - array - array of indices 176345b6f7e9SBarry Smith 176445b6f7e9SBarry Smith Level: advanced 176545b6f7e9SBarry Smith 1766cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingGetIndices()` 176745b6f7e9SBarry Smith @*/ 1768d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog, const PetscInt **array) 1769d71ae5a4SJacob Faibussowitsch { 177045b6f7e9SBarry Smith PetscFunctionBegin; 177145b6f7e9SBarry Smith PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1); 177245b6f7e9SBarry Smith PetscValidPointer(array, 2); 177308401ef6SPierre Jolivet PetscCheck(*array == ltog->indices, PETSC_COMM_SELF, PETSC_ERR_ARG_BADPTR, "Trying to return mismatched pointer"); 17740298fd71SBarry Smith *array = NULL; 177586994e45SJed Brown PetscFunctionReturn(0); 177686994e45SJed Brown } 1777f7efa3c7SJed Brown 1778f7efa3c7SJed Brown /*@C 1779f7efa3c7SJed Brown ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings 1780f7efa3c7SJed Brown 1781f7efa3c7SJed Brown Not Collective 1782f7efa3c7SJed Brown 17834165533cSJose E. Roman Input Parameters: 1784f7efa3c7SJed Brown + comm - communicator for the new mapping, must contain the communicator of every mapping to concatenate 1785f7efa3c7SJed Brown . n - number of mappings to concatenate 1786f7efa3c7SJed Brown - ltogs - local to global mappings 1787f7efa3c7SJed Brown 17884165533cSJose E. Roman Output Parameter: 1789f7efa3c7SJed Brown . ltogcat - new mapping 1790f7efa3c7SJed Brown 1791f7efa3c7SJed Brown Level: advanced 1792f7efa3c7SJed Brown 1793cab54364SBarry Smith Note: 1794cab54364SBarry Smith This currently always returns a mapping with block size of 1 1795cab54364SBarry Smith 1796cab54364SBarry Smith Developer Note: 1797cab54364SBarry Smith If all the input mapping have the same block size we could easily handle that as a special case 1798cab54364SBarry Smith 1799cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingCreate()` 1800f7efa3c7SJed Brown @*/ 1801d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm, PetscInt n, const ISLocalToGlobalMapping ltogs[], ISLocalToGlobalMapping *ltogcat) 1802d71ae5a4SJacob Faibussowitsch { 1803f7efa3c7SJed Brown PetscInt i, cnt, m, *idx; 1804f7efa3c7SJed Brown 1805f7efa3c7SJed Brown PetscFunctionBegin; 180608401ef6SPierre Jolivet PetscCheck(n >= 0, comm, PETSC_ERR_ARG_OUTOFRANGE, "Must have a non-negative number of mappings, given %" PetscInt_FMT, n); 1807f7efa3c7SJed Brown if (n > 0) PetscValidPointer(ltogs, 3); 1808f7efa3c7SJed Brown for (i = 0; i < n; i++) PetscValidHeaderSpecific(ltogs[i], IS_LTOGM_CLASSID, 3); 1809f7efa3c7SJed Brown PetscValidPointer(ltogcat, 4); 1810f7efa3c7SJed Brown for (cnt = 0, i = 0; i < n; i++) { 18119566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(ltogs[i], &m)); 1812f7efa3c7SJed Brown cnt += m; 1813f7efa3c7SJed Brown } 18149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(cnt, &idx)); 1815f7efa3c7SJed Brown for (cnt = 0, i = 0; i < n; i++) { 1816f7efa3c7SJed Brown const PetscInt *subidx; 18179566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(ltogs[i], &m)); 18189566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(ltogs[i], &subidx)); 18199566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(&idx[cnt], subidx, m)); 18209566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(ltogs[i], &subidx)); 1821f7efa3c7SJed Brown cnt += m; 1822f7efa3c7SJed Brown } 18239566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm, 1, cnt, idx, PETSC_OWN_POINTER, ltogcat)); 1824f7efa3c7SJed Brown PetscFunctionReturn(0); 1825f7efa3c7SJed Brown } 182604a59952SBarry Smith 1827413f72f0SBarry Smith /*MC 1828cab54364SBarry Smith ISLOCALTOGLOBALMAPPINGBASIC - basic implementation of the `ISLocalToGlobalMapping` object. When `ISGlobalToLocalMappingApply()` is 1829413f72f0SBarry Smith used this is good for only small and moderate size problems. 1830413f72f0SBarry Smith 1831413f72f0SBarry Smith Options Database Keys: 1832a2b725a8SWilliam Gropp . -islocaltoglobalmapping_type basic - select this method 1833413f72f0SBarry Smith 1834413f72f0SBarry Smith Level: beginner 1835413f72f0SBarry Smith 1836cab54364SBarry Smith Developer Note: 1837cab54364SBarry Smith This stores all the mapping information on each MPI rank. 1838cab54364SBarry Smith 1839cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`, `ISLOCALTOGLOBALMAPPINGHASH` 1840413f72f0SBarry Smith M*/ 1841d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Basic(ISLocalToGlobalMapping ltog) 1842d71ae5a4SJacob Faibussowitsch { 1843413f72f0SBarry Smith PetscFunctionBegin; 1844413f72f0SBarry Smith ltog->ops->globaltolocalmappingapply = ISGlobalToLocalMappingApply_Basic; 1845413f72f0SBarry Smith ltog->ops->globaltolocalmappingsetup = ISGlobalToLocalMappingSetUp_Basic; 1846413f72f0SBarry Smith ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Basic; 1847413f72f0SBarry Smith ltog->ops->destroy = ISLocalToGlobalMappingDestroy_Basic; 1848413f72f0SBarry Smith PetscFunctionReturn(0); 1849413f72f0SBarry Smith } 1850413f72f0SBarry Smith 1851413f72f0SBarry Smith /*MC 1852cab54364SBarry Smith ISLOCALTOGLOBALMAPPINGHASH - hash implementation of the `ISLocalToGlobalMapping` object. When `ISGlobalToLocalMappingApply()` is 1853ed56e8eaSBarry Smith used this is good for large memory problems. 1854413f72f0SBarry Smith 1855413f72f0SBarry Smith Options Database Keys: 1856a2b725a8SWilliam Gropp . -islocaltoglobalmapping_type hash - select this method 1857413f72f0SBarry Smith 1858413f72f0SBarry Smith Level: beginner 1859413f72f0SBarry Smith 1860cab54364SBarry Smith Note: 1861cab54364SBarry Smith This is selected automatically for large problems if the user does not set the type. 1862cab54364SBarry Smith 1863cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`, `ISLOCALTOGLOBALMAPPINGBASIC` 1864413f72f0SBarry Smith M*/ 1865d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Hash(ISLocalToGlobalMapping ltog) 1866d71ae5a4SJacob Faibussowitsch { 1867413f72f0SBarry Smith PetscFunctionBegin; 1868413f72f0SBarry Smith ltog->ops->globaltolocalmappingapply = ISGlobalToLocalMappingApply_Hash; 1869413f72f0SBarry Smith ltog->ops->globaltolocalmappingsetup = ISGlobalToLocalMappingSetUp_Hash; 1870413f72f0SBarry Smith ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Hash; 1871413f72f0SBarry Smith ltog->ops->destroy = ISLocalToGlobalMappingDestroy_Hash; 1872413f72f0SBarry Smith PetscFunctionReturn(0); 1873413f72f0SBarry Smith } 1874413f72f0SBarry Smith 1875413f72f0SBarry Smith /*@C 1876cab54364SBarry Smith ISLocalToGlobalMappingRegister - Registers a method for applying a global to local mapping with an `ISLocalToGlobalMapping` 1877413f72f0SBarry Smith 1878413f72f0SBarry Smith Not Collective 1879413f72f0SBarry Smith 1880413f72f0SBarry Smith Input Parameters: 1881413f72f0SBarry Smith + sname - name of a new method 1882413f72f0SBarry Smith - routine_create - routine to create method context 1883413f72f0SBarry Smith 1884413f72f0SBarry Smith Sample usage: 1885413f72f0SBarry Smith .vb 1886413f72f0SBarry Smith ISLocalToGlobalMappingRegister("my_mapper",MyCreate); 1887413f72f0SBarry Smith .ve 1888413f72f0SBarry Smith 1889ed56e8eaSBarry Smith Then, your mapping can be chosen with the procedural interface via 1890413f72f0SBarry Smith $ ISLocalToGlobalMappingSetType(ltog,"my_mapper") 1891413f72f0SBarry Smith or at runtime via the option 1892ed56e8eaSBarry Smith $ -islocaltoglobalmapping_type my_mapper 1893413f72f0SBarry Smith 1894413f72f0SBarry Smith Level: advanced 1895413f72f0SBarry Smith 1896cab54364SBarry Smith Note: 1897cab54364SBarry Smith `ISLocalToGlobalMappingRegister()` may be called multiple times to add several user-defined mappings. 1898413f72f0SBarry Smith 1899cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingRegisterAll()`, `ISLocalToGlobalMappingRegisterDestroy()`, `ISLOCALTOGLOBALMAPPINGBASIC`, 1900cab54364SBarry Smith `ISLOCALTOGLOBALMAPPINGHASH`, `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApply()` 1901413f72f0SBarry Smith @*/ 1902d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRegister(const char sname[], PetscErrorCode (*function)(ISLocalToGlobalMapping)) 1903d71ae5a4SJacob Faibussowitsch { 1904413f72f0SBarry Smith PetscFunctionBegin; 19059566063dSJacob Faibussowitsch PetscCall(ISInitializePackage()); 19069566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&ISLocalToGlobalMappingList, sname, function)); 1907413f72f0SBarry Smith PetscFunctionReturn(0); 1908413f72f0SBarry Smith } 1909413f72f0SBarry Smith 1910413f72f0SBarry Smith /*@C 1911cab54364SBarry Smith ISLocalToGlobalMappingSetType - Sets the implementation type `ISLocalToGlobalMapping` will use 1912413f72f0SBarry Smith 1913c3339decSBarry Smith Logically Collective 1914413f72f0SBarry Smith 1915413f72f0SBarry Smith Input Parameters: 1916cab54364SBarry Smith + ltog - the `ISLocalToGlobalMapping` object 1917413f72f0SBarry Smith - type - a known method 1918413f72f0SBarry Smith 1919413f72f0SBarry Smith Options Database Key: 1920cab54364SBarry Smith . -islocaltoglobalmapping_type <method> - Sets the method; use -help for a list of available methods (for instance, basic or hash) 1921cab54364SBarry Smith 1922cab54364SBarry Smith Level: intermediate 1923413f72f0SBarry Smith 1924413f72f0SBarry Smith Notes: 1925413f72f0SBarry Smith See "petsc/include/petscis.h" for available methods 1926413f72f0SBarry Smith 1927cab54364SBarry Smith Normally, it is best to use the `ISLocalToGlobalMappingSetFromOptions()` command and 1928cab54364SBarry Smith then set the `ISLocalToGlobalMappingType` from the options database rather than by using 1929413f72f0SBarry Smith this routine. 1930413f72f0SBarry Smith 1931cab54364SBarry Smith Developer Note: 1932cab54364SBarry Smith `ISLocalToGlobalMappingRegister()` is used to add new types to `ISLocalToGlobalMappingList` from which they 1933cab54364SBarry Smith are accessed by `ISLocalToGlobalMappingSetType()`. 1934413f72f0SBarry Smith 1935cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingRegister()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingGetType()` 1936413f72f0SBarry Smith @*/ 1937d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType type) 1938d71ae5a4SJacob Faibussowitsch { 1939413f72f0SBarry Smith PetscBool match; 19405f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(ISLocalToGlobalMapping) = NULL; 1941413f72f0SBarry Smith 1942413f72f0SBarry Smith PetscFunctionBegin; 1943413f72f0SBarry Smith PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1); 1944a0d79125SStefano Zampini if (type) PetscValidCharPointer(type, 2); 1945413f72f0SBarry Smith 19469566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)ltog, type, &match)); 1947413f72f0SBarry Smith if (match) PetscFunctionReturn(0); 1948413f72f0SBarry Smith 1949a0d79125SStefano Zampini /* L2G maps defer type setup at globaltolocal calls, allow passing NULL here */ 1950a0d79125SStefano Zampini if (type) { 19519566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(ISLocalToGlobalMappingList, type, &r)); 1952a0d79125SStefano Zampini PetscCheck(r, PetscObjectComm((PetscObject)ltog), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested ISLocalToGlobalMapping type %s", type); 1953a0d79125SStefano Zampini } 1954413f72f0SBarry Smith /* Destroy the previous private LTOG context */ 1955dbbe0bcdSBarry Smith PetscTryTypeMethod(ltog, destroy); 1956413f72f0SBarry Smith ltog->ops->destroy = NULL; 1957dbbe0bcdSBarry Smith 19589566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)ltog, type)); 19599566063dSJacob Faibussowitsch if (r) PetscCall((*r)(ltog)); 1960a0d79125SStefano Zampini PetscFunctionReturn(0); 1961a0d79125SStefano Zampini } 1962a0d79125SStefano Zampini 1963a0d79125SStefano Zampini /*@C 1964cab54364SBarry Smith ISLocalToGlobalMappingGetType - Get the type of the `ISLocalToGlobalMapping` 1965a0d79125SStefano Zampini 1966a0d79125SStefano Zampini Not Collective 1967a0d79125SStefano Zampini 1968a0d79125SStefano Zampini Input Parameter: 1969cab54364SBarry Smith . ltog - the `ISLocalToGlobalMapping` object 1970a0d79125SStefano Zampini 1971a0d79125SStefano Zampini Output Parameter: 1972a0d79125SStefano Zampini . type - the type 1973a0d79125SStefano Zampini 197449762cbcSSatish Balay Level: intermediate 197549762cbcSSatish Balay 1976cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingRegister()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()` 1977a0d79125SStefano Zampini @*/ 1978d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType *type) 1979d71ae5a4SJacob Faibussowitsch { 1980a0d79125SStefano Zampini PetscFunctionBegin; 1981a0d79125SStefano Zampini PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1); 1982a0d79125SStefano Zampini PetscValidPointer(type, 2); 1983a0d79125SStefano Zampini *type = ((PetscObject)ltog)->type_name; 1984413f72f0SBarry Smith PetscFunctionReturn(0); 1985413f72f0SBarry Smith } 1986413f72f0SBarry Smith 1987413f72f0SBarry Smith PetscBool ISLocalToGlobalMappingRegisterAllCalled = PETSC_FALSE; 1988413f72f0SBarry Smith 1989413f72f0SBarry Smith /*@C 1990cab54364SBarry Smith ISLocalToGlobalMappingRegisterAll - Registers all of the local to global mapping components in the `IS` package. 1991413f72f0SBarry Smith 1992413f72f0SBarry Smith Not Collective 1993413f72f0SBarry Smith 1994413f72f0SBarry Smith Level: advanced 1995413f72f0SBarry Smith 1996cab54364SBarry Smith .seealso: [](sec_scatter), `ISRegister()`, `ISLocalToGlobalRegister()` 1997413f72f0SBarry Smith @*/ 1998d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRegisterAll(void) 1999d71ae5a4SJacob Faibussowitsch { 2000413f72f0SBarry Smith PetscFunctionBegin; 2001413f72f0SBarry Smith if (ISLocalToGlobalMappingRegisterAllCalled) PetscFunctionReturn(0); 2002413f72f0SBarry Smith ISLocalToGlobalMappingRegisterAllCalled = PETSC_TRUE; 20039566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGBASIC, ISLocalToGlobalMappingCreate_Basic)); 20049566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGHASH, ISLocalToGlobalMappingCreate_Hash)); 2005413f72f0SBarry Smith PetscFunctionReturn(0); 2006413f72f0SBarry Smith } 2007