1af0996ceSBarry Smith #include <petsc/private/isimpl.h> /*I "petscis.h" I*/ 2e8f14785SLisandro Dalcin #include <petsc/private/hashmapi.h> 30c312b8eSJed Brown #include <petscsf.h> 4665c2dedSJed Brown #include <petscviewer.h> 52362add9SBarry Smith 67087cfbeSBarry Smith PetscClassId IS_LTOGM_CLASSID; 7268a049cSStefano Zampini static PetscErrorCode ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping, PetscInt *, PetscInt **, PetscInt **, PetscInt ***); 88e58c17dSMatthew Knepley 9413f72f0SBarry Smith typedef struct { 10413f72f0SBarry Smith PetscInt *globals; 11413f72f0SBarry Smith } ISLocalToGlobalMapping_Basic; 12413f72f0SBarry Smith 13413f72f0SBarry Smith typedef struct { 14e8f14785SLisandro Dalcin PetscHMapI globalht; 15413f72f0SBarry Smith } ISLocalToGlobalMapping_Hash; 16413f72f0SBarry Smith 176528b96dSMatthew G. Knepley /*@C 18cab54364SBarry Smith ISGetPointRange - Returns a description of the points in an `IS` suitable for traversal 19413f72f0SBarry Smith 2020662ed9SBarry Smith Not Collective 216528b96dSMatthew G. Knepley 226528b96dSMatthew G. Knepley Input Parameter: 23cab54364SBarry Smith . pointIS - The `IS` object 246528b96dSMatthew G. Knepley 256528b96dSMatthew G. Knepley Output Parameters: 266528b96dSMatthew G. Knepley + pStart - The first index, see notes 276528b96dSMatthew G. Knepley . pEnd - One past the last index, see notes 286528b96dSMatthew G. Knepley - points - The indices, see notes 296528b96dSMatthew G. Knepley 306528b96dSMatthew G. Knepley Level: intermediate 316528b96dSMatthew G. Knepley 32cab54364SBarry Smith Notes: 3320662ed9SBarry Smith If the `IS` contains contiguous indices in an `ISSTRIDE`, then the indices are contained in [pStart, pEnd) and points = `NULL`. 3420662ed9SBarry Smith 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; 3920662ed9SBarry Smith // use point 40cab54364SBarry Smith } 41cab54364SBarry Smith ISRestorePointRange(is, &pstart, &pEnd, &points); 42cab54364SBarry Smith .ve 4320662ed9SBarry Smith Hence the same code can be written for `pointIS` being a `ISSTRIDE` or `ISGENERAL` 44cab54364SBarry Smith 45cab54364SBarry Smith .seealso: [](sec_scatter), `IS`, `ISRestorePointRange()`, `ISGetPointSubrange()`, `ISGetIndices()`, `ISCreateStride()` 466528b96dSMatthew G. Knepley @*/ 47d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGetPointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points) 48d71ae5a4SJacob Faibussowitsch { 499305a4c7SMatthew G. Knepley PetscInt numCells, step = 1; 509305a4c7SMatthew G. Knepley PetscBool isStride; 519305a4c7SMatthew G. Knepley 529305a4c7SMatthew G. Knepley PetscFunctionBeginHot; 539305a4c7SMatthew G. Knepley *pStart = 0; 549305a4c7SMatthew G. Knepley *points = NULL; 559566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(pointIS, &numCells)); 569566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pointIS, ISSTRIDE, &isStride)); 579566063dSJacob Faibussowitsch if (isStride) PetscCall(ISStrideGetInfo(pointIS, pStart, &step)); 589305a4c7SMatthew G. Knepley *pEnd = *pStart + numCells; 599566063dSJacob Faibussowitsch if (!isStride || step != 1) PetscCall(ISGetIndices(pointIS, points)); 603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 619305a4c7SMatthew G. Knepley } 629305a4c7SMatthew G. Knepley 636528b96dSMatthew G. Knepley /*@C 64cab54364SBarry Smith ISRestorePointRange - Destroys the traversal description created with `ISGetPointRange()` 656528b96dSMatthew G. Knepley 6620662ed9SBarry Smith Not Collective 676528b96dSMatthew G. Knepley 686528b96dSMatthew G. Knepley Input Parameters: 69cab54364SBarry Smith + pointIS - The `IS` object 70cab54364SBarry Smith . pStart - The first index, from `ISGetPointRange()` 71cab54364SBarry Smith . pEnd - One past the last index, from `ISGetPointRange()` 72cab54364SBarry Smith - points - The indices, from `ISGetPointRange()` 736528b96dSMatthew G. Knepley 746528b96dSMatthew G. Knepley Level: intermediate 756528b96dSMatthew G. Knepley 76cab54364SBarry Smith .seealso: [](sec_scatter), `IS`, `ISGetPointRange()`, `ISGetPointSubrange()`, `ISGetIndices()`, `ISCreateStride()` 776528b96dSMatthew G. Knepley @*/ 78d71ae5a4SJacob Faibussowitsch PetscErrorCode ISRestorePointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points) 79d71ae5a4SJacob Faibussowitsch { 809305a4c7SMatthew G. Knepley PetscInt step = 1; 819305a4c7SMatthew G. Knepley PetscBool isStride; 829305a4c7SMatthew G. Knepley 839305a4c7SMatthew G. Knepley PetscFunctionBeginHot; 849566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pointIS, ISSTRIDE, &isStride)); 859566063dSJacob Faibussowitsch if (isStride) PetscCall(ISStrideGetInfo(pointIS, pStart, &step)); 869566063dSJacob Faibussowitsch if (!isStride || step != 1) PetscCall(ISGetIndices(pointIS, points)); 873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 889305a4c7SMatthew G. Knepley } 899305a4c7SMatthew G. Knepley 906528b96dSMatthew G. Knepley /*@C 91cab54364SBarry Smith ISGetPointSubrange - Configures the input `IS` to be a subrange for the traversal information given 926528b96dSMatthew G. Knepley 9320662ed9SBarry Smith Not Collective 946528b96dSMatthew G. Knepley 956528b96dSMatthew G. Knepley Input Parameters: 96cab54364SBarry Smith + subpointIS - The `IS` object to be configured 976528b96dSMatthew G. Knepley . pStart - The first index of the subrange 986528b96dSMatthew G. Knepley . pEnd - One past the last index for the subrange 99cab54364SBarry Smith - points - The indices for the entire range, from `ISGetPointRange()` 1006528b96dSMatthew G. Knepley 1016528b96dSMatthew G. Knepley Output Parameters: 102cab54364SBarry Smith . subpointIS - The `IS` object now configured to be a subrange 1036528b96dSMatthew G. Knepley 1046528b96dSMatthew G. Knepley Level: intermediate 1056528b96dSMatthew G. Knepley 106cab54364SBarry Smith Note: 107cab54364SBarry Smith The input `IS` will now respond properly to calls to `ISGetPointRange()` and return the subrange. 108cab54364SBarry Smith 109cab54364SBarry Smith .seealso: [](sec_scatter), `IS`, `ISGetPointRange()`, `ISRestorePointRange()`, `ISGetIndices()`, `ISCreateStride()` 1106528b96dSMatthew G. Knepley @*/ 111d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGetPointSubrange(IS subpointIS, PetscInt pStart, PetscInt pEnd, const PetscInt *points) 112d71ae5a4SJacob Faibussowitsch { 1139305a4c7SMatthew G. Knepley PetscFunctionBeginHot; 1149305a4c7SMatthew G. Knepley if (points) { 1159566063dSJacob Faibussowitsch PetscCall(ISSetType(subpointIS, ISGENERAL)); 1169566063dSJacob Faibussowitsch PetscCall(ISGeneralSetIndices(subpointIS, pEnd - pStart, &points[pStart], PETSC_USE_POINTER)); 1179305a4c7SMatthew G. Knepley } else { 1189566063dSJacob Faibussowitsch PetscCall(ISSetType(subpointIS, ISSTRIDE)); 1199566063dSJacob Faibussowitsch PetscCall(ISStrideSetStride(subpointIS, pEnd - pStart, pStart, 1)); 1209305a4c7SMatthew G. Knepley } 1213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1229305a4c7SMatthew G. Knepley } 1239305a4c7SMatthew G. Knepley 124413f72f0SBarry Smith /* -----------------------------------------------------------------------------------------*/ 125413f72f0SBarry Smith 126413f72f0SBarry Smith /* 127413f72f0SBarry Smith Creates the global mapping information in the ISLocalToGlobalMapping structure 128413f72f0SBarry Smith 129413f72f0SBarry Smith If the user has not selected how to handle the global to local mapping then use HASH for "large" problems 130413f72f0SBarry Smith */ 131d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISGlobalToLocalMappingSetUp(ISLocalToGlobalMapping mapping) 132d71ae5a4SJacob Faibussowitsch { 133413f72f0SBarry Smith PetscInt i, *idx = mapping->indices, n = mapping->n, end, start; 134413f72f0SBarry Smith 135413f72f0SBarry Smith PetscFunctionBegin; 1363ba16761SJacob Faibussowitsch if (mapping->data) PetscFunctionReturn(PETSC_SUCCESS); 137413f72f0SBarry Smith end = 0; 138413f72f0SBarry Smith start = PETSC_MAX_INT; 139413f72f0SBarry Smith 140413f72f0SBarry Smith for (i = 0; i < n; i++) { 141413f72f0SBarry Smith if (idx[i] < 0) continue; 142413f72f0SBarry Smith if (idx[i] < start) start = idx[i]; 143413f72f0SBarry Smith if (idx[i] > end) end = idx[i]; 144413f72f0SBarry Smith } 1459371c9d4SSatish Balay if (start > end) { 1469371c9d4SSatish Balay start = 0; 1479371c9d4SSatish Balay end = -1; 1489371c9d4SSatish Balay } 149413f72f0SBarry Smith mapping->globalstart = start; 150413f72f0SBarry Smith mapping->globalend = end; 151413f72f0SBarry Smith if (!((PetscObject)mapping)->type_name) { 152413f72f0SBarry Smith if ((end - start) > PetscMax(4 * n, 1000000)) { 1539566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingSetType(mapping, ISLOCALTOGLOBALMAPPINGHASH)); 154413f72f0SBarry Smith } else { 1559566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingSetType(mapping, ISLOCALTOGLOBALMAPPINGBASIC)); 156413f72f0SBarry Smith } 157413f72f0SBarry Smith } 158dbbe0bcdSBarry Smith PetscTryTypeMethod(mapping, globaltolocalmappingsetup); 1593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 160413f72f0SBarry Smith } 161413f72f0SBarry Smith 162d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISGlobalToLocalMappingSetUp_Basic(ISLocalToGlobalMapping mapping) 163d71ae5a4SJacob Faibussowitsch { 164413f72f0SBarry Smith PetscInt i, *idx = mapping->indices, n = mapping->n, end, start, *globals; 165413f72f0SBarry Smith ISLocalToGlobalMapping_Basic *map; 166413f72f0SBarry Smith 167413f72f0SBarry Smith PetscFunctionBegin; 168413f72f0SBarry Smith start = mapping->globalstart; 169413f72f0SBarry Smith end = mapping->globalend; 1709566063dSJacob Faibussowitsch PetscCall(PetscNew(&map)); 1719566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(end - start + 2, &globals)); 172413f72f0SBarry Smith map->globals = globals; 173413f72f0SBarry Smith for (i = 0; i < end - start + 1; i++) globals[i] = -1; 174413f72f0SBarry Smith for (i = 0; i < n; i++) { 175413f72f0SBarry Smith if (idx[i] < 0) continue; 176413f72f0SBarry Smith globals[idx[i] - start] = i; 177413f72f0SBarry Smith } 178413f72f0SBarry Smith mapping->data = (void *)map; 1793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 180413f72f0SBarry Smith } 181413f72f0SBarry Smith 182d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISGlobalToLocalMappingSetUp_Hash(ISLocalToGlobalMapping mapping) 183d71ae5a4SJacob Faibussowitsch { 184413f72f0SBarry Smith PetscInt i, *idx = mapping->indices, n = mapping->n; 185413f72f0SBarry Smith ISLocalToGlobalMapping_Hash *map; 186413f72f0SBarry Smith 187413f72f0SBarry Smith PetscFunctionBegin; 1889566063dSJacob Faibussowitsch PetscCall(PetscNew(&map)); 1899566063dSJacob Faibussowitsch PetscCall(PetscHMapICreate(&map->globalht)); 190413f72f0SBarry Smith for (i = 0; i < n; i++) { 191413f72f0SBarry Smith if (idx[i] < 0) continue; 1929566063dSJacob Faibussowitsch PetscCall(PetscHMapISet(map->globalht, idx[i], i)); 193413f72f0SBarry Smith } 194413f72f0SBarry Smith mapping->data = (void *)map; 1953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 196413f72f0SBarry Smith } 197413f72f0SBarry Smith 198d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISLocalToGlobalMappingDestroy_Basic(ISLocalToGlobalMapping mapping) 199d71ae5a4SJacob Faibussowitsch { 200413f72f0SBarry Smith ISLocalToGlobalMapping_Basic *map = (ISLocalToGlobalMapping_Basic *)mapping->data; 201413f72f0SBarry Smith 202413f72f0SBarry Smith PetscFunctionBegin; 2033ba16761SJacob Faibussowitsch if (!map) PetscFunctionReturn(PETSC_SUCCESS); 2049566063dSJacob Faibussowitsch PetscCall(PetscFree(map->globals)); 2059566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->data)); 2063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 207413f72f0SBarry Smith } 208413f72f0SBarry Smith 209d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISLocalToGlobalMappingDestroy_Hash(ISLocalToGlobalMapping mapping) 210d71ae5a4SJacob Faibussowitsch { 211413f72f0SBarry Smith ISLocalToGlobalMapping_Hash *map = (ISLocalToGlobalMapping_Hash *)mapping->data; 212413f72f0SBarry Smith 213413f72f0SBarry Smith PetscFunctionBegin; 2143ba16761SJacob Faibussowitsch if (!map) PetscFunctionReturn(PETSC_SUCCESS); 2159566063dSJacob Faibussowitsch PetscCall(PetscHMapIDestroy(&map->globalht)); 2169566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->data)); 2173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 218413f72f0SBarry Smith } 219413f72f0SBarry Smith 220413f72f0SBarry Smith #define GTOLTYPE _Basic 221413f72f0SBarry Smith #define GTOLNAME _Basic 222541bf97eSAdrian Croucher #define GTOLBS mapping->bs 2239371c9d4SSatish Balay #define GTOL(g, local) \ 2249371c9d4SSatish Balay do { \ 225413f72f0SBarry Smith local = map->globals[g / bs - start]; \ 2260040bde1SJunchao Zhang if (local >= 0) local = bs * local + (g % bs); \ 227413f72f0SBarry Smith } while (0) 228541bf97eSAdrian Croucher 229413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h> 230413f72f0SBarry Smith 231413f72f0SBarry Smith #define GTOLTYPE _Basic 232413f72f0SBarry Smith #define GTOLNAME Block_Basic 233541bf97eSAdrian Croucher #define GTOLBS 1 2349371c9d4SSatish Balay #define GTOL(g, local) \ 235d71ae5a4SJacob Faibussowitsch do { \ 236d71ae5a4SJacob Faibussowitsch local = map->globals[g - start]; \ 237d71ae5a4SJacob Faibussowitsch } while (0) 238413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h> 239413f72f0SBarry Smith 240413f72f0SBarry Smith #define GTOLTYPE _Hash 241413f72f0SBarry Smith #define GTOLNAME _Hash 242541bf97eSAdrian Croucher #define GTOLBS mapping->bs 2439371c9d4SSatish Balay #define GTOL(g, local) \ 2449371c9d4SSatish Balay do { \ 245e8f14785SLisandro Dalcin (void)PetscHMapIGet(map->globalht, g / bs, &local); \ 2460040bde1SJunchao Zhang if (local >= 0) local = bs * local + (g % bs); \ 247413f72f0SBarry Smith } while (0) 248413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h> 249413f72f0SBarry Smith 250413f72f0SBarry Smith #define GTOLTYPE _Hash 251413f72f0SBarry Smith #define GTOLNAME Block_Hash 252541bf97eSAdrian Croucher #define GTOLBS 1 2539371c9d4SSatish Balay #define GTOL(g, local) \ 254d71ae5a4SJacob Faibussowitsch do { \ 255d71ae5a4SJacob Faibussowitsch (void)PetscHMapIGet(map->globalht, g, &local); \ 256d71ae5a4SJacob Faibussowitsch } while (0) 257413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h> 258413f72f0SBarry Smith 2596658fb44Sstefano_zampini /*@ 2606658fb44Sstefano_zampini ISLocalToGlobalMappingDuplicate - Duplicates the local to global mapping object 2616658fb44Sstefano_zampini 2626658fb44Sstefano_zampini Not Collective 2636658fb44Sstefano_zampini 2646658fb44Sstefano_zampini Input Parameter: 2656658fb44Sstefano_zampini . ltog - local to global mapping 2666658fb44Sstefano_zampini 2676658fb44Sstefano_zampini Output Parameter: 2686658fb44Sstefano_zampini . nltog - the duplicated local to global mapping 2696658fb44Sstefano_zampini 2706658fb44Sstefano_zampini Level: advanced 2716658fb44Sstefano_zampini 272cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()` 2736658fb44Sstefano_zampini @*/ 274d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingDuplicate(ISLocalToGlobalMapping ltog, ISLocalToGlobalMapping *nltog) 275d71ae5a4SJacob Faibussowitsch { 276a0d79125SStefano Zampini ISLocalToGlobalMappingType l2gtype; 2776658fb44Sstefano_zampini 2786658fb44Sstefano_zampini PetscFunctionBegin; 2796658fb44Sstefano_zampini PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1); 2809566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)ltog), ltog->bs, ltog->n, ltog->indices, PETSC_COPY_VALUES, nltog)); 2819566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetType(ltog, &l2gtype)); 2829566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingSetType(*nltog, l2gtype)); 2833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2846658fb44Sstefano_zampini } 2856658fb44Sstefano_zampini 286565245c5SBarry Smith /*@ 287107e9a97SBarry Smith ISLocalToGlobalMappingGetSize - Gets the local size of a local to global mapping 2883b9aefa3SBarry Smith 2893b9aefa3SBarry Smith Not Collective 2903b9aefa3SBarry Smith 2913b9aefa3SBarry Smith Input Parameter: 29238b5cf2dSJacob Faibussowitsch . mapping - local to global mapping 2933b9aefa3SBarry Smith 2943b9aefa3SBarry Smith Output Parameter: 295cab54364SBarry Smith . n - the number of entries in the local mapping, `ISLocalToGlobalMappingGetIndices()` returns an array of this length 2963b9aefa3SBarry Smith 2973b9aefa3SBarry Smith Level: advanced 2983b9aefa3SBarry Smith 299cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()` 3003b9aefa3SBarry Smith @*/ 301d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping mapping, PetscInt *n) 302d71ae5a4SJacob Faibussowitsch { 3033b9aefa3SBarry Smith PetscFunctionBegin; 3040700a824SBarry Smith PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 3054f572ea9SToby Isaac PetscAssertPointer(n, 2); 306107e9a97SBarry Smith *n = mapping->bs * mapping->n; 3073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3083b9aefa3SBarry Smith } 3093b9aefa3SBarry Smith 3105a5d4f66SBarry Smith /*@C 311cab54364SBarry Smith ISLocalToGlobalMappingViewFromOptions - View an `ISLocalToGlobalMapping` based on values in the options database 312fe2efc57SMark 313c3339decSBarry Smith Collective 314fe2efc57SMark 315fe2efc57SMark Input Parameters: 316fe2efc57SMark + A - the local to global mapping object 31720662ed9SBarry Smith . obj - Optional object that provides the options prefix used for the options database query 318736c3998SJose E. Roman - name - command line option 319fe2efc57SMark 320fe2efc57SMark Level: intermediate 321cab54364SBarry Smith 32220662ed9SBarry Smith Note: 32320662ed9SBarry Smith See `PetscObjectViewFromOptions()` for the available `PetscViewer` and `PetscViewerFormat` 32420662ed9SBarry Smith 32520662ed9SBarry Smith .seealso: [](sec_scatter), `PetscViewer`, ``ISLocalToGlobalMapping`, `ISLocalToGlobalMappingView`, `PetscObjectViewFromOptions()`, `ISLocalToGlobalMappingCreate()` 326fe2efc57SMark @*/ 327d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingViewFromOptions(ISLocalToGlobalMapping A, PetscObject obj, const char name[]) 328d71ae5a4SJacob Faibussowitsch { 329fe2efc57SMark PetscFunctionBegin; 330fe2efc57SMark PetscValidHeaderSpecific(A, IS_LTOGM_CLASSID, 1); 3319566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name)); 3323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 333fe2efc57SMark } 334fe2efc57SMark 335fe2efc57SMark /*@C 3365a5d4f66SBarry Smith ISLocalToGlobalMappingView - View a local to global mapping 3375a5d4f66SBarry Smith 338b9cd556bSLois Curfman McInnes Not Collective 339b9cd556bSLois Curfman McInnes 3405a5d4f66SBarry Smith Input Parameters: 34138b5cf2dSJacob Faibussowitsch + mapping - local to global mapping 3423b9aefa3SBarry Smith - viewer - viewer 3435a5d4f66SBarry Smith 344a997ad1aSLois Curfman McInnes Level: advanced 345a997ad1aSLois Curfman McInnes 34620662ed9SBarry Smith .seealso: [](sec_scatter), `PetscViewer`, `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()` 3475a5d4f66SBarry Smith @*/ 348d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingView(ISLocalToGlobalMapping mapping, PetscViewer viewer) 349d71ae5a4SJacob Faibussowitsch { 35032dcc486SBarry Smith PetscInt i; 35132dcc486SBarry Smith PetscMPIInt rank; 352ace3abfcSBarry Smith PetscBool iascii; 3535a5d4f66SBarry Smith 3545a5d4f66SBarry Smith PetscFunctionBegin; 3550700a824SBarry Smith PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 35648a46eb9SPierre Jolivet if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mapping), &viewer)); 3570700a824SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 3585a5d4f66SBarry Smith 3599566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mapping), &rank)); 3609566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 36132077d6dSBarry Smith if (iascii) { 3629566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)mapping, viewer)); 3639566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 364f2c6b1a2SJed Brown for (i = 0; i < mapping->n; i++) { 365f2c6b1a2SJed Brown PetscInt bs = mapping->bs, g = mapping->indices[i]; 366f2c6b1a2SJed Brown if (bs == 1) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %" PetscInt_FMT " %" PetscInt_FMT "\n", rank, i, g)); 367f2c6b1a2SJed Brown else PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %" PetscInt_FMT ":%" PetscInt_FMT " %" PetscInt_FMT ":%" PetscInt_FMT "\n", rank, i * bs, (i + 1) * bs, g * bs, (g + 1) * bs)); 368f2c6b1a2SJed Brown } 3699566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 3709566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 3711575c14dSBarry Smith } 3723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3735a5d4f66SBarry Smith } 3745a5d4f66SBarry Smith 3751f428162SBarry Smith /*@ 3762bdab257SBarry Smith ISLocalToGlobalMappingCreateIS - Creates a mapping between a local (0 to n) 3772bdab257SBarry Smith ordering and a global parallel ordering. 3782bdab257SBarry Smith 37920662ed9SBarry Smith Not Collective 380b9cd556bSLois Curfman McInnes 381a997ad1aSLois Curfman McInnes Input Parameter: 3828c03b21aSDmitry Karpeev . is - index set containing the global numbers for each local number 3832bdab257SBarry Smith 384a997ad1aSLois Curfman McInnes Output Parameter: 3852bdab257SBarry Smith . mapping - new mapping data structure 3862bdab257SBarry Smith 387a997ad1aSLois Curfman McInnes Level: advanced 388a997ad1aSLois Curfman McInnes 389cab54364SBarry Smith Note: 390cab54364SBarry Smith the block size of the `IS` determines the block size of the mapping 391cab54364SBarry Smith 392cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetFromOptions()` 3932bdab257SBarry Smith @*/ 394d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingCreateIS(IS is, ISLocalToGlobalMapping *mapping) 395d71ae5a4SJacob Faibussowitsch { 3963bbf0e92SBarry Smith PetscInt n, bs; 3975d0c19d7SBarry Smith const PetscInt *indices; 3982bdab257SBarry Smith MPI_Comm comm; 3993bbf0e92SBarry Smith PetscBool isblock; 4003a40ed3dSBarry Smith 4013a40ed3dSBarry Smith PetscFunctionBegin; 4020700a824SBarry Smith PetscValidHeaderSpecific(is, IS_CLASSID, 1); 4034f572ea9SToby Isaac PetscAssertPointer(mapping, 2); 4042bdab257SBarry Smith 4059566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)is, &comm)); 4069566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &n)); 4079566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)is, ISBLOCK, &isblock)); 4086006e8d2SBarry Smith if (!isblock) { 4099566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &indices)); 4109566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm, 1, n, indices, PETSC_COPY_VALUES, mapping)); 4119566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &indices)); 4126006e8d2SBarry Smith } else { 4139566063dSJacob Faibussowitsch PetscCall(ISGetBlockSize(is, &bs)); 4149566063dSJacob Faibussowitsch PetscCall(ISBlockGetIndices(is, &indices)); 4159566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm, bs, n / bs, indices, PETSC_COPY_VALUES, mapping)); 4169566063dSJacob Faibussowitsch PetscCall(ISBlockRestoreIndices(is, &indices)); 4176006e8d2SBarry Smith } 4183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4192bdab257SBarry Smith } 4205a5d4f66SBarry Smith 421a4d96a55SJed Brown /*@C 422a4d96a55SJed Brown ISLocalToGlobalMappingCreateSF - Creates a mapping between a local (0 to n) 423a4d96a55SJed Brown ordering and a global parallel ordering. 424a4d96a55SJed Brown 425a4d96a55SJed Brown Collective 426a4d96a55SJed Brown 427d8d19677SJose E. Roman Input Parameters: 428a4d96a55SJed Brown + sf - star forest mapping contiguous local indices to (rank, offset) 429cab54364SBarry Smith - start - first global index on this process, or `PETSC_DECIDE` to compute contiguous global numbering automatically 430a4d96a55SJed Brown 431a4d96a55SJed Brown Output Parameter: 432a4d96a55SJed Brown . mapping - new mapping data structure 433a4d96a55SJed Brown 434a4d96a55SJed Brown Level: advanced 435a4d96a55SJed Brown 43620662ed9SBarry Smith Note: 43720662ed9SBarry Smith If any processor calls this with `start` = `PETSC_DECIDE` then all processors must, otherwise the program will hang. 4389a535bafSVaclav Hapla 439cab54364SBarry Smith .seealso: [](sec_scatter), `PetscSF`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingSetFromOptions()` 440a4d96a55SJed Brown @*/ 441d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF sf, PetscInt start, ISLocalToGlobalMapping *mapping) 442d71ae5a4SJacob Faibussowitsch { 443a4d96a55SJed Brown PetscInt i, maxlocal, nroots, nleaves, *globals, *ltog; 444a4d96a55SJed Brown MPI_Comm comm; 445a4d96a55SJed Brown 446a4d96a55SJed Brown PetscFunctionBegin; 447a4d96a55SJed Brown PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 4484f572ea9SToby Isaac PetscAssertPointer(mapping, 3); 4499566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf, &comm)); 45041f4c31fSVaclav Hapla PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, NULL, NULL)); 4519a535bafSVaclav Hapla if (start == PETSC_DECIDE) { 4529a535bafSVaclav Hapla start = 0; 4539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Exscan(&nroots, &start, 1, MPIU_INT, MPI_SUM, comm)); 45441f4c31fSVaclav Hapla } else PetscCheck(start >= 0, comm, PETSC_ERR_ARG_OUTOFRANGE, "start must be nonnegative or PETSC_DECIDE"); 45541f4c31fSVaclav Hapla PetscCall(PetscSFGetLeafRange(sf, NULL, &maxlocal)); 45641f4c31fSVaclav Hapla ++maxlocal; 4579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nroots, &globals)); 4589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxlocal, <og)); 459a4d96a55SJed Brown for (i = 0; i < nroots; i++) globals[i] = start + i; 460a4d96a55SJed Brown for (i = 0; i < maxlocal; i++) ltog[i] = -1; 4619566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf, MPIU_INT, globals, ltog, MPI_REPLACE)); 4629566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf, MPIU_INT, globals, ltog, MPI_REPLACE)); 4639566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm, 1, maxlocal, ltog, PETSC_OWN_POINTER, mapping)); 4649566063dSJacob Faibussowitsch PetscCall(PetscFree(globals)); 4653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 466a4d96a55SJed Brown } 467b46b645bSBarry Smith 46863fa5c83Sstefano_zampini /*@ 46963fa5c83Sstefano_zampini ISLocalToGlobalMappingSetBlockSize - Sets the blocksize of the mapping 47063fa5c83Sstefano_zampini 47120662ed9SBarry Smith Not Collective 47263fa5c83Sstefano_zampini 47363fa5c83Sstefano_zampini Input Parameters: 474a2b725a8SWilliam Gropp + mapping - mapping data structure 475a2b725a8SWilliam Gropp - bs - the blocksize 47663fa5c83Sstefano_zampini 47763fa5c83Sstefano_zampini Level: advanced 47863fa5c83Sstefano_zampini 479cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()` 48063fa5c83Sstefano_zampini @*/ 481d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingSetBlockSize(ISLocalToGlobalMapping mapping, PetscInt bs) 482d71ae5a4SJacob Faibussowitsch { 483a59f3c4dSstefano_zampini PetscInt *nid; 484a59f3c4dSstefano_zampini const PetscInt *oid; 485a59f3c4dSstefano_zampini PetscInt i, cn, on, obs, nn; 48663fa5c83Sstefano_zampini 48763fa5c83Sstefano_zampini PetscFunctionBegin; 48863fa5c83Sstefano_zampini PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 48908401ef6SPierre Jolivet PetscCheck(bs >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid block size %" PetscInt_FMT, bs); 4903ba16761SJacob Faibussowitsch if (bs == mapping->bs) PetscFunctionReturn(PETSC_SUCCESS); 49163fa5c83Sstefano_zampini on = mapping->n; 49263fa5c83Sstefano_zampini obs = mapping->bs; 49363fa5c83Sstefano_zampini oid = mapping->indices; 49463fa5c83Sstefano_zampini nn = (on * obs) / bs; 49508401ef6SPierre 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); 496a59f3c4dSstefano_zampini 4979566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nn, &nid)); 4989566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(mapping, &oid)); 499a59f3c4dSstefano_zampini for (i = 0; i < nn; i++) { 500a59f3c4dSstefano_zampini PetscInt j; 501a59f3c4dSstefano_zampini for (j = 0, cn = 0; j < bs - 1; j++) { 5029371c9d4SSatish Balay if (oid[i * bs + j] < 0) { 5039371c9d4SSatish Balay cn++; 5049371c9d4SSatish Balay continue; 5059371c9d4SSatish Balay } 50608401ef6SPierre 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]); 507a59f3c4dSstefano_zampini } 508a59f3c4dSstefano_zampini if (oid[i * bs + j] < 0) cn++; 5098b7cb0e6Sstefano_zampini if (cn) { 51008401ef6SPierre 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); 511a59f3c4dSstefano_zampini nid[i] = -1; 5128b7cb0e6Sstefano_zampini } else { 513a59f3c4dSstefano_zampini nid[i] = oid[i * bs] / bs; 51463fa5c83Sstefano_zampini } 51563fa5c83Sstefano_zampini } 5169566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(mapping, &oid)); 517a59f3c4dSstefano_zampini 51863fa5c83Sstefano_zampini mapping->n = nn; 51963fa5c83Sstefano_zampini mapping->bs = bs; 5209566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->indices)); 52163fa5c83Sstefano_zampini mapping->indices = nid; 522c9345713Sstefano_zampini mapping->globalstart = 0; 523c9345713Sstefano_zampini mapping->globalend = 0; 5241bd0b88eSStefano Zampini 5251bd0b88eSStefano Zampini /* reset the cached information */ 5269566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->info_procs)); 5279566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->info_numprocs)); 5281bd0b88eSStefano Zampini if (mapping->info_indices) { 5291bd0b88eSStefano Zampini PetscInt i; 5301bd0b88eSStefano Zampini 531*f4f49eeaSPierre Jolivet PetscCall(PetscFree(mapping->info_indices[0])); 53248a46eb9SPierre Jolivet for (i = 1; i < mapping->info_nproc; i++) PetscCall(PetscFree(mapping->info_indices[i])); 5339566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->info_indices)); 5341bd0b88eSStefano Zampini } 5351bd0b88eSStefano Zampini mapping->info_cached = PETSC_FALSE; 5361bd0b88eSStefano Zampini 537dbbe0bcdSBarry Smith PetscTryTypeMethod(mapping, destroy); 5383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 53963fa5c83Sstefano_zampini } 54063fa5c83Sstefano_zampini 54145b6f7e9SBarry Smith /*@ 54245b6f7e9SBarry Smith ISLocalToGlobalMappingGetBlockSize - Gets the blocksize of the mapping 54345b6f7e9SBarry Smith ordering and a global parallel ordering. 54445b6f7e9SBarry Smith 54545b6f7e9SBarry Smith Not Collective 54645b6f7e9SBarry Smith 5472fe279fdSBarry Smith Input Parameter: 54845b6f7e9SBarry Smith . mapping - mapping data structure 54945b6f7e9SBarry Smith 55045b6f7e9SBarry Smith Output Parameter: 55145b6f7e9SBarry Smith . bs - the blocksize 55245b6f7e9SBarry Smith 55345b6f7e9SBarry Smith Level: advanced 55445b6f7e9SBarry Smith 555cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()` 55645b6f7e9SBarry Smith @*/ 557d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping mapping, PetscInt *bs) 558d71ae5a4SJacob Faibussowitsch { 55945b6f7e9SBarry Smith PetscFunctionBegin; 560cbc1caf0SMatthew G. Knepley PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 56145b6f7e9SBarry Smith *bs = mapping->bs; 5623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 56345b6f7e9SBarry Smith } 56445b6f7e9SBarry Smith 565ba5bb76aSSatish Balay /*@ 56690f02eecSBarry Smith ISLocalToGlobalMappingCreate - Creates a mapping between a local (0 to n) 56790f02eecSBarry Smith ordering and a global parallel ordering. 5682362add9SBarry Smith 56989d82c54SBarry Smith Not Collective, but communicator may have more than one process 570b9cd556bSLois Curfman McInnes 5712362add9SBarry Smith Input Parameters: 57289d82c54SBarry Smith + comm - MPI communicator 573f0413b6fSBarry Smith . bs - the block size 57428bc9809SBarry Smith . n - the number of local elements divided by the block size, or equivalently the number of block indices 57528bc9809SBarry 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 576d5ad8652SBarry Smith - mode - see PetscCopyMode 5772362add9SBarry Smith 578a997ad1aSLois Curfman McInnes Output Parameter: 57990f02eecSBarry Smith . mapping - new mapping data structure 5802362add9SBarry Smith 581cab54364SBarry Smith Level: advanced 582cab54364SBarry Smith 58395452b02SPatrick Sanan Notes: 58495452b02SPatrick Sanan There is one integer value in indices per block and it represents the actual indices bs*idx + j, where j=0,..,bs-1 585413f72f0SBarry Smith 586cab54364SBarry Smith For "small" problems when using `ISGlobalToLocalMappingApply()` and `ISGlobalToLocalMappingApplyBlock()`, the `ISLocalToGlobalMappingType` 587cab54364SBarry Smith of `ISLOCALTOGLOBALMAPPINGBASIC` will be used; this uses more memory but is faster; this approach is not scalable for extremely large mappings. 588413f72f0SBarry Smith 589cab54364SBarry Smith For large problems `ISLOCALTOGLOBALMAPPINGHASH` is used, this is scalable. 59020662ed9SBarry Smith Use `ISLocalToGlobalMappingSetType()` or call `ISLocalToGlobalMappingSetFromOptions()` with the option 59120662ed9SBarry Smith `-islocaltoglobalmapping_type` <`basic`,`hash`> to control which is used. 592a997ad1aSLois Curfman McInnes 59320662ed9SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingSetFromOptions()`, 59420662ed9SBarry Smith `ISLOCALTOGLOBALMAPPINGBASIC`, `ISLOCALTOGLOBALMAPPINGHASH` 595db781477SPatrick Sanan `ISLocalToGlobalMappingSetType()`, `ISLocalToGlobalMappingType` 5962362add9SBarry Smith @*/ 597d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingCreate(MPI_Comm comm, PetscInt bs, PetscInt n, const PetscInt indices[], PetscCopyMode mode, ISLocalToGlobalMapping *mapping) 598d71ae5a4SJacob Faibussowitsch { 59932dcc486SBarry Smith PetscInt *in; 600b46b645bSBarry Smith 601b46b645bSBarry Smith PetscFunctionBegin; 6024f572ea9SToby Isaac if (n) PetscAssertPointer(indices, 4); 6034f572ea9SToby Isaac PetscAssertPointer(mapping, 6); 604b46b645bSBarry Smith 6050298fd71SBarry Smith *mapping = NULL; 6069566063dSJacob Faibussowitsch PetscCall(ISInitializePackage()); 6072362add9SBarry Smith 6089566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(*mapping, IS_LTOGM_CLASSID, "ISLocalToGlobalMapping", "Local to global mapping", "IS", comm, ISLocalToGlobalMappingDestroy, ISLocalToGlobalMappingView)); 609d4bb536fSBarry Smith (*mapping)->n = n; 610f0413b6fSBarry Smith (*mapping)->bs = bs; 611d5ad8652SBarry Smith if (mode == PETSC_COPY_VALUES) { 6129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &in)); 6139566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(in, indices, n)); 614d5ad8652SBarry Smith (*mapping)->indices = in; 61571910c26SVaclav Hapla (*mapping)->dealloc_indices = PETSC_TRUE; 6166389a1a1SBarry Smith } else if (mode == PETSC_OWN_POINTER) { 6176389a1a1SBarry Smith (*mapping)->indices = (PetscInt *)indices; 61871910c26SVaclav Hapla (*mapping)->dealloc_indices = PETSC_TRUE; 61971910c26SVaclav Hapla } else if (mode == PETSC_USE_POINTER) { 62071910c26SVaclav Hapla (*mapping)->indices = (PetscInt *)indices; 6219371c9d4SSatish Balay } else SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid mode %d", mode); 6223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6232362add9SBarry Smith } 6242362add9SBarry Smith 625413f72f0SBarry Smith PetscFunctionList ISLocalToGlobalMappingList = NULL; 626413f72f0SBarry Smith 62790f02eecSBarry Smith /*@ 6287e99dc12SLawrence Mitchell ISLocalToGlobalMappingSetFromOptions - Set mapping options from the options database. 6297e99dc12SLawrence Mitchell 63020662ed9SBarry Smith Not Collective 6317e99dc12SLawrence Mitchell 6322fe279fdSBarry Smith Input Parameter: 6337e99dc12SLawrence Mitchell . mapping - mapping data structure 6347e99dc12SLawrence Mitchell 63520662ed9SBarry Smith Options Database Key: 63620662ed9SBarry Smith . -islocaltoglobalmapping_type - <basic,hash> nonscalable and scalable versions 63720662ed9SBarry Smith 6387e99dc12SLawrence Mitchell Level: advanced 6397e99dc12SLawrence Mitchell 64042747ad1SJacob Faibussowitsch .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingDestroy()`, 64142747ad1SJacob Faibussowitsch `ISLocalToGlobalMappingCreateIS()`, `ISLOCALTOGLOBALMAPPINGBASIC`, 64242747ad1SJacob Faibussowitsch `ISLOCALTOGLOBALMAPPINGHASH`, `ISLocalToGlobalMappingSetType()`, `ISLocalToGlobalMappingType` 6437e99dc12SLawrence Mitchell @*/ 644d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingSetFromOptions(ISLocalToGlobalMapping mapping) 645d71ae5a4SJacob Faibussowitsch { 646413f72f0SBarry Smith char type[256]; 647413f72f0SBarry Smith ISLocalToGlobalMappingType defaulttype = "Not set"; 6487e99dc12SLawrence Mitchell PetscBool flg; 6497e99dc12SLawrence Mitchell 6507e99dc12SLawrence Mitchell PetscFunctionBegin; 6517e99dc12SLawrence Mitchell PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 6529566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRegisterAll()); 653d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)mapping); 6549566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-islocaltoglobalmapping_type", "ISLocalToGlobalMapping method", "ISLocalToGlobalMappingSetType", ISLocalToGlobalMappingList, (char *)(((PetscObject)mapping)->type_name) ? ((PetscObject)mapping)->type_name : defaulttype, type, 256, &flg)); 6551baa6e33SBarry Smith if (flg) PetscCall(ISLocalToGlobalMappingSetType(mapping, type)); 656d0609cedSBarry Smith PetscOptionsEnd(); 6573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6587e99dc12SLawrence Mitchell } 6597e99dc12SLawrence Mitchell 6607e99dc12SLawrence Mitchell /*@ 66190f02eecSBarry Smith ISLocalToGlobalMappingDestroy - Destroys a mapping between a local (0 to n) 66290f02eecSBarry Smith ordering and a global parallel ordering. 66390f02eecSBarry Smith 66420662ed9SBarry Smith Not Collective 665b9cd556bSLois Curfman McInnes 6662fe279fdSBarry Smith Input Parameter: 66790f02eecSBarry Smith . mapping - mapping data structure 66890f02eecSBarry Smith 669a997ad1aSLois Curfman McInnes Level: advanced 670a997ad1aSLois Curfman McInnes 671cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingCreate()` 67290f02eecSBarry Smith @*/ 673d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping *mapping) 674d71ae5a4SJacob Faibussowitsch { 6753a40ed3dSBarry Smith PetscFunctionBegin; 6763ba16761SJacob Faibussowitsch if (!*mapping) PetscFunctionReturn(PETSC_SUCCESS); 677*f4f49eeaSPierre Jolivet PetscValidHeaderSpecific(*mapping, IS_LTOGM_CLASSID, 1); 678*f4f49eeaSPierre Jolivet if (--((PetscObject)*mapping)->refct > 0) { 6799371c9d4SSatish Balay *mapping = NULL; 6803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 68171910c26SVaclav Hapla } 68248a46eb9SPierre Jolivet if ((*mapping)->dealloc_indices) PetscCall(PetscFree((*mapping)->indices)); 6839566063dSJacob Faibussowitsch PetscCall(PetscFree((*mapping)->info_procs)); 6849566063dSJacob Faibussowitsch PetscCall(PetscFree((*mapping)->info_numprocs)); 685268a049cSStefano Zampini if ((*mapping)->info_indices) { 686268a049cSStefano Zampini PetscInt i; 687268a049cSStefano Zampini 6889566063dSJacob Faibussowitsch PetscCall(PetscFree(((*mapping)->info_indices)[0])); 68948a46eb9SPierre Jolivet for (i = 1; i < (*mapping)->info_nproc; i++) PetscCall(PetscFree(((*mapping)->info_indices)[i])); 6909566063dSJacob Faibussowitsch PetscCall(PetscFree((*mapping)->info_indices)); 691268a049cSStefano Zampini } 69248a46eb9SPierre Jolivet if ((*mapping)->info_nodei) PetscCall(PetscFree(((*mapping)->info_nodei)[0])); 6939566063dSJacob Faibussowitsch PetscCall(PetscFree2((*mapping)->info_nodec, (*mapping)->info_nodei)); 694213acdd3SPierre Jolivet PetscTryTypeMethod(*mapping, destroy); 6959566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(mapping)); 6964c8fdceaSLisandro Dalcin *mapping = NULL; 6973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 69890f02eecSBarry Smith } 69990f02eecSBarry Smith 70090f02eecSBarry Smith /*@ 701cab54364SBarry Smith ISLocalToGlobalMappingApplyIS - Creates from an `IS` in the local numbering 702cab54364SBarry Smith a new index set using the global numbering defined in an `ISLocalToGlobalMapping` 7033acfe500SLois Curfman McInnes context. 70490f02eecSBarry Smith 705c3339decSBarry Smith Collective 706b9cd556bSLois Curfman McInnes 70790f02eecSBarry Smith Input Parameters: 708b9cd556bSLois Curfman McInnes + mapping - mapping between local and global numbering 709b9cd556bSLois Curfman McInnes - is - index set in local numbering 71090f02eecSBarry Smith 711cab54364SBarry Smith Output Parameter: 71290f02eecSBarry Smith . newis - index set in global numbering 71390f02eecSBarry Smith 714a997ad1aSLois Curfman McInnes Level: advanced 715a997ad1aSLois Curfman McInnes 716cab54364SBarry Smith Note: 71720662ed9SBarry Smith The output `IS` will have the same communicator as the input `IS`. 718cab54364SBarry Smith 719cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingCreate()`, 720db781477SPatrick Sanan `ISLocalToGlobalMappingDestroy()`, `ISGlobalToLocalMappingApply()` 72190f02eecSBarry Smith @*/ 722d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping, IS is, IS *newis) 723d71ae5a4SJacob Faibussowitsch { 724e24637baSBarry Smith PetscInt n, *idxout; 7255d0c19d7SBarry Smith const PetscInt *idxin; 7263a40ed3dSBarry Smith 7273a40ed3dSBarry Smith PetscFunctionBegin; 7280700a824SBarry Smith PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 7290700a824SBarry Smith PetscValidHeaderSpecific(is, IS_CLASSID, 2); 7304f572ea9SToby Isaac PetscAssertPointer(newis, 3); 73190f02eecSBarry Smith 7329566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &n)); 7339566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &idxin)); 7349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &idxout)); 7359566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(mapping, n, idxin, idxout)); 7369566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &idxin)); 7379566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)is), n, idxout, PETSC_OWN_POINTER, newis)); 7383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 73990f02eecSBarry Smith } 74090f02eecSBarry Smith 741b89cb25eSSatish Balay /*@ 7423acfe500SLois Curfman McInnes ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering 7433acfe500SLois Curfman McInnes and converts them to the global numbering. 74490f02eecSBarry Smith 74520662ed9SBarry Smith Not Collective 746b9cd556bSLois Curfman McInnes 747bb25748dSBarry Smith Input Parameters: 748b9cd556bSLois Curfman McInnes + mapping - the local to global mapping context 749bb25748dSBarry Smith . N - number of integers 750b9cd556bSLois Curfman McInnes - in - input indices in local numbering 751bb25748dSBarry Smith 752bb25748dSBarry Smith Output Parameter: 753bb25748dSBarry Smith . out - indices in global numbering 754bb25748dSBarry Smith 755a997ad1aSLois Curfman McInnes Level: advanced 756a997ad1aSLois Curfman McInnes 757cab54364SBarry Smith Note: 75820662ed9SBarry Smith The `in` and `out` array parameters may be identical. 759cab54364SBarry Smith 760cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApplyBlock()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingDestroy()`, 761c2e3fba1SPatrick Sanan `ISLocalToGlobalMappingApplyIS()`, `AOCreateBasic()`, `AOApplicationToPetsc()`, 762db781477SPatrick Sanan `AOPetscToApplication()`, `ISGlobalToLocalMappingApply()` 763afcb2eb5SJed Brown @*/ 764d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping, PetscInt N, const PetscInt in[], PetscInt out[]) 765d71ae5a4SJacob Faibussowitsch { 766cbc1caf0SMatthew G. Knepley PetscInt i, bs, Nmax; 76745b6f7e9SBarry Smith 76845b6f7e9SBarry Smith PetscFunctionBegin; 769cbc1caf0SMatthew G. Knepley PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 770cbc1caf0SMatthew G. Knepley bs = mapping->bs; 771cbc1caf0SMatthew G. Knepley Nmax = bs * mapping->n; 77245b6f7e9SBarry Smith if (bs == 1) { 773cbc1caf0SMatthew G. Knepley const PetscInt *idx = mapping->indices; 77445b6f7e9SBarry Smith for (i = 0; i < N; i++) { 77545b6f7e9SBarry Smith if (in[i] < 0) { 77645b6f7e9SBarry Smith out[i] = in[i]; 77745b6f7e9SBarry Smith continue; 77845b6f7e9SBarry Smith } 77908401ef6SPierre 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); 78045b6f7e9SBarry Smith out[i] = idx[in[i]]; 78145b6f7e9SBarry Smith } 78245b6f7e9SBarry Smith } else { 783cbc1caf0SMatthew G. Knepley const PetscInt *idx = mapping->indices; 78445b6f7e9SBarry Smith for (i = 0; i < N; i++) { 78545b6f7e9SBarry Smith if (in[i] < 0) { 78645b6f7e9SBarry Smith out[i] = in[i]; 78745b6f7e9SBarry Smith continue; 78845b6f7e9SBarry Smith } 78908401ef6SPierre 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); 79045b6f7e9SBarry Smith out[i] = idx[in[i] / bs] * bs + (in[i] % bs); 79145b6f7e9SBarry Smith } 79245b6f7e9SBarry Smith } 7933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 79445b6f7e9SBarry Smith } 79545b6f7e9SBarry Smith 79645b6f7e9SBarry Smith /*@ 7976006e8d2SBarry Smith ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local block numbering and converts them to the global block numbering 79845b6f7e9SBarry Smith 79920662ed9SBarry Smith Not Collective 80045b6f7e9SBarry Smith 80145b6f7e9SBarry Smith Input Parameters: 80245b6f7e9SBarry Smith + mapping - the local to global mapping context 80345b6f7e9SBarry Smith . N - number of integers 8046006e8d2SBarry Smith - in - input indices in local block numbering 80545b6f7e9SBarry Smith 80645b6f7e9SBarry Smith Output Parameter: 8076006e8d2SBarry Smith . out - indices in global block numbering 80845b6f7e9SBarry Smith 8096006e8d2SBarry Smith Example: 810cab54364SBarry Smith If the index values are {0,1,6,7} set with a call to `ISLocalToGlobalMappingCreate`(`PETSC_COMM_SELF`,2,2,{0,3}) then the mapping applied to 0 8116006e8d2SBarry Smith (the first block) would produce 0 and the mapping applied to 1 (the second block) would produce 3. 8126006e8d2SBarry Smith 81320662ed9SBarry Smith Level: advanced 81420662ed9SBarry Smith 81520662ed9SBarry Smith Note: 81620662ed9SBarry Smith The `in` and `out` array parameters may be identical. 81720662ed9SBarry Smith 818cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingDestroy()`, 819c2e3fba1SPatrick Sanan `ISLocalToGlobalMappingApplyIS()`, `AOCreateBasic()`, `AOApplicationToPetsc()`, 820db781477SPatrick Sanan `AOPetscToApplication()`, `ISGlobalToLocalMappingApply()` 82145b6f7e9SBarry Smith @*/ 822d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping, PetscInt N, const PetscInt in[], PetscInt out[]) 823d71ae5a4SJacob Faibussowitsch { 8248a1f772fSStefano Zampini PetscInt i, Nmax; 8258a1f772fSStefano Zampini const PetscInt *idx; 826d4bb536fSBarry Smith 827a0d79125SStefano Zampini PetscFunctionBegin; 828a0d79125SStefano Zampini PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 8298a1f772fSStefano Zampini Nmax = mapping->n; 8308a1f772fSStefano Zampini idx = mapping->indices; 831afcb2eb5SJed Brown for (i = 0; i < N; i++) { 832afcb2eb5SJed Brown if (in[i] < 0) { 833afcb2eb5SJed Brown out[i] = in[i]; 834afcb2eb5SJed Brown continue; 835afcb2eb5SJed Brown } 83608401ef6SPierre 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); 837afcb2eb5SJed Brown out[i] = idx[in[i]]; 838afcb2eb5SJed Brown } 8393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 840afcb2eb5SJed Brown } 841d4bb536fSBarry Smith 8427e99dc12SLawrence Mitchell /*@ 843a997ad1aSLois Curfman McInnes ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers 844a997ad1aSLois Curfman McInnes specified with a global numbering. 845d4bb536fSBarry Smith 84620662ed9SBarry Smith Not Collective 847b9cd556bSLois Curfman McInnes 848d4bb536fSBarry Smith Input Parameters: 849b9cd556bSLois Curfman McInnes + mapping - mapping between local and global numbering 850cab54364SBarry Smith . type - `IS_GTOLM_MASK` - maps global indices with no local value to -1 in the output list (i.e., mask them) 851cab54364SBarry Smith `IS_GTOLM_DROP` - drops the indices with no local value from the output list 852d4bb536fSBarry Smith . n - number of global indices to map 853b9cd556bSLois Curfman McInnes - idx - global indices to map 854d4bb536fSBarry Smith 855d4bb536fSBarry Smith Output Parameters: 856cab54364SBarry Smith + nout - number of indices in output array (if type == `IS_GTOLM_MASK` then nout = n) 857b9cd556bSLois Curfman McInnes - idxout - local index of each global index, one must pass in an array long enough 858cab54364SBarry Smith to hold all the indices. You can call `ISGlobalToLocalMappingApply()` with 8590298fd71SBarry Smith idxout == NULL to determine the required length (returned in nout) 860cab54364SBarry Smith and then allocate the required space and call `ISGlobalToLocalMappingApply()` 861e182c471SBarry Smith a second time to set the values. 862d4bb536fSBarry Smith 863cab54364SBarry Smith Level: advanced 864cab54364SBarry Smith 865b9cd556bSLois Curfman McInnes Notes: 86620662ed9SBarry Smith Either `nout` or `idxout` may be `NULL`. `idx` and `idxout` may be identical. 867d4bb536fSBarry Smith 86820662ed9SBarry Smith For "small" problems when using `ISGlobalToLocalMappingApply()` and `ISGlobalToLocalMappingApplyBlock()`, the `ISLocalToGlobalMappingType` of 86920662ed9SBarry Smith `ISLOCALTOGLOBALMAPPINGBASIC` will be used; 870cab54364SBarry 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. 871cab54364SBarry Smith Use `ISLocalToGlobalMappingSetType()` or call `ISLocalToGlobalMappingSetFromOptions()` with the option -islocaltoglobalmapping_type <basic,hash> to control which is used. 8720f5bd95cSBarry Smith 87338b5cf2dSJacob Faibussowitsch Developer Notes: 87420662ed9SBarry Smith The manual page states that `idx` and `idxout` may be identical but the calling 87520662ed9SBarry Smith sequence declares `idx` as const so it cannot be the same as `idxout`. 87632fd6b96SBarry Smith 877cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApply()`, `ISGlobalToLocalMappingApplyBlock()`, `ISLocalToGlobalMappingCreate()`, 878db781477SPatrick Sanan `ISLocalToGlobalMappingDestroy()` 879d4bb536fSBarry Smith @*/ 880d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, PetscInt n, const PetscInt idx[], PetscInt *nout, PetscInt idxout[]) 881d71ae5a4SJacob Faibussowitsch { 8829d90f715SBarry Smith PetscFunctionBegin; 8839d90f715SBarry Smith PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 88448a46eb9SPierre Jolivet if (!mapping->data) PetscCall(ISGlobalToLocalMappingSetUp(mapping)); 885dbbe0bcdSBarry Smith PetscUseTypeMethod(mapping, globaltolocalmappingapply, type, n, idx, nout, idxout); 8863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8879d90f715SBarry Smith } 8889d90f715SBarry Smith 889d4fe737eSStefano Zampini /*@ 890cab54364SBarry Smith ISGlobalToLocalMappingApplyIS - Creates from an `IS` in the global numbering 891cab54364SBarry Smith a new index set using the local numbering defined in an `ISLocalToGlobalMapping` 892d4fe737eSStefano Zampini context. 893d4fe737eSStefano Zampini 89420662ed9SBarry Smith Not Collective 895d4fe737eSStefano Zampini 896d4fe737eSStefano Zampini Input Parameters: 897d4fe737eSStefano Zampini + mapping - mapping between local and global numbering 898cab54364SBarry Smith . type - `IS_GTOLM_MASK` - maps global indices with no local value to -1 in the output list (i.e., mask them) 899cab54364SBarry Smith `IS_GTOLM_DROP` - drops the indices with no local value from the output list 900d4fe737eSStefano Zampini - is - index set in global numbering 901d4fe737eSStefano Zampini 9022fe279fdSBarry Smith Output Parameter: 903d4fe737eSStefano Zampini . newis - index set in local numbering 904d4fe737eSStefano Zampini 905d4fe737eSStefano Zampini Level: advanced 906d4fe737eSStefano Zampini 907cab54364SBarry Smith Note: 908cab54364SBarry Smith The output `IS` will be sequential, as it encodes a purely local operation 909cab54364SBarry Smith 910cab54364SBarry Smith .seealso: [](sec_scatter), `ISGlobalToLocalMapping`, `ISGlobalToLocalMappingApply()`, `ISLocalToGlobalMappingCreate()`, 911db781477SPatrick Sanan `ISLocalToGlobalMappingDestroy()` 912d4fe737eSStefano Zampini @*/ 913d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, IS is, IS *newis) 914d71ae5a4SJacob Faibussowitsch { 915d4fe737eSStefano Zampini PetscInt n, nout, *idxout; 916d4fe737eSStefano Zampini const PetscInt *idxin; 917d4fe737eSStefano Zampini 918d4fe737eSStefano Zampini PetscFunctionBegin; 919d4fe737eSStefano Zampini PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 920d4fe737eSStefano Zampini PetscValidHeaderSpecific(is, IS_CLASSID, 3); 9214f572ea9SToby Isaac PetscAssertPointer(newis, 4); 922d4fe737eSStefano Zampini 9239566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &n)); 9249566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &idxin)); 925d4fe737eSStefano Zampini if (type == IS_GTOLM_MASK) { 9269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &idxout)); 927d4fe737eSStefano Zampini } else { 9289566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(mapping, type, n, idxin, &nout, NULL)); 9299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nout, &idxout)); 930d4fe737eSStefano Zampini } 9319566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(mapping, type, n, idxin, &nout, idxout)); 9329566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &idxin)); 9339566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nout, idxout, PETSC_OWN_POINTER, newis)); 9343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 935d4fe737eSStefano Zampini } 936d4fe737eSStefano Zampini 9379d90f715SBarry Smith /*@ 9389d90f715SBarry Smith ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers 9399d90f715SBarry Smith specified with a block global numbering. 9409d90f715SBarry Smith 94120662ed9SBarry Smith Not Collective 9429d90f715SBarry Smith 9439d90f715SBarry Smith Input Parameters: 9449d90f715SBarry Smith + mapping - mapping between local and global numbering 945cab54364SBarry Smith . type - `IS_GTOLM_MASK` - maps global indices with no local value to -1 in the output list (i.e., mask them) 946cab54364SBarry Smith `IS_GTOLM_DROP` - drops the indices with no local value from the output list 9479d90f715SBarry Smith . n - number of global indices to map 9489d90f715SBarry Smith - idx - global indices to map 9499d90f715SBarry Smith 9509d90f715SBarry Smith Output Parameters: 951cab54364SBarry Smith + nout - number of indices in output array (if type == `IS_GTOLM_MASK` then nout = n) 9529d90f715SBarry Smith - idxout - local index of each global index, one must pass in an array long enough 953cab54364SBarry Smith to hold all the indices. You can call `ISGlobalToLocalMappingApplyBlock()` with 9549d90f715SBarry Smith idxout == NULL to determine the required length (returned in nout) 955cab54364SBarry Smith and then allocate the required space and call `ISGlobalToLocalMappingApplyBlock()` 9569d90f715SBarry Smith a second time to set the values. 9579d90f715SBarry Smith 958cab54364SBarry Smith Level: advanced 959cab54364SBarry Smith 9609d90f715SBarry Smith Notes: 96120662ed9SBarry Smith Either `nout` or `idxout` may be `NULL`. `idx` and `idxout` may be identical. 9629d90f715SBarry Smith 96320662ed9SBarry Smith For "small" problems when using `ISGlobalToLocalMappingApply()` and `ISGlobalToLocalMappingApplyBlock()`, the `ISLocalToGlobalMappingType` of 96420662ed9SBarry Smith `ISLOCALTOGLOBALMAPPINGBASIC` will be used; 965cab54364SBarry 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. 966cab54364SBarry Smith Use `ISLocalToGlobalMappingSetType()` or call `ISLocalToGlobalMappingSetFromOptions()` with the option -islocaltoglobalmapping_type <basic,hash> to control which is used. 9679d90f715SBarry Smith 96838b5cf2dSJacob Faibussowitsch Developer Notes: 96920662ed9SBarry Smith The manual page states that `idx` and `idxout` may be identical but the calling 97020662ed9SBarry Smith sequence declares `idx` as const so it cannot be the same as `idxout`. 9719d90f715SBarry Smith 972cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApply()`, `ISGlobalToLocalMappingApply()`, `ISLocalToGlobalMappingCreate()`, 973db781477SPatrick Sanan `ISLocalToGlobalMappingDestroy()` 9749d90f715SBarry Smith @*/ 975d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping, ISGlobalToLocalMappingMode type, PetscInt n, const PetscInt idx[], PetscInt *nout, PetscInt idxout[]) 976d71ae5a4SJacob Faibussowitsch { 9773a40ed3dSBarry Smith PetscFunctionBegin; 9780700a824SBarry Smith PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 97948a46eb9SPierre Jolivet if (!mapping->data) PetscCall(ISGlobalToLocalMappingSetUp(mapping)); 980dbbe0bcdSBarry Smith PetscUseTypeMethod(mapping, globaltolocalmappingapplyblock, type, n, idx, nout, idxout); 9813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 982d4bb536fSBarry Smith } 98390f02eecSBarry Smith 98489d82c54SBarry Smith /*@C 9856a818285SBarry Smith ISLocalToGlobalMappingGetBlockInfo - Gets the neighbor information for each processor and 98689d82c54SBarry Smith each index shared by more than one processor 98789d82c54SBarry Smith 988c3339decSBarry Smith Collective 98989d82c54SBarry Smith 990f899ff85SJose E. Roman Input Parameter: 99189d82c54SBarry Smith . mapping - the mapping from local to global indexing 99289d82c54SBarry Smith 993d8d19677SJose E. Roman Output Parameters: 99489d82c54SBarry Smith + nproc - number of processors that are connected to this one 99538b5cf2dSJacob Faibussowitsch . procs - neighboring processors 99638b5cf2dSJacob Faibussowitsch . numprocs - number of indices for each subdomain (processor) 9973463a7baSJed Brown - indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering) 99889d82c54SBarry Smith 99989d82c54SBarry Smith Level: advanced 100089d82c54SBarry Smith 100138b5cf2dSJacob Faibussowitsch Fortran Notes: 1002e33f79d8SJacob Faibussowitsch There is no `ISLocalToGlobalMappingRestoreInfo()` in Fortran. You must make sure that 1003e33f79d8SJacob Faibussowitsch `procs`[], `numprocs`[] and `indices`[][] are large enough arrays, either by allocating them 1004e33f79d8SJacob Faibussowitsch dynamically or defining static ones large enough. 1005cab54364SBarry Smith .vb 10062cfcea29SBarry Smith PetscInt indices[nproc][numprocmax],ierr) 1007cab54364SBarry Smith ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by 1008cab54364SBarry Smith ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc], 1009cab54364SBarry Smith .ve 1010cab54364SBarry Smith 1011cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 1012db781477SPatrick Sanan `ISLocalToGlobalMappingRestoreInfo()` 101389d82c54SBarry Smith @*/ 1014d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[]) 1015d71ae5a4SJacob Faibussowitsch { 1016268a049cSStefano Zampini PetscFunctionBegin; 1017268a049cSStefano Zampini PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 1018268a049cSStefano Zampini if (mapping->info_cached) { 1019268a049cSStefano Zampini *nproc = mapping->info_nproc; 1020268a049cSStefano Zampini *procs = mapping->info_procs; 1021268a049cSStefano Zampini *numprocs = mapping->info_numprocs; 1022268a049cSStefano Zampini *indices = mapping->info_indices; 1023268a049cSStefano Zampini } else { 10249566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockInfo_Private(mapping, nproc, procs, numprocs, indices)); 1025268a049cSStefano Zampini } 10263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1027268a049cSStefano Zampini } 1028268a049cSStefano Zampini 1029d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[]) 1030d71ae5a4SJacob Faibussowitsch { 103197f1f81fSBarry Smith PetscMPIInt size, rank, tag1, tag2, tag3, *len, *source, imdex; 103232dcc486SBarry Smith PetscInt i, n = mapping->n, Ng, ng, max = 0, *lindices = mapping->indices; 103332dcc486SBarry Smith PetscInt *nprocs, *owner, nsends, *sends, j, *starts, nmax, nrecvs, *recvs, proc; 1034c599c493SJunchao Zhang PetscInt cnt, scale, *ownedsenders, *nownedsenders, rstart; 103532dcc486SBarry Smith PetscInt node, nownedm, nt, *sends2, nsends2, *starts2, *lens2, *dest, nrecvs2, *starts3, *recvs2, k, *bprocs, *tmp; 103632dcc486SBarry Smith PetscInt first_procs, first_numprocs, *first_indices; 103789d82c54SBarry Smith MPI_Request *recv_waits, *send_waits; 103830dcb7c9SBarry Smith MPI_Status recv_status, *send_status, *recv_statuses; 1039ce94432eSBarry Smith MPI_Comm comm; 1040ace3abfcSBarry Smith PetscBool debug = PETSC_FALSE; 104189d82c54SBarry Smith 104289d82c54SBarry Smith PetscFunctionBegin; 10439566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)mapping, &comm)); 10449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 10459566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 104624cf384cSBarry Smith if (size == 1) { 104724cf384cSBarry Smith *nproc = 0; 10480298fd71SBarry Smith *procs = NULL; 10499566063dSJacob Faibussowitsch PetscCall(PetscNew(numprocs)); 10501e2105dcSBarry Smith (*numprocs)[0] = 0; 10519566063dSJacob Faibussowitsch PetscCall(PetscNew(indices)); 10520298fd71SBarry Smith (*indices)[0] = NULL; 1053268a049cSStefano Zampini /* save info for reuse */ 1054268a049cSStefano Zampini mapping->info_nproc = *nproc; 1055268a049cSStefano Zampini mapping->info_procs = *procs; 1056268a049cSStefano Zampini mapping->info_numprocs = *numprocs; 1057268a049cSStefano Zampini mapping->info_indices = *indices; 1058268a049cSStefano Zampini mapping->info_cached = PETSC_TRUE; 10593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 106024cf384cSBarry Smith } 106124cf384cSBarry Smith 10629566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)mapping)->options, NULL, "-islocaltoglobalmappinggetinfo_debug", &debug, NULL)); 106307b52d57SBarry Smith 10643677ff5aSBarry Smith /* 10656a818285SBarry Smith Notes on ISLocalToGlobalMappingGetBlockInfo 10663677ff5aSBarry Smith 10673677ff5aSBarry Smith globally owned node - the nodes that have been assigned to this processor in global 10683677ff5aSBarry Smith numbering, just for this routine. 10693677ff5aSBarry Smith 10703677ff5aSBarry Smith nontrivial globally owned node - node assigned to this processor that is on a subdomain 10713677ff5aSBarry Smith boundary (i.e. is has more than one local owner) 10723677ff5aSBarry Smith 10733677ff5aSBarry Smith locally owned node - node that exists on this processors subdomain 10743677ff5aSBarry Smith 10753677ff5aSBarry Smith nontrivial locally owned node - node that is not in the interior (i.e. has more than one 10763677ff5aSBarry Smith local subdomain 10773677ff5aSBarry Smith */ 10789566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)mapping, &tag1)); 10799566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)mapping, &tag2)); 10809566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)mapping, &tag3)); 108189d82c54SBarry Smith 108289d82c54SBarry Smith for (i = 0; i < n; i++) { 108389d82c54SBarry Smith if (lindices[i] > max) max = lindices[i]; 108489d82c54SBarry Smith } 10851c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&max, &Ng, 1, MPIU_INT, MPI_MAX, comm)); 108678058e43SBarry Smith Ng++; 10879566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 10889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1089bc8ff85bSBarry Smith scale = Ng / size + 1; 10909371c9d4SSatish Balay ng = scale; 10919371c9d4SSatish Balay if (rank == size - 1) ng = Ng - scale * (size - 1); 10929371c9d4SSatish Balay ng = PetscMax(1, ng); 1093caba0dd0SBarry Smith rstart = scale * rank; 109489d82c54SBarry Smith 109589d82c54SBarry Smith /* determine ownership ranges of global indices */ 10969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * size, &nprocs)); 10979566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(nprocs, 2 * size)); 109889d82c54SBarry Smith 109989d82c54SBarry Smith /* determine owners of each local node */ 11009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &owner)); 110189d82c54SBarry Smith for (i = 0; i < n; i++) { 11023677ff5aSBarry Smith proc = lindices[i] / scale; /* processor that globally owns this index */ 110327c402fcSBarry Smith nprocs[2 * proc + 1] = 1; /* processor globally owns at least one of ours */ 11043677ff5aSBarry Smith owner[i] = proc; 110527c402fcSBarry Smith nprocs[2 * proc]++; /* count of how many that processor globally owns of ours */ 110689d82c54SBarry Smith } 11079371c9d4SSatish Balay nsends = 0; 11089371c9d4SSatish Balay for (i = 0; i < size; i++) nsends += nprocs[2 * i + 1]; 11099566063dSJacob Faibussowitsch PetscCall(PetscInfo(mapping, "Number of global owners for my local data %" PetscInt_FMT "\n", nsends)); 111089d82c54SBarry Smith 111189d82c54SBarry Smith /* inform other processors of number of messages and max length*/ 11129566063dSJacob Faibussowitsch PetscCall(PetscMaxSum(comm, nprocs, &nmax, &nrecvs)); 11139566063dSJacob Faibussowitsch PetscCall(PetscInfo(mapping, "Number of local owners for my global data %" PetscInt_FMT "\n", nrecvs)); 111489d82c54SBarry Smith 111589d82c54SBarry Smith /* post receives for owned rows */ 11169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((2 * nrecvs + 1) * (nmax + 1), &recvs)); 11179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs + 1, &recv_waits)); 111848a46eb9SPierre 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)); 111989d82c54SBarry Smith 112089d82c54SBarry Smith /* pack messages containing lists of local nodes to owners */ 11219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * n + 1, &sends)); 11229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size + 1, &starts)); 112389d82c54SBarry Smith starts[0] = 0; 1124f6e5521dSKarl Rupp for (i = 1; i < size; i++) starts[i] = starts[i - 1] + 2 * nprocs[2 * i - 2]; 112589d82c54SBarry Smith for (i = 0; i < n; i++) { 112689d82c54SBarry Smith sends[starts[owner[i]]++] = lindices[i]; 112730dcb7c9SBarry Smith sends[starts[owner[i]]++] = i; 112889d82c54SBarry Smith } 11299566063dSJacob Faibussowitsch PetscCall(PetscFree(owner)); 113089d82c54SBarry Smith starts[0] = 0; 1131f6e5521dSKarl Rupp for (i = 1; i < size; i++) starts[i] = starts[i - 1] + 2 * nprocs[2 * i - 2]; 113289d82c54SBarry Smith 113389d82c54SBarry Smith /* send the messages */ 11349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends + 1, &send_waits)); 11359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends + 1, &dest)); 113689d82c54SBarry Smith cnt = 0; 113789d82c54SBarry Smith for (i = 0; i < size; i++) { 113827c402fcSBarry Smith if (nprocs[2 * i]) { 11399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Isend(sends + starts[i], 2 * nprocs[2 * i], MPIU_INT, i, tag1, comm, send_waits + cnt)); 114030dcb7c9SBarry Smith dest[cnt] = i; 114189d82c54SBarry Smith cnt++; 114289d82c54SBarry Smith } 114389d82c54SBarry Smith } 11449566063dSJacob Faibussowitsch PetscCall(PetscFree(starts)); 114589d82c54SBarry Smith 114689d82c54SBarry Smith /* wait on receives */ 11479566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs + 1, &source)); 11489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs + 1, &len)); 114989d82c54SBarry Smith cnt = nrecvs; 11509566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(ng + 1, &nownedsenders)); 115189d82c54SBarry Smith while (cnt) { 11529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitany(nrecvs, recv_waits, &imdex, &recv_status)); 115389d82c54SBarry Smith /* unpack receives into our local space */ 11549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&recv_status, MPIU_INT, &len[imdex])); 115589d82c54SBarry Smith source[imdex] = recv_status.MPI_SOURCE; 115630dcb7c9SBarry Smith len[imdex] = len[imdex] / 2; 1157caba0dd0SBarry Smith /* count how many local owners for each of my global owned indices */ 115830dcb7c9SBarry Smith for (i = 0; i < len[imdex]; i++) nownedsenders[recvs[2 * imdex * nmax + 2 * i] - rstart]++; 115989d82c54SBarry Smith cnt--; 116089d82c54SBarry Smith } 11619566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_waits)); 116289d82c54SBarry Smith 116330dcb7c9SBarry Smith /* count how many globally owned indices are on an edge multiplied by how many processors own them. */ 1164bc8ff85bSBarry Smith nownedm = 0; 1165bc8ff85bSBarry Smith for (i = 0; i < ng; i++) { 1166c599c493SJunchao Zhang if (nownedsenders[i] > 1) nownedm += nownedsenders[i]; 1167bc8ff85bSBarry Smith } 1168bc8ff85bSBarry Smith 1169bc8ff85bSBarry Smith /* create single array to contain rank of all local owners of each globally owned index */ 11709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nownedm + 1, &ownedsenders)); 11719566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ng + 1, &starts)); 1172bc8ff85bSBarry Smith starts[0] = 0; 1173bc8ff85bSBarry Smith for (i = 1; i < ng; i++) { 1174bc8ff85bSBarry Smith if (nownedsenders[i - 1] > 1) starts[i] = starts[i - 1] + nownedsenders[i - 1]; 1175bc8ff85bSBarry Smith else starts[i] = starts[i - 1]; 1176bc8ff85bSBarry Smith } 1177bc8ff85bSBarry Smith 11786aad120cSJose E. Roman /* for each nontrivial globally owned node list all arriving processors */ 1179bc8ff85bSBarry Smith for (i = 0; i < nrecvs; i++) { 1180bc8ff85bSBarry Smith for (j = 0; j < len[i]; j++) { 118130dcb7c9SBarry Smith node = recvs[2 * i * nmax + 2 * j] - rstart; 1182f6e5521dSKarl Rupp if (nownedsenders[node] > 1) ownedsenders[starts[node]++] = source[i]; 1183bc8ff85bSBarry Smith } 1184bc8ff85bSBarry Smith } 1185bc8ff85bSBarry Smith 118607b52d57SBarry Smith if (debug) { /* ----------------------------------- */ 118730dcb7c9SBarry Smith starts[0] = 0; 118830dcb7c9SBarry Smith for (i = 1; i < ng; i++) { 118930dcb7c9SBarry Smith if (nownedsenders[i - 1] > 1) starts[i] = starts[i - 1] + nownedsenders[i - 1]; 119030dcb7c9SBarry Smith else starts[i] = starts[i - 1]; 119130dcb7c9SBarry Smith } 119230dcb7c9SBarry Smith for (i = 0; i < ng; i++) { 119330dcb7c9SBarry Smith if (nownedsenders[i] > 1) { 11949566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm, "[%d] global node %" PetscInt_FMT " local owner processors: ", rank, i + rstart)); 119548a46eb9SPierre Jolivet for (j = 0; j < nownedsenders[i]; j++) PetscCall(PetscSynchronizedPrintf(comm, "%" PetscInt_FMT " ", ownedsenders[starts[i] + j])); 11969566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm, "\n")); 119730dcb7c9SBarry Smith } 119830dcb7c9SBarry Smith } 11999566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 120007b52d57SBarry Smith } /* ----------------------------------- */ 120130dcb7c9SBarry Smith 12023677ff5aSBarry Smith /* wait on original sends */ 12033a96401aSBarry Smith if (nsends) { 12049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends, &send_status)); 12059566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nsends, send_waits, send_status)); 12069566063dSJacob Faibussowitsch PetscCall(PetscFree(send_status)); 12073a96401aSBarry Smith } 12089566063dSJacob Faibussowitsch PetscCall(PetscFree(send_waits)); 12099566063dSJacob Faibussowitsch PetscCall(PetscFree(sends)); 12109566063dSJacob Faibussowitsch PetscCall(PetscFree(nprocs)); 12113677ff5aSBarry Smith 12123677ff5aSBarry Smith /* pack messages to send back to local owners */ 121330dcb7c9SBarry Smith starts[0] = 0; 121430dcb7c9SBarry Smith for (i = 1; i < ng; i++) { 121530dcb7c9SBarry Smith if (nownedsenders[i - 1] > 1) starts[i] = starts[i - 1] + nownedsenders[i - 1]; 121630dcb7c9SBarry Smith else starts[i] = starts[i - 1]; 121730dcb7c9SBarry Smith } 121830dcb7c9SBarry Smith nsends2 = nrecvs; 12199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends2 + 1, &nprocs)); /* length of each message */ 122030dcb7c9SBarry Smith for (i = 0; i < nrecvs; i++) { 122130dcb7c9SBarry Smith nprocs[i] = 1; 122230dcb7c9SBarry Smith for (j = 0; j < len[i]; j++) { 122330dcb7c9SBarry Smith node = recvs[2 * i * nmax + 2 * j] - rstart; 1224f6e5521dSKarl Rupp if (nownedsenders[node] > 1) nprocs[i] += 2 + nownedsenders[node]; 122530dcb7c9SBarry Smith } 122630dcb7c9SBarry Smith } 1227f6e5521dSKarl Rupp nt = 0; 1228f6e5521dSKarl Rupp for (i = 0; i < nsends2; i++) nt += nprocs[i]; 1229f6e5521dSKarl Rupp 12309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nt + 1, &sends2)); 12319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends2 + 1, &starts2)); 1232f6e5521dSKarl Rupp 1233f6e5521dSKarl Rupp starts2[0] = 0; 1234f6e5521dSKarl Rupp for (i = 1; i < nsends2; i++) starts2[i] = starts2[i - 1] + nprocs[i - 1]; 123530dcb7c9SBarry Smith /* 123630dcb7c9SBarry Smith Each message is 1 + nprocs[i] long, and consists of 123730dcb7c9SBarry Smith (0) the number of nodes being sent back 123830dcb7c9SBarry Smith (1) the local node number, 123930dcb7c9SBarry Smith (2) the number of processors sharing it, 124030dcb7c9SBarry Smith (3) the processors sharing it 124130dcb7c9SBarry Smith */ 124230dcb7c9SBarry Smith for (i = 0; i < nsends2; i++) { 124330dcb7c9SBarry Smith cnt = 1; 124430dcb7c9SBarry Smith sends2[starts2[i]] = 0; 124530dcb7c9SBarry Smith for (j = 0; j < len[i]; j++) { 124630dcb7c9SBarry Smith node = recvs[2 * i * nmax + 2 * j] - rstart; 124730dcb7c9SBarry Smith if (nownedsenders[node] > 1) { 124830dcb7c9SBarry Smith sends2[starts2[i]]++; 124930dcb7c9SBarry Smith sends2[starts2[i] + cnt++] = recvs[2 * i * nmax + 2 * j + 1]; 125030dcb7c9SBarry Smith sends2[starts2[i] + cnt++] = nownedsenders[node]; 12519566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(&sends2[starts2[i] + cnt], &ownedsenders[starts[node]], nownedsenders[node])); 125230dcb7c9SBarry Smith cnt += nownedsenders[node]; 125330dcb7c9SBarry Smith } 125430dcb7c9SBarry Smith } 125530dcb7c9SBarry Smith } 125630dcb7c9SBarry Smith 125730dcb7c9SBarry Smith /* receive the message lengths */ 125830dcb7c9SBarry Smith nrecvs2 = nsends; 12599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2 + 1, &lens2)); 12609566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2 + 1, &starts3)); 12619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2 + 1, &recv_waits)); 126248a46eb9SPierre Jolivet for (i = 0; i < nrecvs2; i++) PetscCallMPI(MPI_Irecv(&lens2[i], 1, MPIU_INT, dest[i], tag2, comm, recv_waits + i)); 1263d44834fbSBarry Smith 12648a8e0b3aSBarry Smith /* send the message lengths */ 126548a46eb9SPierre Jolivet for (i = 0; i < nsends2; i++) PetscCallMPI(MPI_Send(&nprocs[i], 1, MPIU_INT, source[i], tag2, comm)); 12668a8e0b3aSBarry Smith 1267d44834fbSBarry Smith /* wait on receives of lens */ 12680c468ba9SBarry Smith if (nrecvs2) { 12699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2, &recv_statuses)); 12709566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrecvs2, recv_waits, recv_statuses)); 12719566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_statuses)); 12720c468ba9SBarry Smith } 12739566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_waits)); 1274d44834fbSBarry Smith 127530dcb7c9SBarry Smith starts3[0] = 0; 1276d44834fbSBarry Smith nt = 0; 127730dcb7c9SBarry Smith for (i = 0; i < nrecvs2 - 1; i++) { 127830dcb7c9SBarry Smith starts3[i + 1] = starts3[i] + lens2[i]; 1279d44834fbSBarry Smith nt += lens2[i]; 128030dcb7c9SBarry Smith } 128176466f69SStefano Zampini if (nrecvs2) nt += lens2[nrecvs2 - 1]; 1282d44834fbSBarry Smith 12839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nt + 1, &recvs2)); 12849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2 + 1, &recv_waits)); 128548a46eb9SPierre Jolivet for (i = 0; i < nrecvs2; i++) PetscCallMPI(MPI_Irecv(recvs2 + starts3[i], lens2[i], MPIU_INT, dest[i], tag3, comm, recv_waits + i)); 128630dcb7c9SBarry Smith 128730dcb7c9SBarry Smith /* send the messages */ 12889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends2 + 1, &send_waits)); 128948a46eb9SPierre Jolivet for (i = 0; i < nsends2; i++) PetscCallMPI(MPI_Isend(sends2 + starts2[i], nprocs[i], MPIU_INT, source[i], tag3, comm, send_waits + i)); 129030dcb7c9SBarry Smith 129130dcb7c9SBarry Smith /* wait on receives */ 12920c468ba9SBarry Smith if (nrecvs2) { 12939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2, &recv_statuses)); 12949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrecvs2, recv_waits, recv_statuses)); 12959566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_statuses)); 12960c468ba9SBarry Smith } 12979566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_waits)); 12989566063dSJacob Faibussowitsch PetscCall(PetscFree(nprocs)); 129930dcb7c9SBarry Smith 130007b52d57SBarry Smith if (debug) { /* ----------------------------------- */ 130130dcb7c9SBarry Smith cnt = 0; 130230dcb7c9SBarry Smith for (i = 0; i < nrecvs2; i++) { 130330dcb7c9SBarry Smith nt = recvs2[cnt++]; 130430dcb7c9SBarry Smith for (j = 0; j < nt; j++) { 13059566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm, "[%d] local node %" PetscInt_FMT " number of subdomains %" PetscInt_FMT ": ", rank, recvs2[cnt], recvs2[cnt + 1])); 130648a46eb9SPierre Jolivet for (k = 0; k < recvs2[cnt + 1]; k++) PetscCall(PetscSynchronizedPrintf(comm, "%" PetscInt_FMT " ", recvs2[cnt + 2 + k])); 130730dcb7c9SBarry Smith cnt += 2 + recvs2[cnt + 1]; 13089566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm, "\n")); 130930dcb7c9SBarry Smith } 131030dcb7c9SBarry Smith } 13119566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 131207b52d57SBarry Smith } /* ----------------------------------- */ 131330dcb7c9SBarry Smith 131430dcb7c9SBarry Smith /* count number subdomains for each local node */ 13159566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(size, &nprocs)); 131630dcb7c9SBarry Smith cnt = 0; 131730dcb7c9SBarry Smith for (i = 0; i < nrecvs2; i++) { 131830dcb7c9SBarry Smith nt = recvs2[cnt++]; 131930dcb7c9SBarry Smith for (j = 0; j < nt; j++) { 1320f6e5521dSKarl Rupp for (k = 0; k < recvs2[cnt + 1]; k++) nprocs[recvs2[cnt + 2 + k]]++; 132130dcb7c9SBarry Smith cnt += 2 + recvs2[cnt + 1]; 132230dcb7c9SBarry Smith } 132330dcb7c9SBarry Smith } 13249371c9d4SSatish Balay nt = 0; 13259371c9d4SSatish Balay for (i = 0; i < size; i++) nt += (nprocs[i] > 0); 132630dcb7c9SBarry Smith *nproc = nt; 13279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nt + 1, procs)); 13289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nt + 1, numprocs)); 13299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nt + 1, indices)); 13300298fd71SBarry Smith for (i = 0; i < nt + 1; i++) (*indices)[i] = NULL; 13319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &bprocs)); 133230dcb7c9SBarry Smith cnt = 0; 133330dcb7c9SBarry Smith for (i = 0; i < size; i++) { 133430dcb7c9SBarry Smith if (nprocs[i] > 0) { 133530dcb7c9SBarry Smith bprocs[i] = cnt; 133630dcb7c9SBarry Smith (*procs)[cnt] = i; 133730dcb7c9SBarry Smith (*numprocs)[cnt] = nprocs[i]; 13389566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nprocs[i], &(*indices)[cnt])); 133930dcb7c9SBarry Smith cnt++; 134030dcb7c9SBarry Smith } 134130dcb7c9SBarry Smith } 134230dcb7c9SBarry Smith 134330dcb7c9SBarry Smith /* make the list of subdomains for each nontrivial local node */ 13449566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(*numprocs, nt)); 134530dcb7c9SBarry Smith cnt = 0; 134630dcb7c9SBarry Smith for (i = 0; i < nrecvs2; i++) { 134730dcb7c9SBarry Smith nt = recvs2[cnt++]; 134830dcb7c9SBarry Smith for (j = 0; j < nt; j++) { 1349f6e5521dSKarl Rupp for (k = 0; k < recvs2[cnt + 1]; k++) (*indices)[bprocs[recvs2[cnt + 2 + k]]][(*numprocs)[bprocs[recvs2[cnt + 2 + k]]]++] = recvs2[cnt]; 135030dcb7c9SBarry Smith cnt += 2 + recvs2[cnt + 1]; 135130dcb7c9SBarry Smith } 135230dcb7c9SBarry Smith } 13539566063dSJacob Faibussowitsch PetscCall(PetscFree(bprocs)); 13549566063dSJacob Faibussowitsch PetscCall(PetscFree(recvs2)); 135530dcb7c9SBarry Smith 135607b52d57SBarry Smith /* sort the node indexing by their global numbers */ 135707b52d57SBarry Smith nt = *nproc; 135807b52d57SBarry Smith for (i = 0; i < nt; i++) { 13599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((*numprocs)[i], &tmp)); 1360f6e5521dSKarl Rupp for (j = 0; j < (*numprocs)[i]; j++) tmp[j] = lindices[(*indices)[i][j]]; 13619566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithArray((*numprocs)[i], tmp, (*indices)[i])); 13629566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp)); 136307b52d57SBarry Smith } 136407b52d57SBarry Smith 136507b52d57SBarry Smith if (debug) { /* ----------------------------------- */ 136630dcb7c9SBarry Smith nt = *nproc; 136730dcb7c9SBarry Smith for (i = 0; i < nt; i++) { 13689566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm, "[%d] subdomain %" PetscInt_FMT " number of indices %" PetscInt_FMT ": ", rank, (*procs)[i], (*numprocs)[i])); 136948a46eb9SPierre Jolivet for (j = 0; j < (*numprocs)[i]; j++) PetscCall(PetscSynchronizedPrintf(comm, "%" PetscInt_FMT " ", (*indices)[i][j])); 13709566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm, "\n")); 137130dcb7c9SBarry Smith } 13729566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 137307b52d57SBarry Smith } /* ----------------------------------- */ 137430dcb7c9SBarry Smith 137530dcb7c9SBarry Smith /* wait on sends */ 137630dcb7c9SBarry Smith if (nsends2) { 13779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends2, &send_status)); 13789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nsends2, send_waits, send_status)); 13799566063dSJacob Faibussowitsch PetscCall(PetscFree(send_status)); 138030dcb7c9SBarry Smith } 138130dcb7c9SBarry Smith 13829566063dSJacob Faibussowitsch PetscCall(PetscFree(starts3)); 13839566063dSJacob Faibussowitsch PetscCall(PetscFree(dest)); 13849566063dSJacob Faibussowitsch PetscCall(PetscFree(send_waits)); 13853677ff5aSBarry Smith 13869566063dSJacob Faibussowitsch PetscCall(PetscFree(nownedsenders)); 13879566063dSJacob Faibussowitsch PetscCall(PetscFree(ownedsenders)); 13889566063dSJacob Faibussowitsch PetscCall(PetscFree(starts)); 13899566063dSJacob Faibussowitsch PetscCall(PetscFree(starts2)); 13909566063dSJacob Faibussowitsch PetscCall(PetscFree(lens2)); 139189d82c54SBarry Smith 13929566063dSJacob Faibussowitsch PetscCall(PetscFree(source)); 13939566063dSJacob Faibussowitsch PetscCall(PetscFree(len)); 13949566063dSJacob Faibussowitsch PetscCall(PetscFree(recvs)); 13959566063dSJacob Faibussowitsch PetscCall(PetscFree(nprocs)); 13969566063dSJacob Faibussowitsch PetscCall(PetscFree(sends2)); 139724cf384cSBarry Smith 139824cf384cSBarry Smith /* put the information about myself as the first entry in the list */ 139924cf384cSBarry Smith first_procs = (*procs)[0]; 140024cf384cSBarry Smith first_numprocs = (*numprocs)[0]; 140124cf384cSBarry Smith first_indices = (*indices)[0]; 140224cf384cSBarry Smith for (i = 0; i < *nproc; i++) { 140324cf384cSBarry Smith if ((*procs)[i] == rank) { 140424cf384cSBarry Smith (*procs)[0] = (*procs)[i]; 140524cf384cSBarry Smith (*numprocs)[0] = (*numprocs)[i]; 140624cf384cSBarry Smith (*indices)[0] = (*indices)[i]; 140724cf384cSBarry Smith (*procs)[i] = first_procs; 140824cf384cSBarry Smith (*numprocs)[i] = first_numprocs; 140924cf384cSBarry Smith (*indices)[i] = first_indices; 141024cf384cSBarry Smith break; 141124cf384cSBarry Smith } 141224cf384cSBarry Smith } 1413268a049cSStefano Zampini 1414268a049cSStefano Zampini /* save info for reuse */ 1415268a049cSStefano Zampini mapping->info_nproc = *nproc; 1416268a049cSStefano Zampini mapping->info_procs = *procs; 1417268a049cSStefano Zampini mapping->info_numprocs = *numprocs; 1418268a049cSStefano Zampini mapping->info_indices = *indices; 1419268a049cSStefano Zampini mapping->info_cached = PETSC_TRUE; 14203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 142189d82c54SBarry Smith } 142289d82c54SBarry Smith 14236a818285SBarry Smith /*@C 1424cab54364SBarry Smith ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by `ISLocalToGlobalMappingGetBlockInfo()` 14256a818285SBarry Smith 1426c3339decSBarry Smith Collective 14276a818285SBarry Smith 1428f899ff85SJose E. Roman Input Parameter: 14296a818285SBarry Smith . mapping - the mapping from local to global indexing 14306a818285SBarry Smith 1431d8d19677SJose E. Roman Output Parameters: 14326a818285SBarry Smith + nproc - number of processors that are connected to this one 143338b5cf2dSJacob Faibussowitsch . procs - neighboring processors 143438b5cf2dSJacob Faibussowitsch . numprocs - number of indices for each processor 14356a818285SBarry Smith - indices - indices of local nodes shared with neighbor (sorted by global numbering) 14366a818285SBarry Smith 14376a818285SBarry Smith Level: advanced 14386a818285SBarry Smith 1439cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 1440db781477SPatrick Sanan `ISLocalToGlobalMappingGetInfo()` 14416a818285SBarry Smith @*/ 1442d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[]) 1443d71ae5a4SJacob Faibussowitsch { 14446a818285SBarry Smith PetscFunctionBegin; 1445cbc1caf0SMatthew G. Knepley PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 1446268a049cSStefano Zampini if (mapping->info_free) { 14479566063dSJacob Faibussowitsch PetscCall(PetscFree(*numprocs)); 14486a818285SBarry Smith if (*indices) { 1449268a049cSStefano Zampini PetscInt i; 1450268a049cSStefano Zampini 14519566063dSJacob Faibussowitsch PetscCall(PetscFree((*indices)[0])); 145248a46eb9SPierre Jolivet for (i = 1; i < *nproc; i++) PetscCall(PetscFree((*indices)[i])); 14539566063dSJacob Faibussowitsch PetscCall(PetscFree(*indices)); 14546a818285SBarry Smith } 1455268a049cSStefano Zampini } 1456268a049cSStefano Zampini *nproc = 0; 1457268a049cSStefano Zampini *procs = NULL; 1458268a049cSStefano Zampini *numprocs = NULL; 1459268a049cSStefano Zampini *indices = NULL; 14603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14616a818285SBarry Smith } 14626a818285SBarry Smith 14636a818285SBarry Smith /*@C 14646a818285SBarry Smith ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each processor and 14656a818285SBarry Smith each index shared by more than one processor 14666a818285SBarry Smith 1467c3339decSBarry Smith Collective 14686a818285SBarry Smith 1469f899ff85SJose E. Roman Input Parameter: 14706a818285SBarry Smith . mapping - the mapping from local to global indexing 14716a818285SBarry Smith 1472d8d19677SJose E. Roman Output Parameters: 14736a818285SBarry Smith + nproc - number of processors that are connected to this one 147438b5cf2dSJacob Faibussowitsch . procs - neighboring processors 147538b5cf2dSJacob Faibussowitsch . numprocs - number of indices for each subdomain (processor) 14766a818285SBarry Smith - indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering) 14776a818285SBarry Smith 14786a818285SBarry Smith Level: advanced 14796a818285SBarry Smith 1480cab54364SBarry Smith Note: 1481cab54364SBarry Smith The user needs to call `ISLocalToGlobalMappingRestoreInfo()` when the data is no longer needed. 14821bd0b88eSStefano Zampini 148338b5cf2dSJacob Faibussowitsch Fortran Notes: 1484e33f79d8SJacob Faibussowitsch There is no `ISLocalToGlobalMappingRestoreInfo()` in Fortran. You must make sure that 1485e33f79d8SJacob Faibussowitsch `procs`[], `numprocs`[] and `indices`[][] are large enough arrays, either by allocating them 1486e33f79d8SJacob Faibussowitsch dynamically or defining static ones large enough. 1487e33f79d8SJacob Faibussowitsch 1488cab54364SBarry Smith .vb 14896a818285SBarry Smith PetscInt indices[nproc][numprocmax],ierr) 1490cab54364SBarry Smith ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by 1491cab54364SBarry Smith ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc], 1492cab54364SBarry Smith .ve 1493cab54364SBarry Smith 1494cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 1495db781477SPatrick Sanan `ISLocalToGlobalMappingRestoreInfo()` 14966a818285SBarry Smith @*/ 1497d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[]) 1498d71ae5a4SJacob Faibussowitsch { 14998a1f772fSStefano Zampini PetscInt **bindices = NULL, *bnumprocs = NULL, bs, i, j, k; 15006a818285SBarry Smith 15016a818285SBarry Smith PetscFunctionBegin; 15026a818285SBarry Smith PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 15038a1f772fSStefano Zampini bs = mapping->bs; 15049566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockInfo(mapping, nproc, procs, &bnumprocs, &bindices)); 1505268a049cSStefano Zampini if (bs > 1) { /* we need to expand the cached info */ 15069566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(*nproc, &*indices)); 15079566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(*nproc, &*numprocs)); 15086a818285SBarry Smith for (i = 0; i < *nproc; i++) { 15099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * bnumprocs[i], &(*indices)[i])); 1510268a049cSStefano Zampini for (j = 0; j < bnumprocs[i]; j++) { 1511ad540459SPierre Jolivet for (k = 0; k < bs; k++) (*indices)[i][j * bs + k] = bs * bindices[i][j] + k; 15126a818285SBarry Smith } 1513268a049cSStefano Zampini (*numprocs)[i] = bnumprocs[i] * bs; 15146a818285SBarry Smith } 1515268a049cSStefano Zampini mapping->info_free = PETSC_TRUE; 1516268a049cSStefano Zampini } else { 1517268a049cSStefano Zampini *numprocs = bnumprocs; 1518268a049cSStefano Zampini *indices = bindices; 15196a818285SBarry Smith } 15203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15216a818285SBarry Smith } 15226a818285SBarry Smith 152307b52d57SBarry Smith /*@C 1524cab54364SBarry Smith ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by `ISLocalToGlobalMappingGetInfo()` 152589d82c54SBarry Smith 1526c3339decSBarry Smith Collective 152707b52d57SBarry Smith 1528f899ff85SJose E. Roman Input Parameter: 152907b52d57SBarry Smith . mapping - the mapping from local to global indexing 153007b52d57SBarry Smith 1531d8d19677SJose E. Roman Output Parameters: 153207b52d57SBarry Smith + nproc - number of processors that are connected to this one 153338b5cf2dSJacob Faibussowitsch . procs - neighboring processors 153438b5cf2dSJacob Faibussowitsch . numprocs - number of indices for each processor 153507b52d57SBarry Smith - indices - indices of local nodes shared with neighbor (sorted by global numbering) 153607b52d57SBarry Smith 153707b52d57SBarry Smith Level: advanced 153807b52d57SBarry Smith 1539cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 1540db781477SPatrick Sanan `ISLocalToGlobalMappingGetInfo()` 154107b52d57SBarry Smith @*/ 1542d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping, PetscInt *nproc, PetscInt *procs[], PetscInt *numprocs[], PetscInt **indices[]) 1543d71ae5a4SJacob Faibussowitsch { 154407b52d57SBarry Smith PetscFunctionBegin; 15459566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockInfo(mapping, nproc, procs, numprocs, indices)); 15463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 154707b52d57SBarry Smith } 154886994e45SJed Brown 154986994e45SJed Brown /*@C 1550cab54364SBarry Smith ISLocalToGlobalMappingGetNodeInfo - Gets the neighbor information for each MPI rank 15511bd0b88eSStefano Zampini 1552c3339decSBarry Smith Collective 15531bd0b88eSStefano Zampini 1554f899ff85SJose E. Roman Input Parameter: 15551bd0b88eSStefano Zampini . mapping - the mapping from local to global indexing 15561bd0b88eSStefano Zampini 1557d8d19677SJose E. Roman Output Parameters: 1558cab54364SBarry Smith + nnodes - number of local nodes (same `ISLocalToGlobalMappingGetSize()`) 15591bd0b88eSStefano Zampini . count - number of neighboring processors per node 15601bd0b88eSStefano Zampini - indices - indices of processes sharing the node (sorted) 15611bd0b88eSStefano Zampini 15621bd0b88eSStefano Zampini Level: advanced 15631bd0b88eSStefano Zampini 1564cab54364SBarry Smith Note: 1565cab54364SBarry Smith The user needs to call `ISLocalToGlobalMappingRestoreInfo()` when the data is no longer needed. 15661bd0b88eSStefano Zampini 1567cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 1568db781477SPatrick Sanan `ISLocalToGlobalMappingGetInfo()`, `ISLocalToGlobalMappingRestoreNodeInfo()` 15691bd0b88eSStefano Zampini @*/ 1570d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetNodeInfo(ISLocalToGlobalMapping mapping, PetscInt *nnodes, PetscInt *count[], PetscInt **indices[]) 1571d71ae5a4SJacob Faibussowitsch { 15721bd0b88eSStefano Zampini PetscInt n; 15731bd0b88eSStefano Zampini 15741bd0b88eSStefano Zampini PetscFunctionBegin; 15751bd0b88eSStefano Zampini PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 15769566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(mapping, &n)); 15771bd0b88eSStefano Zampini if (!mapping->info_nodec) { 15781bd0b88eSStefano Zampini PetscInt i, m, n_neigh, *neigh, *n_shared, **shared; 15791bd0b88eSStefano Zampini 15809566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n + 1, &mapping->info_nodec, n, &mapping->info_nodei)); 15819566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetInfo(mapping, &n_neigh, &neigh, &n_shared, &shared)); 1582ad540459SPierre Jolivet for (i = 0; i < n; i++) mapping->info_nodec[i] = 1; 1583071fcb05SBarry Smith m = n; 1584071fcb05SBarry Smith mapping->info_nodec[n] = 0; 15851bd0b88eSStefano Zampini for (i = 1; i < n_neigh; i++) { 15861bd0b88eSStefano Zampini PetscInt j; 15871bd0b88eSStefano Zampini 15881bd0b88eSStefano Zampini m += n_shared[i]; 15891bd0b88eSStefano Zampini for (j = 0; j < n_shared[i]; j++) mapping->info_nodec[shared[i][j]] += 1; 15901bd0b88eSStefano Zampini } 15919566063dSJacob Faibussowitsch if (n) PetscCall(PetscMalloc1(m, &mapping->info_nodei[0])); 15921bd0b88eSStefano Zampini for (i = 1; i < n; i++) mapping->info_nodei[i] = mapping->info_nodei[i - 1] + mapping->info_nodec[i - 1]; 15939566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(mapping->info_nodec, n)); 15949371c9d4SSatish Balay for (i = 0; i < n; i++) { 15959371c9d4SSatish Balay mapping->info_nodec[i] = 1; 15969371c9d4SSatish Balay mapping->info_nodei[i][0] = neigh[0]; 15979371c9d4SSatish Balay } 15981bd0b88eSStefano Zampini for (i = 1; i < n_neigh; i++) { 15991bd0b88eSStefano Zampini PetscInt j; 16001bd0b88eSStefano Zampini 16011bd0b88eSStefano Zampini for (j = 0; j < n_shared[i]; j++) { 16021bd0b88eSStefano Zampini PetscInt k = shared[i][j]; 16031bd0b88eSStefano Zampini 16041bd0b88eSStefano Zampini mapping->info_nodei[k][mapping->info_nodec[k]] = neigh[i]; 16051bd0b88eSStefano Zampini mapping->info_nodec[k] += 1; 16061bd0b88eSStefano Zampini } 16071bd0b88eSStefano Zampini } 16089566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(PetscSortRemoveDupsInt(&mapping->info_nodec[i], mapping->info_nodei[i])); 16099566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreInfo(mapping, &n_neigh, &neigh, &n_shared, &shared)); 16101bd0b88eSStefano Zampini } 16111bd0b88eSStefano Zampini if (nnodes) *nnodes = n; 16121bd0b88eSStefano Zampini if (count) *count = mapping->info_nodec; 16131bd0b88eSStefano Zampini if (indices) *indices = mapping->info_nodei; 16143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16151bd0b88eSStefano Zampini } 16161bd0b88eSStefano Zampini 16171bd0b88eSStefano Zampini /*@C 1618cab54364SBarry Smith ISLocalToGlobalMappingRestoreNodeInfo - Frees the memory allocated by `ISLocalToGlobalMappingGetNodeInfo()` 16191bd0b88eSStefano Zampini 1620c3339decSBarry Smith Collective 16211bd0b88eSStefano Zampini 1622f899ff85SJose E. Roman Input Parameter: 16231bd0b88eSStefano Zampini . mapping - the mapping from local to global indexing 16241bd0b88eSStefano Zampini 1625d8d19677SJose E. Roman Output Parameters: 16261bd0b88eSStefano Zampini + nnodes - number of local nodes 16271bd0b88eSStefano Zampini . count - number of neighboring processors per node 16281bd0b88eSStefano Zampini - indices - indices of processes sharing the node (sorted) 16291bd0b88eSStefano Zampini 16301bd0b88eSStefano Zampini Level: advanced 16311bd0b88eSStefano Zampini 1632cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 1633db781477SPatrick Sanan `ISLocalToGlobalMappingGetInfo()` 16341bd0b88eSStefano Zampini @*/ 1635d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreNodeInfo(ISLocalToGlobalMapping mapping, PetscInt *nnodes, PetscInt *count[], PetscInt **indices[]) 1636d71ae5a4SJacob Faibussowitsch { 16371bd0b88eSStefano Zampini PetscFunctionBegin; 16381bd0b88eSStefano Zampini PetscValidHeaderSpecific(mapping, IS_LTOGM_CLASSID, 1); 16391bd0b88eSStefano Zampini if (nnodes) *nnodes = 0; 16401bd0b88eSStefano Zampini if (count) *count = NULL; 16411bd0b88eSStefano Zampini if (indices) *indices = NULL; 16423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16431bd0b88eSStefano Zampini } 16441bd0b88eSStefano Zampini 16456ce40d10SBarry Smith /*MC 16466ce40d10SBarry Smith ISLocalToGlobalMappingGetIndicesF90 - Get global indices for every local point that is mapped 16476ce40d10SBarry Smith 16486ce40d10SBarry Smith Synopsis: 16496ce40d10SBarry Smith ISLocalToGlobalMappingGetIndicesF90(ISLocalToGlobalMapping ltog, PetscInt, pointer :: array(:)}, integer ierr) 16506ce40d10SBarry Smith 16516ce40d10SBarry Smith Not Collective 16526ce40d10SBarry Smith 16536ce40d10SBarry Smith Input Parameter: 16546ce40d10SBarry Smith . A - the matrix 16556ce40d10SBarry Smith 16562fe279fdSBarry Smith Output Parameter: 16576ce40d10SBarry Smith . array - array of indices, the length of this array may be obtained with `ISLocalToGlobalMappingGetSize()` 16586ce40d10SBarry Smith 16596ce40d10SBarry Smith Level: advanced 16606ce40d10SBarry Smith 16616ce40d10SBarry Smith Note: 16626ce40d10SBarry Smith Use `ISLocalToGlobalMappingGetIndicesF90()` when you no longer need access to the data 16636ce40d10SBarry Smith 166420662ed9SBarry Smith .seealso: [](sec_fortranarrays), [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingGetIndices()`, `ISLocalToGlobalMappingRestoreIndices()`, 166520662ed9SBarry Smith `ISLocalToGlobalMappingRestoreIndicesF90()` 16666ce40d10SBarry Smith M*/ 16676ce40d10SBarry Smith 16686ce40d10SBarry Smith /*MC 16696ce40d10SBarry Smith ISLocalToGlobalMappingRestoreIndicesF90 - restores the global indices for every local point that is mapped obtained with `ISLocalToGlobalMappingGetIndicesF90()` 16706ce40d10SBarry Smith 16716ce40d10SBarry Smith Synopsis: 16726ce40d10SBarry Smith ISLocalToGlobalMappingRestoreIndicesF90(ISLocalToGlobalMapping ltog, PetscInt, pointer :: array(:)}, integer ierr) 16736ce40d10SBarry Smith 16746ce40d10SBarry Smith Not Collective 16756ce40d10SBarry Smith 16766ce40d10SBarry Smith Input Parameters: 16776ce40d10SBarry Smith + A - the matrix 16786ce40d10SBarry Smith - array - array of indices, the length of this array may be obtained with `ISLocalToGlobalMappingGetSize()` 16796ce40d10SBarry Smith 16806ce40d10SBarry Smith Level: advanced 16816ce40d10SBarry Smith 168220662ed9SBarry Smith .seealso: [](sec_fortranarrays), [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingGetIndices()`, `ISLocalToGlobalMappingRestoreIndices()`, 168320662ed9SBarry Smith `ISLocalToGlobalMappingGetIndicesF90()` 16846ce40d10SBarry Smith M*/ 16856ce40d10SBarry Smith 16861bd0b88eSStefano Zampini /*@C 1687107e9a97SBarry Smith ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped 168886994e45SJed Brown 168986994e45SJed Brown Not Collective 169086994e45SJed Brown 16914165533cSJose E. Roman Input Parameter: 169286994e45SJed Brown . ltog - local to global mapping 169386994e45SJed Brown 16944165533cSJose E. Roman Output Parameter: 1695cab54364SBarry Smith . array - array of indices, the length of this array may be obtained with `ISLocalToGlobalMappingGetSize()` 169686994e45SJed Brown 169786994e45SJed Brown Level: advanced 169886994e45SJed Brown 1699cab54364SBarry Smith Note: 1700cab54364SBarry Smith `ISLocalToGlobalMappingGetSize()` returns the length the this array 1701107e9a97SBarry Smith 170220662ed9SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingRestoreIndices()`, 170320662ed9SBarry Smith `ISLocalToGlobalMappingGetBlockIndices()`, `ISLocalToGlobalMappingRestoreBlockIndices()` 170486994e45SJed Brown @*/ 1705d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog, const PetscInt **array) 1706d71ae5a4SJacob Faibussowitsch { 170786994e45SJed Brown PetscFunctionBegin; 170886994e45SJed Brown PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1); 17094f572ea9SToby Isaac PetscAssertPointer(array, 2); 171045b6f7e9SBarry Smith if (ltog->bs == 1) { 171186994e45SJed Brown *array = ltog->indices; 171245b6f7e9SBarry Smith } else { 171345b6f7e9SBarry Smith PetscInt *jj, k, i, j, n = ltog->n, bs = ltog->bs; 171445b6f7e9SBarry Smith const PetscInt *ii; 171545b6f7e9SBarry Smith 17169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * n, &jj)); 171745b6f7e9SBarry Smith *array = jj; 171845b6f7e9SBarry Smith k = 0; 171945b6f7e9SBarry Smith ii = ltog->indices; 172045b6f7e9SBarry Smith for (i = 0; i < n; i++) 17219371c9d4SSatish Balay for (j = 0; j < bs; j++) jj[k++] = bs * ii[i] + j; 172245b6f7e9SBarry Smith } 17233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 172486994e45SJed Brown } 172586994e45SJed Brown 172686994e45SJed Brown /*@C 1727cab54364SBarry Smith ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with `ISLocalToGlobalMappingGetIndices()` 172886994e45SJed Brown 172986994e45SJed Brown Not Collective 173086994e45SJed Brown 17314165533cSJose E. Roman Input Parameters: 173286994e45SJed Brown + ltog - local to global mapping 173386994e45SJed Brown - array - array of indices 173486994e45SJed Brown 173586994e45SJed Brown Level: advanced 173686994e45SJed Brown 1737cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingGetIndices()` 173886994e45SJed Brown @*/ 1739d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog, const PetscInt **array) 1740d71ae5a4SJacob Faibussowitsch { 174186994e45SJed Brown PetscFunctionBegin; 174286994e45SJed Brown PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1); 17434f572ea9SToby Isaac PetscAssertPointer(array, 2); 1744c9cc58a2SBarry Smith PetscCheck(ltog->bs != 1 || *array == ltog->indices, PETSC_COMM_SELF, PETSC_ERR_ARG_BADPTR, "Trying to return mismatched pointer"); 174545b6f7e9SBarry Smith 174648a46eb9SPierre Jolivet if (ltog->bs > 1) PetscCall(PetscFree(*(void **)array)); 17473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 174845b6f7e9SBarry Smith } 174945b6f7e9SBarry Smith 175045b6f7e9SBarry Smith /*@C 175145b6f7e9SBarry Smith ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block 175245b6f7e9SBarry Smith 175345b6f7e9SBarry Smith Not Collective 175445b6f7e9SBarry Smith 17554165533cSJose E. Roman Input Parameter: 175645b6f7e9SBarry Smith . ltog - local to global mapping 175745b6f7e9SBarry Smith 17584165533cSJose E. Roman Output Parameter: 175945b6f7e9SBarry Smith . array - array of indices 176045b6f7e9SBarry Smith 176145b6f7e9SBarry Smith Level: advanced 176245b6f7e9SBarry Smith 1763cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingRestoreBlockIndices()` 176445b6f7e9SBarry Smith @*/ 1765d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog, const PetscInt **array) 1766d71ae5a4SJacob Faibussowitsch { 176745b6f7e9SBarry Smith PetscFunctionBegin; 176845b6f7e9SBarry Smith PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1); 17694f572ea9SToby Isaac PetscAssertPointer(array, 2); 177045b6f7e9SBarry Smith *array = ltog->indices; 17713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 177245b6f7e9SBarry Smith } 177345b6f7e9SBarry Smith 177445b6f7e9SBarry Smith /*@C 1775cab54364SBarry Smith ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with `ISLocalToGlobalMappingGetBlockIndices()` 177645b6f7e9SBarry Smith 177745b6f7e9SBarry Smith Not Collective 177845b6f7e9SBarry Smith 17794165533cSJose E. Roman Input Parameters: 178045b6f7e9SBarry Smith + ltog - local to global mapping 178145b6f7e9SBarry Smith - array - array of indices 178245b6f7e9SBarry Smith 178345b6f7e9SBarry Smith Level: advanced 178445b6f7e9SBarry Smith 1785cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingGetIndices()` 178645b6f7e9SBarry Smith @*/ 1787d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog, const PetscInt **array) 1788d71ae5a4SJacob Faibussowitsch { 178945b6f7e9SBarry Smith PetscFunctionBegin; 179045b6f7e9SBarry Smith PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1); 17914f572ea9SToby Isaac PetscAssertPointer(array, 2); 179208401ef6SPierre Jolivet PetscCheck(*array == ltog->indices, PETSC_COMM_SELF, PETSC_ERR_ARG_BADPTR, "Trying to return mismatched pointer"); 17930298fd71SBarry Smith *array = NULL; 17943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 179586994e45SJed Brown } 1796f7efa3c7SJed Brown 1797f7efa3c7SJed Brown /*@C 1798f7efa3c7SJed Brown ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings 1799f7efa3c7SJed Brown 1800f7efa3c7SJed Brown Not Collective 1801f7efa3c7SJed Brown 18024165533cSJose E. Roman Input Parameters: 1803f7efa3c7SJed Brown + comm - communicator for the new mapping, must contain the communicator of every mapping to concatenate 1804f7efa3c7SJed Brown . n - number of mappings to concatenate 1805f7efa3c7SJed Brown - ltogs - local to global mappings 1806f7efa3c7SJed Brown 18074165533cSJose E. Roman Output Parameter: 1808f7efa3c7SJed Brown . ltogcat - new mapping 1809f7efa3c7SJed Brown 1810f7efa3c7SJed Brown Level: advanced 1811f7efa3c7SJed Brown 1812cab54364SBarry Smith Note: 1813cab54364SBarry Smith This currently always returns a mapping with block size of 1 1814cab54364SBarry Smith 181538b5cf2dSJacob Faibussowitsch Developer Notes: 1816cab54364SBarry Smith If all the input mapping have the same block size we could easily handle that as a special case 1817cab54364SBarry Smith 1818cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingCreate()` 1819f7efa3c7SJed Brown @*/ 1820d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm, PetscInt n, const ISLocalToGlobalMapping ltogs[], ISLocalToGlobalMapping *ltogcat) 1821d71ae5a4SJacob Faibussowitsch { 1822f7efa3c7SJed Brown PetscInt i, cnt, m, *idx; 1823f7efa3c7SJed Brown 1824f7efa3c7SJed Brown PetscFunctionBegin; 182508401ef6SPierre Jolivet PetscCheck(n >= 0, comm, PETSC_ERR_ARG_OUTOFRANGE, "Must have a non-negative number of mappings, given %" PetscInt_FMT, n); 18264f572ea9SToby Isaac if (n > 0) PetscAssertPointer(ltogs, 3); 1827f7efa3c7SJed Brown for (i = 0; i < n; i++) PetscValidHeaderSpecific(ltogs[i], IS_LTOGM_CLASSID, 3); 18284f572ea9SToby Isaac PetscAssertPointer(ltogcat, 4); 1829f7efa3c7SJed Brown for (cnt = 0, i = 0; i < n; i++) { 18309566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(ltogs[i], &m)); 1831f7efa3c7SJed Brown cnt += m; 1832f7efa3c7SJed Brown } 18339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(cnt, &idx)); 1834f7efa3c7SJed Brown for (cnt = 0, i = 0; i < n; i++) { 1835f7efa3c7SJed Brown const PetscInt *subidx; 18369566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(ltogs[i], &m)); 18379566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(ltogs[i], &subidx)); 18389566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(&idx[cnt], subidx, m)); 18399566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(ltogs[i], &subidx)); 1840f7efa3c7SJed Brown cnt += m; 1841f7efa3c7SJed Brown } 18429566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm, 1, cnt, idx, PETSC_OWN_POINTER, ltogcat)); 18433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1844f7efa3c7SJed Brown } 184504a59952SBarry Smith 1846413f72f0SBarry Smith /*MC 1847cab54364SBarry Smith ISLOCALTOGLOBALMAPPINGBASIC - basic implementation of the `ISLocalToGlobalMapping` object. When `ISGlobalToLocalMappingApply()` is 1848413f72f0SBarry Smith used this is good for only small and moderate size problems. 1849413f72f0SBarry Smith 185020662ed9SBarry Smith Options Database Key: 1851a2b725a8SWilliam Gropp . -islocaltoglobalmapping_type basic - select this method 1852413f72f0SBarry Smith 1853413f72f0SBarry Smith Level: beginner 1854413f72f0SBarry Smith 1855cab54364SBarry Smith Developer Note: 1856cab54364SBarry Smith This stores all the mapping information on each MPI rank. 1857cab54364SBarry Smith 1858cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`, `ISLOCALTOGLOBALMAPPINGHASH` 1859413f72f0SBarry Smith M*/ 1860d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Basic(ISLocalToGlobalMapping ltog) 1861d71ae5a4SJacob Faibussowitsch { 1862413f72f0SBarry Smith PetscFunctionBegin; 1863413f72f0SBarry Smith ltog->ops->globaltolocalmappingapply = ISGlobalToLocalMappingApply_Basic; 1864413f72f0SBarry Smith ltog->ops->globaltolocalmappingsetup = ISGlobalToLocalMappingSetUp_Basic; 1865413f72f0SBarry Smith ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Basic; 1866413f72f0SBarry Smith ltog->ops->destroy = ISLocalToGlobalMappingDestroy_Basic; 18673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1868413f72f0SBarry Smith } 1869413f72f0SBarry Smith 1870413f72f0SBarry Smith /*MC 1871cab54364SBarry Smith ISLOCALTOGLOBALMAPPINGHASH - hash implementation of the `ISLocalToGlobalMapping` object. When `ISGlobalToLocalMappingApply()` is 1872ed56e8eaSBarry Smith used this is good for large memory problems. 1873413f72f0SBarry Smith 187420662ed9SBarry Smith Options Database Key: 1875a2b725a8SWilliam Gropp . -islocaltoglobalmapping_type hash - select this method 1876413f72f0SBarry Smith 1877413f72f0SBarry Smith Level: beginner 1878413f72f0SBarry Smith 1879cab54364SBarry Smith Note: 1880cab54364SBarry Smith This is selected automatically for large problems if the user does not set the type. 1881cab54364SBarry Smith 1882cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`, `ISLOCALTOGLOBALMAPPINGBASIC` 1883413f72f0SBarry Smith M*/ 1884d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Hash(ISLocalToGlobalMapping ltog) 1885d71ae5a4SJacob Faibussowitsch { 1886413f72f0SBarry Smith PetscFunctionBegin; 1887413f72f0SBarry Smith ltog->ops->globaltolocalmappingapply = ISGlobalToLocalMappingApply_Hash; 1888413f72f0SBarry Smith ltog->ops->globaltolocalmappingsetup = ISGlobalToLocalMappingSetUp_Hash; 1889413f72f0SBarry Smith ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Hash; 1890413f72f0SBarry Smith ltog->ops->destroy = ISLocalToGlobalMappingDestroy_Hash; 18913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1892413f72f0SBarry Smith } 1893413f72f0SBarry Smith 1894413f72f0SBarry Smith /*@C 1895cab54364SBarry Smith ISLocalToGlobalMappingRegister - Registers a method for applying a global to local mapping with an `ISLocalToGlobalMapping` 1896413f72f0SBarry Smith 1897413f72f0SBarry Smith Not Collective 1898413f72f0SBarry Smith 1899413f72f0SBarry Smith Input Parameters: 1900413f72f0SBarry Smith + sname - name of a new method 19012fe279fdSBarry Smith - function - routine to create method context 1902413f72f0SBarry Smith 190338b5cf2dSJacob Faibussowitsch Example Usage: 1904413f72f0SBarry Smith .vb 1905413f72f0SBarry Smith ISLocalToGlobalMappingRegister("my_mapper", MyCreate); 1906413f72f0SBarry Smith .ve 1907413f72f0SBarry Smith 1908ed56e8eaSBarry Smith Then, your mapping can be chosen with the procedural interface via 1909413f72f0SBarry Smith $ ISLocalToGlobalMappingSetType(ltog, "my_mapper") 1910413f72f0SBarry Smith or at runtime via the option 1911ed56e8eaSBarry Smith $ -islocaltoglobalmapping_type my_mapper 1912413f72f0SBarry Smith 1913413f72f0SBarry Smith Level: advanced 1914413f72f0SBarry Smith 1915cab54364SBarry Smith Note: 1916cab54364SBarry Smith `ISLocalToGlobalMappingRegister()` may be called multiple times to add several user-defined mappings. 1917413f72f0SBarry Smith 1918cab54364SBarry Smith .seealso: [](sec_scatter), `ISLocalToGlobalMappingRegisterAll()`, `ISLocalToGlobalMappingRegisterDestroy()`, `ISLOCALTOGLOBALMAPPINGBASIC`, 1919cab54364SBarry Smith `ISLOCALTOGLOBALMAPPINGHASH`, `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingApply()` 1920413f72f0SBarry Smith @*/ 1921d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRegister(const char sname[], PetscErrorCode (*function)(ISLocalToGlobalMapping)) 1922d71ae5a4SJacob Faibussowitsch { 1923413f72f0SBarry Smith PetscFunctionBegin; 19249566063dSJacob Faibussowitsch PetscCall(ISInitializePackage()); 19259566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&ISLocalToGlobalMappingList, sname, function)); 19263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1927413f72f0SBarry Smith } 1928413f72f0SBarry Smith 1929413f72f0SBarry Smith /*@C 1930cab54364SBarry Smith ISLocalToGlobalMappingSetType - Sets the implementation type `ISLocalToGlobalMapping` will use 1931413f72f0SBarry Smith 1932c3339decSBarry Smith Logically Collective 1933413f72f0SBarry Smith 1934413f72f0SBarry Smith Input Parameters: 1935cab54364SBarry Smith + ltog - the `ISLocalToGlobalMapping` object 1936413f72f0SBarry Smith - type - a known method 1937413f72f0SBarry Smith 1938413f72f0SBarry Smith Options Database Key: 1939cab54364SBarry Smith . -islocaltoglobalmapping_type <method> - Sets the method; use -help for a list of available methods (for instance, basic or hash) 1940cab54364SBarry Smith 1941cab54364SBarry Smith Level: intermediate 1942413f72f0SBarry Smith 1943413f72f0SBarry Smith Notes: 194420662ed9SBarry Smith See `ISLocalToGlobalMappingType` for available methods 1945413f72f0SBarry Smith 1946cab54364SBarry Smith Normally, it is best to use the `ISLocalToGlobalMappingSetFromOptions()` command and 1947cab54364SBarry Smith then set the `ISLocalToGlobalMappingType` from the options database rather than by using 1948413f72f0SBarry Smith this routine. 1949413f72f0SBarry Smith 195038b5cf2dSJacob Faibussowitsch Developer Notes: 1951cab54364SBarry Smith `ISLocalToGlobalMappingRegister()` is used to add new types to `ISLocalToGlobalMappingList` from which they 1952cab54364SBarry Smith are accessed by `ISLocalToGlobalMappingSetType()`. 1953413f72f0SBarry Smith 195438b5cf2dSJacob Faibussowitsch .seealso: [](sec_scatter), `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingRegister()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingGetType()` 1955413f72f0SBarry Smith @*/ 1956d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType type) 1957d71ae5a4SJacob Faibussowitsch { 1958413f72f0SBarry Smith PetscBool match; 19595f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(ISLocalToGlobalMapping) = NULL; 1960413f72f0SBarry Smith 1961413f72f0SBarry Smith PetscFunctionBegin; 1962413f72f0SBarry Smith PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1); 19634f572ea9SToby Isaac if (type) PetscAssertPointer(type, 2); 1964413f72f0SBarry Smith 19659566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)ltog, type, &match)); 19663ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS); 1967413f72f0SBarry Smith 1968a0d79125SStefano Zampini /* L2G maps defer type setup at globaltolocal calls, allow passing NULL here */ 1969a0d79125SStefano Zampini if (type) { 19709566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(ISLocalToGlobalMappingList, type, &r)); 1971a0d79125SStefano Zampini PetscCheck(r, PetscObjectComm((PetscObject)ltog), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested ISLocalToGlobalMapping type %s", type); 1972a0d79125SStefano Zampini } 1973413f72f0SBarry Smith /* Destroy the previous private LTOG context */ 1974dbbe0bcdSBarry Smith PetscTryTypeMethod(ltog, destroy); 1975413f72f0SBarry Smith ltog->ops->destroy = NULL; 1976dbbe0bcdSBarry Smith 19779566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)ltog, type)); 19789566063dSJacob Faibussowitsch if (r) PetscCall((*r)(ltog)); 19793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1980a0d79125SStefano Zampini } 1981a0d79125SStefano Zampini 1982a0d79125SStefano Zampini /*@C 1983cab54364SBarry Smith ISLocalToGlobalMappingGetType - Get the type of the `ISLocalToGlobalMapping` 1984a0d79125SStefano Zampini 1985a0d79125SStefano Zampini Not Collective 1986a0d79125SStefano Zampini 1987a0d79125SStefano Zampini Input Parameter: 1988cab54364SBarry Smith . ltog - the `ISLocalToGlobalMapping` object 1989a0d79125SStefano Zampini 1990a0d79125SStefano Zampini Output Parameter: 1991a0d79125SStefano Zampini . type - the type 1992a0d79125SStefano Zampini 199349762cbcSSatish Balay Level: intermediate 199449762cbcSSatish Balay 199538b5cf2dSJacob Faibussowitsch .seealso: [](sec_scatter), `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingRegister()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()` 1996a0d79125SStefano Zampini @*/ 1997d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingGetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType *type) 1998d71ae5a4SJacob Faibussowitsch { 1999a0d79125SStefano Zampini PetscFunctionBegin; 2000a0d79125SStefano Zampini PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 1); 20014f572ea9SToby Isaac PetscAssertPointer(type, 2); 2002a0d79125SStefano Zampini *type = ((PetscObject)ltog)->type_name; 20033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2004413f72f0SBarry Smith } 2005413f72f0SBarry Smith 2006413f72f0SBarry Smith PetscBool ISLocalToGlobalMappingRegisterAllCalled = PETSC_FALSE; 2007413f72f0SBarry Smith 2008413f72f0SBarry Smith /*@C 2009cab54364SBarry Smith ISLocalToGlobalMappingRegisterAll - Registers all of the local to global mapping components in the `IS` package. 2010413f72f0SBarry Smith 2011413f72f0SBarry Smith Not Collective 2012413f72f0SBarry Smith 2013413f72f0SBarry Smith Level: advanced 2014413f72f0SBarry Smith 2015cab54364SBarry Smith .seealso: [](sec_scatter), `ISRegister()`, `ISLocalToGlobalRegister()` 2016413f72f0SBarry Smith @*/ 2017d71ae5a4SJacob Faibussowitsch PetscErrorCode ISLocalToGlobalMappingRegisterAll(void) 2018d71ae5a4SJacob Faibussowitsch { 2019413f72f0SBarry Smith PetscFunctionBegin; 20203ba16761SJacob Faibussowitsch if (ISLocalToGlobalMappingRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS); 2021413f72f0SBarry Smith ISLocalToGlobalMappingRegisterAllCalled = PETSC_TRUE; 20229566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGBASIC, ISLocalToGlobalMappingCreate_Basic)); 20239566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGHASH, ISLocalToGlobalMappingCreate_Hash)); 20243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2025413f72f0SBarry Smith } 2026