12362add9SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/isimpl.h> /*I "petscis.h" I*/ 3e8f14785SLisandro Dalcin #include <petsc/private/hashmapi.h> 40c312b8eSJed Brown #include <petscsf.h> 5665c2dedSJed Brown #include <petscviewer.h> 62362add9SBarry Smith 77087cfbeSBarry Smith PetscClassId IS_LTOGM_CLASSID; 8268a049cSStefano Zampini static PetscErrorCode ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping,PetscInt*,PetscInt**,PetscInt**,PetscInt***); 98e58c17dSMatthew Knepley 10413f72f0SBarry Smith typedef struct { 11413f72f0SBarry Smith PetscInt *globals; 12413f72f0SBarry Smith } ISLocalToGlobalMapping_Basic; 13413f72f0SBarry Smith 14413f72f0SBarry Smith typedef struct { 15e8f14785SLisandro Dalcin PetscHMapI globalht; 16413f72f0SBarry Smith } ISLocalToGlobalMapping_Hash; 17413f72f0SBarry Smith 186528b96dSMatthew G. Knepley /*@C 196528b96dSMatthew G. Knepley ISGetPointRange - Returns a description of the points in an IS suitable for traversal 20413f72f0SBarry Smith 216528b96dSMatthew G. Knepley Not collective 226528b96dSMatthew G. Knepley 236528b96dSMatthew G. Knepley Input Parameter: 246528b96dSMatthew G. Knepley . pointIS - The IS object 256528b96dSMatthew G. Knepley 266528b96dSMatthew G. Knepley Output Parameters: 276528b96dSMatthew G. Knepley + pStart - The first index, see notes 286528b96dSMatthew G. Knepley . pEnd - One past the last index, see notes 296528b96dSMatthew G. Knepley - points - The indices, see notes 306528b96dSMatthew G. Knepley 316528b96dSMatthew G. Knepley Notes: 326528b96dSMatthew G. Knepley If the IS contains contiguous indices in an ISSTRIDE, then the indices are contained in [pStart, pEnd) and points = NULL. Otherwise, pStart = 0, pEnd = numIndices, and points is an array of the indices. This supports the following pattern 336528b96dSMatthew G. Knepley $ ISGetPointRange(is, &pStart, &pEnd, &points); 346528b96dSMatthew G. Knepley $ for (p = pStart; p < pEnd; ++p) { 356528b96dSMatthew G. Knepley $ const PetscInt point = points ? points[p] : p; 366528b96dSMatthew G. Knepley $ } 376528b96dSMatthew G. Knepley $ ISRestorePointRange(is, &pstart, &pEnd, &points); 386528b96dSMatthew G. Knepley 396528b96dSMatthew G. Knepley Level: intermediate 406528b96dSMatthew G. Knepley 41db781477SPatrick Sanan .seealso: `ISRestorePointRange()`, `ISGetPointSubrange()`, `ISGetIndices()`, `ISCreateStride()` 426528b96dSMatthew G. Knepley @*/ 439305a4c7SMatthew G. Knepley PetscErrorCode ISGetPointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points) 449305a4c7SMatthew G. Knepley { 459305a4c7SMatthew G. Knepley PetscInt numCells, step = 1; 469305a4c7SMatthew G. Knepley PetscBool isStride; 479305a4c7SMatthew G. Knepley 489305a4c7SMatthew G. Knepley PetscFunctionBeginHot; 499305a4c7SMatthew G. Knepley *pStart = 0; 509305a4c7SMatthew G. Knepley *points = NULL; 519566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(pointIS, &numCells)); 529566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject) pointIS, ISSTRIDE, &isStride)); 539566063dSJacob Faibussowitsch if (isStride) PetscCall(ISStrideGetInfo(pointIS, pStart, &step)); 549305a4c7SMatthew G. Knepley *pEnd = *pStart + numCells; 559566063dSJacob Faibussowitsch if (!isStride || step != 1) PetscCall(ISGetIndices(pointIS, points)); 569305a4c7SMatthew G. Knepley PetscFunctionReturn(0); 579305a4c7SMatthew G. Knepley } 589305a4c7SMatthew G. Knepley 596528b96dSMatthew G. Knepley /*@C 606528b96dSMatthew G. Knepley ISRestorePointRange - Destroys the traversal description 616528b96dSMatthew G. Knepley 626528b96dSMatthew G. Knepley Not collective 636528b96dSMatthew G. Knepley 646528b96dSMatthew G. Knepley Input Parameters: 656528b96dSMatthew G. Knepley + pointIS - The IS object 666528b96dSMatthew G. Knepley . pStart - The first index, from ISGetPointRange() 676528b96dSMatthew G. Knepley . pEnd - One past the last index, from ISGetPointRange() 686528b96dSMatthew G. Knepley - points - The indices, from ISGetPointRange() 696528b96dSMatthew G. Knepley 706528b96dSMatthew G. Knepley Notes: 716528b96dSMatthew G. Knepley If the IS contains contiguous indices in an ISSTRIDE, then the indices are contained in [pStart, pEnd) and points = NULL. Otherwise, pStart = 0, pEnd = numIndices, and points is an array of the indices. This supports the following pattern 726528b96dSMatthew G. Knepley $ ISGetPointRange(is, &pStart, &pEnd, &points); 736528b96dSMatthew G. Knepley $ for (p = pStart; p < pEnd; ++p) { 746528b96dSMatthew G. Knepley $ const PetscInt point = points ? points[p] : p; 756528b96dSMatthew G. Knepley $ } 766528b96dSMatthew G. Knepley $ ISRestorePointRange(is, &pstart, &pEnd, &points); 776528b96dSMatthew G. Knepley 786528b96dSMatthew G. Knepley Level: intermediate 796528b96dSMatthew G. Knepley 80db781477SPatrick Sanan .seealso: `ISGetPointRange()`, `ISGetPointSubrange()`, `ISGetIndices()`, `ISCreateStride()` 816528b96dSMatthew G. Knepley @*/ 829305a4c7SMatthew G. Knepley PetscErrorCode ISRestorePointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points) 839305a4c7SMatthew G. Knepley { 849305a4c7SMatthew G. Knepley PetscInt step = 1; 859305a4c7SMatthew G. Knepley PetscBool isStride; 869305a4c7SMatthew G. Knepley 879305a4c7SMatthew G. Knepley PetscFunctionBeginHot; 889566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject) pointIS, ISSTRIDE, &isStride)); 899566063dSJacob Faibussowitsch if (isStride) PetscCall(ISStrideGetInfo(pointIS, pStart, &step)); 909566063dSJacob Faibussowitsch if (!isStride || step != 1) PetscCall(ISGetIndices(pointIS, points)); 919305a4c7SMatthew G. Knepley PetscFunctionReturn(0); 929305a4c7SMatthew G. Knepley } 939305a4c7SMatthew G. Knepley 946528b96dSMatthew G. Knepley /*@C 956528b96dSMatthew G. Knepley ISGetPointSubrange - Configures the input IS to be a subrange for the traversal information given 966528b96dSMatthew G. Knepley 976528b96dSMatthew G. Knepley Not collective 986528b96dSMatthew G. Knepley 996528b96dSMatthew G. Knepley Input Parameters: 1006528b96dSMatthew G. Knepley + subpointIS - The IS object to be configured 1016528b96dSMatthew G. Knepley . pStar t - The first index of the subrange 1026528b96dSMatthew G. Knepley . pEnd - One past the last index for the subrange 1036528b96dSMatthew G. Knepley - points - The indices for the entire range, from ISGetPointRange() 1046528b96dSMatthew G. Knepley 1056528b96dSMatthew G. Knepley Output Parameters: 1066528b96dSMatthew G. Knepley . subpointIS - The IS object now configured to be a subrange 1076528b96dSMatthew G. Knepley 1086528b96dSMatthew G. Knepley Notes: 1096528b96dSMatthew G. Knepley The input IS will now respond properly to calls to ISGetPointRange() and return the subrange. 1106528b96dSMatthew G. Knepley 1116528b96dSMatthew G. Knepley Level: intermediate 1126528b96dSMatthew G. Knepley 113db781477SPatrick Sanan .seealso: `ISGetPointRange()`, `ISRestorePointRange()`, `ISGetIndices()`, `ISCreateStride()` 1146528b96dSMatthew G. Knepley @*/ 1159305a4c7SMatthew G. Knepley PetscErrorCode ISGetPointSubrange(IS subpointIS, PetscInt pStart, PetscInt pEnd, const PetscInt *points) 1169305a4c7SMatthew G. Knepley { 1179305a4c7SMatthew G. Knepley PetscFunctionBeginHot; 1189305a4c7SMatthew G. Knepley if (points) { 1199566063dSJacob Faibussowitsch PetscCall(ISSetType(subpointIS, ISGENERAL)); 1209566063dSJacob Faibussowitsch PetscCall(ISGeneralSetIndices(subpointIS, pEnd-pStart, &points[pStart], PETSC_USE_POINTER)); 1219305a4c7SMatthew G. Knepley } else { 1229566063dSJacob Faibussowitsch PetscCall(ISSetType(subpointIS, ISSTRIDE)); 1239566063dSJacob Faibussowitsch PetscCall(ISStrideSetStride(subpointIS, pEnd-pStart, pStart, 1)); 1249305a4c7SMatthew G. Knepley } 1259305a4c7SMatthew G. Knepley PetscFunctionReturn(0); 1269305a4c7SMatthew G. Knepley } 1279305a4c7SMatthew G. Knepley 128413f72f0SBarry Smith /* -----------------------------------------------------------------------------------------*/ 129413f72f0SBarry Smith 130413f72f0SBarry Smith /* 131413f72f0SBarry Smith Creates the global mapping information in the ISLocalToGlobalMapping structure 132413f72f0SBarry Smith 133413f72f0SBarry Smith If the user has not selected how to handle the global to local mapping then use HASH for "large" problems 134413f72f0SBarry Smith */ 135413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp(ISLocalToGlobalMapping mapping) 136413f72f0SBarry Smith { 137413f72f0SBarry Smith PetscInt i,*idx = mapping->indices,n = mapping->n,end,start; 138413f72f0SBarry Smith 139413f72f0SBarry Smith PetscFunctionBegin; 140413f72f0SBarry Smith if (mapping->data) PetscFunctionReturn(0); 141413f72f0SBarry Smith end = 0; 142413f72f0SBarry Smith start = PETSC_MAX_INT; 143413f72f0SBarry Smith 144413f72f0SBarry Smith for (i=0; i<n; i++) { 145413f72f0SBarry Smith if (idx[i] < 0) continue; 146413f72f0SBarry Smith if (idx[i] < start) start = idx[i]; 147413f72f0SBarry Smith if (idx[i] > end) end = idx[i]; 148413f72f0SBarry Smith } 149413f72f0SBarry Smith if (start > end) {start = 0; end = -1;} 150413f72f0SBarry Smith mapping->globalstart = start; 151413f72f0SBarry Smith mapping->globalend = end; 152413f72f0SBarry Smith if (!((PetscObject)mapping)->type_name) { 153413f72f0SBarry Smith if ((end - start) > PetscMax(4*n,1000000)) { 1549566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingSetType(mapping,ISLOCALTOGLOBALMAPPINGHASH)); 155413f72f0SBarry Smith } else { 1569566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingSetType(mapping,ISLOCALTOGLOBALMAPPINGBASIC)); 157413f72f0SBarry Smith } 158413f72f0SBarry Smith } 1599566063dSJacob Faibussowitsch if (mapping->ops->globaltolocalmappingsetup) PetscCall((*mapping->ops->globaltolocalmappingsetup)(mapping)); 160413f72f0SBarry Smith PetscFunctionReturn(0); 161413f72f0SBarry Smith } 162413f72f0SBarry Smith 163413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp_Basic(ISLocalToGlobalMapping mapping) 164413f72f0SBarry Smith { 165413f72f0SBarry Smith PetscInt i,*idx = mapping->indices,n = mapping->n,end,start,*globals; 166413f72f0SBarry Smith ISLocalToGlobalMapping_Basic *map; 167413f72f0SBarry Smith 168413f72f0SBarry Smith PetscFunctionBegin; 169413f72f0SBarry Smith start = mapping->globalstart; 170413f72f0SBarry Smith end = mapping->globalend; 1719566063dSJacob Faibussowitsch PetscCall(PetscNew(&map)); 1729566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(end-start+2,&globals)); 173413f72f0SBarry Smith map->globals = globals; 174413f72f0SBarry Smith for (i=0; i<end-start+1; i++) globals[i] = -1; 175413f72f0SBarry Smith for (i=0; i<n; i++) { 176413f72f0SBarry Smith if (idx[i] < 0) continue; 177413f72f0SBarry Smith globals[idx[i] - start] = i; 178413f72f0SBarry Smith } 179413f72f0SBarry Smith mapping->data = (void*)map; 1809566063dSJacob Faibussowitsch PetscCall(PetscLogObjectMemory((PetscObject)mapping,(end-start+1)*sizeof(PetscInt))); 181413f72f0SBarry Smith PetscFunctionReturn(0); 182413f72f0SBarry Smith } 183413f72f0SBarry Smith 184413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp_Hash(ISLocalToGlobalMapping mapping) 185413f72f0SBarry Smith { 186413f72f0SBarry Smith PetscInt i,*idx = mapping->indices,n = mapping->n; 187413f72f0SBarry Smith ISLocalToGlobalMapping_Hash *map; 188413f72f0SBarry Smith 189413f72f0SBarry Smith PetscFunctionBegin; 1909566063dSJacob Faibussowitsch PetscCall(PetscNew(&map)); 1919566063dSJacob Faibussowitsch PetscCall(PetscHMapICreate(&map->globalht)); 192413f72f0SBarry Smith for (i=0; i<n; i++) { 193413f72f0SBarry Smith if (idx[i] < 0) continue; 1949566063dSJacob Faibussowitsch PetscCall(PetscHMapISet(map->globalht,idx[i],i)); 195413f72f0SBarry Smith } 196413f72f0SBarry Smith mapping->data = (void*)map; 1979566063dSJacob Faibussowitsch PetscCall(PetscLogObjectMemory((PetscObject)mapping,2*n*sizeof(PetscInt))); 198413f72f0SBarry Smith PetscFunctionReturn(0); 199413f72f0SBarry Smith } 200413f72f0SBarry Smith 201413f72f0SBarry Smith static PetscErrorCode ISLocalToGlobalMappingDestroy_Basic(ISLocalToGlobalMapping mapping) 202413f72f0SBarry Smith { 203413f72f0SBarry Smith ISLocalToGlobalMapping_Basic *map = (ISLocalToGlobalMapping_Basic *)mapping->data; 204413f72f0SBarry Smith 205413f72f0SBarry Smith PetscFunctionBegin; 206413f72f0SBarry Smith if (!map) PetscFunctionReturn(0); 2079566063dSJacob Faibussowitsch PetscCall(PetscFree(map->globals)); 2089566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->data)); 209413f72f0SBarry Smith PetscFunctionReturn(0); 210413f72f0SBarry Smith } 211413f72f0SBarry Smith 212413f72f0SBarry Smith static PetscErrorCode ISLocalToGlobalMappingDestroy_Hash(ISLocalToGlobalMapping mapping) 213413f72f0SBarry Smith { 214413f72f0SBarry Smith ISLocalToGlobalMapping_Hash *map = (ISLocalToGlobalMapping_Hash*)mapping->data; 215413f72f0SBarry Smith 216413f72f0SBarry Smith PetscFunctionBegin; 217413f72f0SBarry Smith if (!map) PetscFunctionReturn(0); 2189566063dSJacob Faibussowitsch PetscCall(PetscHMapIDestroy(&map->globalht)); 2199566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->data)); 220413f72f0SBarry Smith PetscFunctionReturn(0); 221413f72f0SBarry Smith } 222413f72f0SBarry Smith 223413f72f0SBarry Smith #define GTOLTYPE _Basic 224413f72f0SBarry Smith #define GTOLNAME _Basic 225541bf97eSAdrian Croucher #define GTOLBS mapping->bs 226413f72f0SBarry Smith #define GTOL(g, local) do { \ 227413f72f0SBarry Smith local = map->globals[g/bs - start]; \ 2280040bde1SJunchao Zhang if (local >= 0) local = bs*local + (g % bs); \ 229413f72f0SBarry Smith } while (0) 230541bf97eSAdrian Croucher 231413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h> 232413f72f0SBarry Smith 233413f72f0SBarry Smith #define GTOLTYPE _Basic 234413f72f0SBarry Smith #define GTOLNAME Block_Basic 235541bf97eSAdrian Croucher #define GTOLBS 1 236413f72f0SBarry Smith #define GTOL(g, local) do { \ 237413f72f0SBarry Smith local = map->globals[g - start]; \ 238413f72f0SBarry Smith } while (0) 239413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h> 240413f72f0SBarry Smith 241413f72f0SBarry Smith #define GTOLTYPE _Hash 242413f72f0SBarry Smith #define GTOLNAME _Hash 243541bf97eSAdrian Croucher #define GTOLBS mapping->bs 244413f72f0SBarry Smith #define GTOL(g, local) do { \ 245e8f14785SLisandro Dalcin (void)PetscHMapIGet(map->globalht,g/bs,&local); \ 2460040bde1SJunchao Zhang if (local >= 0) local = bs*local + (g % bs); \ 247413f72f0SBarry Smith } while (0) 248413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h> 249413f72f0SBarry Smith 250413f72f0SBarry Smith #define GTOLTYPE _Hash 251413f72f0SBarry Smith #define GTOLNAME Block_Hash 252541bf97eSAdrian Croucher #define GTOLBS 1 253413f72f0SBarry Smith #define GTOL(g, local) do { \ 254e8f14785SLisandro Dalcin (void)PetscHMapIGet(map->globalht,g,&local); \ 255413f72f0SBarry Smith } while (0) 256413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h> 257413f72f0SBarry Smith 2586658fb44Sstefano_zampini /*@ 2596658fb44Sstefano_zampini ISLocalToGlobalMappingDuplicate - Duplicates the local to global mapping object 2606658fb44Sstefano_zampini 2616658fb44Sstefano_zampini Not Collective 2626658fb44Sstefano_zampini 2636658fb44Sstefano_zampini Input Parameter: 2646658fb44Sstefano_zampini . ltog - local to global mapping 2656658fb44Sstefano_zampini 2666658fb44Sstefano_zampini Output Parameter: 2676658fb44Sstefano_zampini . nltog - the duplicated local to global mapping 2686658fb44Sstefano_zampini 2696658fb44Sstefano_zampini Level: advanced 2706658fb44Sstefano_zampini 271db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()` 2726658fb44Sstefano_zampini @*/ 2736658fb44Sstefano_zampini PetscErrorCode ISLocalToGlobalMappingDuplicate(ISLocalToGlobalMapping ltog,ISLocalToGlobalMapping* nltog) 2746658fb44Sstefano_zampini { 275a0d79125SStefano Zampini ISLocalToGlobalMappingType l2gtype; 2766658fb44Sstefano_zampini 2776658fb44Sstefano_zampini PetscFunctionBegin; 2786658fb44Sstefano_zampini PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1); 2799566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)ltog),ltog->bs,ltog->n,ltog->indices,PETSC_COPY_VALUES,nltog)); 2809566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetType(ltog,&l2gtype)); 2819566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingSetType(*nltog,l2gtype)); 2826658fb44Sstefano_zampini PetscFunctionReturn(0); 2836658fb44Sstefano_zampini } 2846658fb44Sstefano_zampini 285565245c5SBarry Smith /*@ 286107e9a97SBarry Smith ISLocalToGlobalMappingGetSize - Gets the local size of a local to global mapping 2873b9aefa3SBarry Smith 2883b9aefa3SBarry Smith Not Collective 2893b9aefa3SBarry Smith 2903b9aefa3SBarry Smith Input Parameter: 2913b9aefa3SBarry Smith . ltog - local to global mapping 2923b9aefa3SBarry Smith 2933b9aefa3SBarry Smith Output Parameter: 294107e9a97SBarry Smith . n - the number of entries in the local mapping, ISLocalToGlobalMappingGetIndices() returns an array of this length 2953b9aefa3SBarry Smith 2963b9aefa3SBarry Smith Level: advanced 2973b9aefa3SBarry Smith 298db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()` 2993b9aefa3SBarry Smith @*/ 3007087cfbeSBarry Smith PetscErrorCode ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping mapping,PetscInt *n) 3013b9aefa3SBarry Smith { 3023b9aefa3SBarry Smith PetscFunctionBegin; 3030700a824SBarry Smith PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 3044482741eSBarry Smith PetscValidIntPointer(n,2); 305107e9a97SBarry Smith *n = mapping->bs*mapping->n; 3063b9aefa3SBarry Smith PetscFunctionReturn(0); 3073b9aefa3SBarry Smith } 3083b9aefa3SBarry Smith 3095a5d4f66SBarry Smith /*@C 310fe2efc57SMark ISLocalToGlobalMappingViewFromOptions - View from Options 311fe2efc57SMark 312fe2efc57SMark Collective on ISLocalToGlobalMapping 313fe2efc57SMark 314fe2efc57SMark Input Parameters: 315fe2efc57SMark + A - the local to global mapping object 316736c3998SJose E. Roman . obj - Optional object 317736c3998SJose E. Roman - name - command line option 318fe2efc57SMark 319fe2efc57SMark Level: intermediate 320db781477SPatrick Sanan .seealso: `ISLocalToGlobalMapping`, `ISLocalToGlobalMappingView`, `PetscObjectViewFromOptions()`, `ISLocalToGlobalMappingCreate()` 321fe2efc57SMark @*/ 322fe2efc57SMark PetscErrorCode ISLocalToGlobalMappingViewFromOptions(ISLocalToGlobalMapping A,PetscObject obj,const char name[]) 323fe2efc57SMark { 324fe2efc57SMark PetscFunctionBegin; 325fe2efc57SMark PetscValidHeaderSpecific(A,IS_LTOGM_CLASSID,1); 3269566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)A,obj,name)); 327fe2efc57SMark PetscFunctionReturn(0); 328fe2efc57SMark } 329fe2efc57SMark 330fe2efc57SMark /*@C 3315a5d4f66SBarry Smith ISLocalToGlobalMappingView - View a local to global mapping 3325a5d4f66SBarry Smith 333b9cd556bSLois Curfman McInnes Not Collective 334b9cd556bSLois Curfman McInnes 3355a5d4f66SBarry Smith Input Parameters: 3363b9aefa3SBarry Smith + ltog - local to global mapping 3373b9aefa3SBarry Smith - viewer - viewer 3385a5d4f66SBarry Smith 339a997ad1aSLois Curfman McInnes Level: advanced 340a997ad1aSLois Curfman McInnes 341db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()` 3425a5d4f66SBarry Smith @*/ 3437087cfbeSBarry Smith PetscErrorCode ISLocalToGlobalMappingView(ISLocalToGlobalMapping mapping,PetscViewer viewer) 3445a5d4f66SBarry Smith { 34532dcc486SBarry Smith PetscInt i; 34632dcc486SBarry Smith PetscMPIInt rank; 347ace3abfcSBarry Smith PetscBool iascii; 3485a5d4f66SBarry Smith 3495a5d4f66SBarry Smith PetscFunctionBegin; 3500700a824SBarry Smith PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 3513050cee2SBarry Smith if (!viewer) { 3529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mapping),&viewer)); 3533050cee2SBarry Smith } 3540700a824SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 3555a5d4f66SBarry Smith 3569566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mapping),&rank)); 3579566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 35832077d6dSBarry Smith if (iascii) { 3599566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)mapping,viewer)); 3609566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 3615a5d4f66SBarry Smith for (i=0; i<mapping->n; i++) { 3629566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %" PetscInt_FMT " %" PetscInt_FMT "\n",rank,i,mapping->indices[i])); 3636831982aSBarry Smith } 3649566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 3659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 3661575c14dSBarry Smith } 3675a5d4f66SBarry Smith PetscFunctionReturn(0); 3685a5d4f66SBarry Smith } 3695a5d4f66SBarry Smith 3701f428162SBarry Smith /*@ 3712bdab257SBarry Smith ISLocalToGlobalMappingCreateIS - Creates a mapping between a local (0 to n) 3722bdab257SBarry Smith ordering and a global parallel ordering. 3732bdab257SBarry Smith 3740f5bd95cSBarry Smith Not collective 375b9cd556bSLois Curfman McInnes 376a997ad1aSLois Curfman McInnes Input Parameter: 3778c03b21aSDmitry Karpeev . is - index set containing the global numbers for each local number 3782bdab257SBarry Smith 379a997ad1aSLois Curfman McInnes Output Parameter: 3802bdab257SBarry Smith . mapping - new mapping data structure 3812bdab257SBarry Smith 38295452b02SPatrick Sanan Notes: 38395452b02SPatrick Sanan the block size of the IS determines the block size of the mapping 384a997ad1aSLois Curfman McInnes Level: advanced 385a997ad1aSLois Curfman McInnes 386db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetFromOptions()` 3872bdab257SBarry Smith @*/ 3887087cfbeSBarry Smith PetscErrorCode ISLocalToGlobalMappingCreateIS(IS is,ISLocalToGlobalMapping *mapping) 3892bdab257SBarry Smith { 3903bbf0e92SBarry Smith PetscInt n,bs; 3915d0c19d7SBarry Smith const PetscInt *indices; 3922bdab257SBarry Smith MPI_Comm comm; 3933bbf0e92SBarry Smith PetscBool isblock; 3943a40ed3dSBarry Smith 3953a40ed3dSBarry Smith PetscFunctionBegin; 3960700a824SBarry Smith PetscValidHeaderSpecific(is,IS_CLASSID,1); 3974482741eSBarry Smith PetscValidPointer(mapping,2); 3982bdab257SBarry Smith 3999566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)is,&comm)); 4009566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is,&n)); 4019566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)is,ISBLOCK,&isblock)); 4026006e8d2SBarry Smith if (!isblock) { 4039566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is,&indices)); 4049566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm,1,n,indices,PETSC_COPY_VALUES,mapping)); 4059566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is,&indices)); 4066006e8d2SBarry Smith } else { 4079566063dSJacob Faibussowitsch PetscCall(ISGetBlockSize(is,&bs)); 4089566063dSJacob Faibussowitsch PetscCall(ISBlockGetIndices(is,&indices)); 4099566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm,bs,n/bs,indices,PETSC_COPY_VALUES,mapping)); 4109566063dSJacob Faibussowitsch PetscCall(ISBlockRestoreIndices(is,&indices)); 4116006e8d2SBarry Smith } 4123a40ed3dSBarry Smith PetscFunctionReturn(0); 4132bdab257SBarry Smith } 4145a5d4f66SBarry Smith 415a4d96a55SJed Brown /*@C 416a4d96a55SJed Brown ISLocalToGlobalMappingCreateSF - Creates a mapping between a local (0 to n) 417a4d96a55SJed Brown ordering and a global parallel ordering. 418a4d96a55SJed Brown 419a4d96a55SJed Brown Collective 420a4d96a55SJed Brown 421d8d19677SJose E. Roman Input Parameters: 422a4d96a55SJed Brown + sf - star forest mapping contiguous local indices to (rank, offset) 4239a535bafSVaclav Hapla - start - first global index on this process, or PETSC_DECIDE to compute contiguous global numbering automatically 424a4d96a55SJed Brown 425a4d96a55SJed Brown Output Parameter: 426a4d96a55SJed Brown . mapping - new mapping data structure 427a4d96a55SJed Brown 428a4d96a55SJed Brown Level: advanced 429a4d96a55SJed Brown 4309a535bafSVaclav Hapla Notes: 4319a535bafSVaclav Hapla If any processor calls this with start = PETSC_DECIDE then all processors must, otherwise the program will hang. 4329a535bafSVaclav Hapla 433db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingSetFromOptions()` 434a4d96a55SJed Brown @*/ 435a4d96a55SJed Brown PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF sf,PetscInt start,ISLocalToGlobalMapping *mapping) 436a4d96a55SJed Brown { 437a4d96a55SJed Brown PetscInt i,maxlocal,nroots,nleaves,*globals,*ltog; 438a4d96a55SJed Brown MPI_Comm comm; 439a4d96a55SJed Brown 440a4d96a55SJed Brown PetscFunctionBegin; 441a4d96a55SJed Brown PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 442a4d96a55SJed Brown PetscValidPointer(mapping,3); 4439566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf,&comm)); 44441f4c31fSVaclav Hapla PetscCall(PetscSFGetGraph(sf,&nroots,&nleaves,NULL,NULL)); 4459a535bafSVaclav Hapla if (start == PETSC_DECIDE) { 4469a535bafSVaclav Hapla start = 0; 4479566063dSJacob Faibussowitsch PetscCallMPI(MPI_Exscan(&nroots,&start,1,MPIU_INT,MPI_SUM,comm)); 44841f4c31fSVaclav Hapla } else PetscCheck(start >= 0,comm, PETSC_ERR_ARG_OUTOFRANGE, "start must be nonnegative or PETSC_DECIDE"); 44941f4c31fSVaclav Hapla PetscCall(PetscSFGetLeafRange(sf, NULL, &maxlocal)); 45041f4c31fSVaclav Hapla ++maxlocal; 4519566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nroots,&globals)); 4529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxlocal,<og)); 453a4d96a55SJed Brown for (i=0; i<nroots; i++) globals[i] = start + i; 454a4d96a55SJed Brown for (i=0; i<maxlocal; i++) ltog[i] = -1; 4559566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf,MPIU_INT,globals,ltog,MPI_REPLACE)); 4569566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf,MPIU_INT,globals,ltog,MPI_REPLACE)); 4579566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm,1,maxlocal,ltog,PETSC_OWN_POINTER,mapping)); 4589566063dSJacob Faibussowitsch PetscCall(PetscFree(globals)); 459a4d96a55SJed Brown PetscFunctionReturn(0); 460a4d96a55SJed Brown } 461b46b645bSBarry Smith 46263fa5c83Sstefano_zampini /*@ 46363fa5c83Sstefano_zampini ISLocalToGlobalMappingSetBlockSize - Sets the blocksize of the mapping 46463fa5c83Sstefano_zampini 46563fa5c83Sstefano_zampini Not collective 46663fa5c83Sstefano_zampini 46763fa5c83Sstefano_zampini Input Parameters: 468a2b725a8SWilliam Gropp + mapping - mapping data structure 469a2b725a8SWilliam Gropp - bs - the blocksize 47063fa5c83Sstefano_zampini 47163fa5c83Sstefano_zampini Level: advanced 47263fa5c83Sstefano_zampini 473db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()` 47463fa5c83Sstefano_zampini @*/ 47563fa5c83Sstefano_zampini PetscErrorCode ISLocalToGlobalMappingSetBlockSize(ISLocalToGlobalMapping mapping,PetscInt bs) 47663fa5c83Sstefano_zampini { 477a59f3c4dSstefano_zampini PetscInt *nid; 478a59f3c4dSstefano_zampini const PetscInt *oid; 479a59f3c4dSstefano_zampini PetscInt i,cn,on,obs,nn; 48063fa5c83Sstefano_zampini 48163fa5c83Sstefano_zampini PetscFunctionBegin; 48263fa5c83Sstefano_zampini PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 48308401ef6SPierre Jolivet PetscCheck(bs >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid block size %" PetscInt_FMT,bs); 48463fa5c83Sstefano_zampini if (bs == mapping->bs) PetscFunctionReturn(0); 48563fa5c83Sstefano_zampini on = mapping->n; 48663fa5c83Sstefano_zampini obs = mapping->bs; 48763fa5c83Sstefano_zampini oid = mapping->indices; 48863fa5c83Sstefano_zampini nn = (on*obs)/bs; 48908401ef6SPierre Jolivet PetscCheck((on*obs)%bs == 0,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Block size %" PetscInt_FMT " is inconsistent with block size %" PetscInt_FMT " and number of block indices %" PetscInt_FMT,bs,obs,on); 490a59f3c4dSstefano_zampini 4919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nn,&nid)); 4929566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(mapping,&oid)); 493a59f3c4dSstefano_zampini for (i=0;i<nn;i++) { 494a59f3c4dSstefano_zampini PetscInt j; 495a59f3c4dSstefano_zampini for (j=0,cn=0;j<bs-1;j++) { 496a59f3c4dSstefano_zampini if (oid[i*bs+j] < 0) { cn++; continue; } 49708401ef6SPierre Jolivet PetscCheck(oid[i*bs+j] == oid[i*bs+j+1]-1,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Block sizes %" PetscInt_FMT " and %" PetscInt_FMT " are incompatible with the block indices: non consecutive indices %" PetscInt_FMT " %" PetscInt_FMT,bs,obs,oid[i*bs+j],oid[i*bs+j+1]); 498a59f3c4dSstefano_zampini } 499a59f3c4dSstefano_zampini if (oid[i*bs+j] < 0) cn++; 5008b7cb0e6Sstefano_zampini if (cn) { 50108401ef6SPierre Jolivet PetscCheck(cn == bs,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Block sizes %" PetscInt_FMT " and %" PetscInt_FMT " are incompatible with the block indices: invalid number of negative entries in block %" PetscInt_FMT,bs,obs,cn); 502a59f3c4dSstefano_zampini nid[i] = -1; 5038b7cb0e6Sstefano_zampini } else { 504a59f3c4dSstefano_zampini nid[i] = oid[i*bs]/bs; 50563fa5c83Sstefano_zampini } 50663fa5c83Sstefano_zampini } 5079566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(mapping,&oid)); 508a59f3c4dSstefano_zampini 50963fa5c83Sstefano_zampini mapping->n = nn; 51063fa5c83Sstefano_zampini mapping->bs = bs; 5119566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->indices)); 51263fa5c83Sstefano_zampini mapping->indices = nid; 513c9345713Sstefano_zampini mapping->globalstart = 0; 514c9345713Sstefano_zampini mapping->globalend = 0; 5151bd0b88eSStefano Zampini 5161bd0b88eSStefano Zampini /* reset the cached information */ 5179566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->info_procs)); 5189566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->info_numprocs)); 5191bd0b88eSStefano Zampini if (mapping->info_indices) { 5201bd0b88eSStefano Zampini PetscInt i; 5211bd0b88eSStefano Zampini 5229566063dSJacob Faibussowitsch PetscCall(PetscFree((mapping->info_indices)[0])); 5231bd0b88eSStefano Zampini for (i=1; i<mapping->info_nproc; i++) { 5249566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->info_indices[i])); 5251bd0b88eSStefano Zampini } 5269566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping->info_indices)); 5271bd0b88eSStefano Zampini } 5281bd0b88eSStefano Zampini mapping->info_cached = PETSC_FALSE; 5291bd0b88eSStefano Zampini 530413f72f0SBarry Smith if (mapping->ops->destroy) { 5319566063dSJacob Faibussowitsch PetscCall((*mapping->ops->destroy)(mapping)); 532413f72f0SBarry Smith } 53363fa5c83Sstefano_zampini PetscFunctionReturn(0); 53463fa5c83Sstefano_zampini } 53563fa5c83Sstefano_zampini 53645b6f7e9SBarry Smith /*@ 53745b6f7e9SBarry Smith ISLocalToGlobalMappingGetBlockSize - Gets the blocksize of the mapping 53845b6f7e9SBarry Smith ordering and a global parallel ordering. 53945b6f7e9SBarry Smith 54045b6f7e9SBarry Smith Not Collective 54145b6f7e9SBarry Smith 54245b6f7e9SBarry Smith Input Parameters: 54345b6f7e9SBarry Smith . mapping - mapping data structure 54445b6f7e9SBarry Smith 54545b6f7e9SBarry Smith Output Parameter: 54645b6f7e9SBarry Smith . bs - the blocksize 54745b6f7e9SBarry Smith 54845b6f7e9SBarry Smith Level: advanced 54945b6f7e9SBarry Smith 550db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()` 55145b6f7e9SBarry Smith @*/ 55245b6f7e9SBarry Smith PetscErrorCode ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping mapping,PetscInt *bs) 55345b6f7e9SBarry Smith { 55445b6f7e9SBarry Smith PetscFunctionBegin; 555cbc1caf0SMatthew G. Knepley PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 55645b6f7e9SBarry Smith *bs = mapping->bs; 55745b6f7e9SBarry Smith PetscFunctionReturn(0); 55845b6f7e9SBarry Smith } 55945b6f7e9SBarry Smith 560ba5bb76aSSatish Balay /*@ 56190f02eecSBarry Smith ISLocalToGlobalMappingCreate - Creates a mapping between a local (0 to n) 56290f02eecSBarry Smith ordering and a global parallel ordering. 5632362add9SBarry Smith 56489d82c54SBarry Smith Not Collective, but communicator may have more than one process 565b9cd556bSLois Curfman McInnes 5662362add9SBarry Smith Input Parameters: 56789d82c54SBarry Smith + comm - MPI communicator 568f0413b6fSBarry Smith . bs - the block size 56928bc9809SBarry Smith . n - the number of local elements divided by the block size, or equivalently the number of block indices 57028bc9809SBarry Smith . indices - the global index for each local element, these do not need to be in increasing order (sorted), these values should not be scaled (i.e. multiplied) by the blocksize bs 571d5ad8652SBarry Smith - mode - see PetscCopyMode 5722362add9SBarry Smith 573a997ad1aSLois Curfman McInnes Output Parameter: 57490f02eecSBarry Smith . mapping - new mapping data structure 5752362add9SBarry Smith 57695452b02SPatrick Sanan Notes: 57795452b02SPatrick Sanan There is one integer value in indices per block and it represents the actual indices bs*idx + j, where j=0,..,bs-1 578413f72f0SBarry Smith 5799a7b7924SJed Brown For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used; 580413f72f0SBarry Smith this uses more memory but is faster; this approach is not scalable for extremely large mappings. For large problems ISLOCALTOGLOBALMAPPINGHASH is used, this is scalable. 581413f72f0SBarry Smith Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used. 582413f72f0SBarry Smith 583a997ad1aSLois Curfman McInnes Level: advanced 584a997ad1aSLois Curfman McInnes 585db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingSetFromOptions()`, `ISLOCALTOGLOBALMAPPINGBASIC`, `ISLOCALTOGLOBALMAPPINGHASH` 586db781477SPatrick Sanan `ISLocalToGlobalMappingSetType()`, `ISLocalToGlobalMappingType` 5872362add9SBarry Smith @*/ 58860c7cefcSBarry Smith PetscErrorCode ISLocalToGlobalMappingCreate(MPI_Comm comm,PetscInt bs,PetscInt n,const PetscInt indices[],PetscCopyMode mode,ISLocalToGlobalMapping *mapping) 5892362add9SBarry Smith { 59032dcc486SBarry Smith PetscInt *in; 591b46b645bSBarry Smith 592b46b645bSBarry Smith PetscFunctionBegin; 593064a246eSJacob Faibussowitsch if (n) PetscValidIntPointer(indices,4); 594064a246eSJacob Faibussowitsch PetscValidPointer(mapping,6); 595b46b645bSBarry Smith 5960298fd71SBarry Smith *mapping = NULL; 5979566063dSJacob Faibussowitsch PetscCall(ISInitializePackage()); 5982362add9SBarry Smith 5999566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(*mapping,IS_LTOGM_CLASSID,"ISLocalToGlobalMapping","Local to global mapping","IS",comm,ISLocalToGlobalMappingDestroy,ISLocalToGlobalMappingView)); 600d4bb536fSBarry Smith (*mapping)->n = n; 601f0413b6fSBarry Smith (*mapping)->bs = bs; 602d5ad8652SBarry Smith if (mode == PETSC_COPY_VALUES) { 6039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&in)); 6049566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(in,indices,n)); 605d5ad8652SBarry Smith (*mapping)->indices = in; 60671910c26SVaclav Hapla (*mapping)->dealloc_indices = PETSC_TRUE; 6079566063dSJacob Faibussowitsch PetscCall(PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt))); 6086389a1a1SBarry Smith } else if (mode == PETSC_OWN_POINTER) { 6096389a1a1SBarry Smith (*mapping)->indices = (PetscInt*)indices; 61071910c26SVaclav Hapla (*mapping)->dealloc_indices = PETSC_TRUE; 6119566063dSJacob Faibussowitsch PetscCall(PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt))); 61271910c26SVaclav Hapla } else if (mode == PETSC_USE_POINTER) { 61371910c26SVaclav Hapla (*mapping)->indices = (PetscInt*)indices; 6146389a1a1SBarry Smith } 61598921bdaSJacob Faibussowitsch else SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid mode %d", mode); 6163a40ed3dSBarry Smith PetscFunctionReturn(0); 6172362add9SBarry Smith } 6182362add9SBarry Smith 619413f72f0SBarry Smith PetscFunctionList ISLocalToGlobalMappingList = NULL; 620413f72f0SBarry Smith 62190f02eecSBarry Smith /*@ 6227e99dc12SLawrence Mitchell ISLocalToGlobalMappingSetFromOptions - Set mapping options from the options database. 6237e99dc12SLawrence Mitchell 6247e99dc12SLawrence Mitchell Not collective 6257e99dc12SLawrence Mitchell 6267e99dc12SLawrence Mitchell Input Parameters: 6277e99dc12SLawrence Mitchell . mapping - mapping data structure 6287e99dc12SLawrence Mitchell 6297e99dc12SLawrence Mitchell Level: advanced 6307e99dc12SLawrence Mitchell 6317e99dc12SLawrence Mitchell @*/ 6327e99dc12SLawrence Mitchell PetscErrorCode ISLocalToGlobalMappingSetFromOptions(ISLocalToGlobalMapping mapping) 6337e99dc12SLawrence Mitchell { 634413f72f0SBarry Smith char type[256]; 635413f72f0SBarry Smith ISLocalToGlobalMappingType defaulttype = "Not set"; 6367e99dc12SLawrence Mitchell PetscBool flg; 6377e99dc12SLawrence Mitchell 6387e99dc12SLawrence Mitchell PetscFunctionBegin; 6397e99dc12SLawrence Mitchell PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 6409566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRegisterAll()); 641d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)mapping); 6429566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-islocaltoglobalmapping_type","ISLocalToGlobalMapping method","ISLocalToGlobalMappingSetType",ISLocalToGlobalMappingList,(char*)(((PetscObject)mapping)->type_name) ? ((PetscObject)mapping)->type_name : defaulttype,type,256,&flg)); 643413f72f0SBarry Smith if (flg) { 6449566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingSetType(mapping,type)); 645413f72f0SBarry Smith } 646d0609cedSBarry Smith PetscOptionsEnd(); 6477e99dc12SLawrence Mitchell PetscFunctionReturn(0); 6487e99dc12SLawrence Mitchell } 6497e99dc12SLawrence Mitchell 6507e99dc12SLawrence Mitchell /*@ 65190f02eecSBarry Smith ISLocalToGlobalMappingDestroy - Destroys a mapping between a local (0 to n) 65290f02eecSBarry Smith ordering and a global parallel ordering. 65390f02eecSBarry Smith 6540f5bd95cSBarry Smith Note Collective 655b9cd556bSLois Curfman McInnes 65690f02eecSBarry Smith Input Parameters: 65790f02eecSBarry Smith . mapping - mapping data structure 65890f02eecSBarry Smith 659a997ad1aSLois Curfman McInnes Level: advanced 660a997ad1aSLois Curfman McInnes 661db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()` 66290f02eecSBarry Smith @*/ 6636bf464f9SBarry Smith PetscErrorCode ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping *mapping) 66490f02eecSBarry Smith { 6653a40ed3dSBarry Smith PetscFunctionBegin; 6666bf464f9SBarry Smith if (!*mapping) PetscFunctionReturn(0); 6676bf464f9SBarry Smith PetscValidHeaderSpecific((*mapping),IS_LTOGM_CLASSID,1); 6684c8fdceaSLisandro Dalcin if (--((PetscObject)(*mapping))->refct > 0) {*mapping = NULL;PetscFunctionReturn(0);} 66971910c26SVaclav Hapla if ((*mapping)->dealloc_indices) { 6709566063dSJacob Faibussowitsch PetscCall(PetscFree((*mapping)->indices)); 67171910c26SVaclav Hapla } 6729566063dSJacob Faibussowitsch PetscCall(PetscFree((*mapping)->info_procs)); 6739566063dSJacob Faibussowitsch PetscCall(PetscFree((*mapping)->info_numprocs)); 674268a049cSStefano Zampini if ((*mapping)->info_indices) { 675268a049cSStefano Zampini PetscInt i; 676268a049cSStefano Zampini 6779566063dSJacob Faibussowitsch PetscCall(PetscFree(((*mapping)->info_indices)[0])); 678268a049cSStefano Zampini for (i=1; i<(*mapping)->info_nproc; i++) { 6799566063dSJacob Faibussowitsch PetscCall(PetscFree(((*mapping)->info_indices)[i])); 680268a049cSStefano Zampini } 6819566063dSJacob Faibussowitsch PetscCall(PetscFree((*mapping)->info_indices)); 682268a049cSStefano Zampini } 6831bd0b88eSStefano Zampini if ((*mapping)->info_nodei) { 6849566063dSJacob Faibussowitsch PetscCall(PetscFree(((*mapping)->info_nodei)[0])); 6851bd0b88eSStefano Zampini } 6869566063dSJacob Faibussowitsch PetscCall(PetscFree2((*mapping)->info_nodec,(*mapping)->info_nodei)); 687413f72f0SBarry Smith if ((*mapping)->ops->destroy) { 6889566063dSJacob Faibussowitsch PetscCall((*(*mapping)->ops->destroy)(*mapping)); 689413f72f0SBarry Smith } 6909566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(mapping)); 6914c8fdceaSLisandro Dalcin *mapping = NULL; 6923a40ed3dSBarry Smith PetscFunctionReturn(0); 69390f02eecSBarry Smith } 69490f02eecSBarry Smith 69590f02eecSBarry Smith /*@ 6963acfe500SLois Curfman McInnes ISLocalToGlobalMappingApplyIS - Creates from an IS in the local numbering 6973acfe500SLois Curfman McInnes a new index set using the global numbering defined in an ISLocalToGlobalMapping 6983acfe500SLois Curfman McInnes context. 69990f02eecSBarry Smith 7004cb36875SStefano Zampini Collective on is 701b9cd556bSLois Curfman McInnes 70290f02eecSBarry Smith Input Parameters: 703b9cd556bSLois Curfman McInnes + mapping - mapping between local and global numbering 704b9cd556bSLois Curfman McInnes - is - index set in local numbering 70590f02eecSBarry Smith 70690f02eecSBarry Smith Output Parameters: 70790f02eecSBarry Smith . newis - index set in global numbering 70890f02eecSBarry Smith 7094cb36875SStefano Zampini Notes: 7104cb36875SStefano Zampini The output IS will have the same communicator of the input IS. 7114cb36875SStefano Zampini 712a997ad1aSLois Curfman McInnes Level: advanced 713a997ad1aSLois Curfman McInnes 714db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingCreate()`, 715db781477SPatrick Sanan `ISLocalToGlobalMappingDestroy()`, `ISGlobalToLocalMappingApply()` 71690f02eecSBarry Smith @*/ 7177087cfbeSBarry Smith PetscErrorCode ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping,IS is,IS *newis) 71890f02eecSBarry Smith { 719e24637baSBarry Smith PetscInt n,*idxout; 7205d0c19d7SBarry Smith const PetscInt *idxin; 7213a40ed3dSBarry Smith 7223a40ed3dSBarry Smith PetscFunctionBegin; 7230700a824SBarry Smith PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 7240700a824SBarry Smith PetscValidHeaderSpecific(is,IS_CLASSID,2); 7254482741eSBarry Smith PetscValidPointer(newis,3); 72690f02eecSBarry Smith 7279566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is,&n)); 7289566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is,&idxin)); 7299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&idxout)); 7309566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(mapping,n,idxin,idxout)); 7319566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is,&idxin)); 7329566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)is),n,idxout,PETSC_OWN_POINTER,newis)); 7333a40ed3dSBarry Smith PetscFunctionReturn(0); 73490f02eecSBarry Smith } 73590f02eecSBarry Smith 736b89cb25eSSatish Balay /*@ 7373acfe500SLois Curfman McInnes ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering 7383acfe500SLois Curfman McInnes and converts them to the global numbering. 73990f02eecSBarry Smith 740b9cd556bSLois Curfman McInnes Not collective 741b9cd556bSLois Curfman McInnes 742bb25748dSBarry Smith Input Parameters: 743b9cd556bSLois Curfman McInnes + mapping - the local to global mapping context 744bb25748dSBarry Smith . N - number of integers 745b9cd556bSLois Curfman McInnes - in - input indices in local numbering 746bb25748dSBarry Smith 747bb25748dSBarry Smith Output Parameter: 748bb25748dSBarry Smith . out - indices in global numbering 749bb25748dSBarry Smith 750b9cd556bSLois Curfman McInnes Notes: 751b9cd556bSLois Curfman McInnes The in and out array parameters may be identical. 752d4bb536fSBarry Smith 753a997ad1aSLois Curfman McInnes Level: advanced 754a997ad1aSLois Curfman McInnes 755*c2e3fba1SPatrick Sanan .seealso: `ISLocalToGlobalMappingApplyBlock()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingDestroy()`, 756*c2e3fba1SPatrick Sanan `ISLocalToGlobalMappingApplyIS()`, `AOCreateBasic()`, `AOApplicationToPetsc()`, 757db781477SPatrick Sanan `AOPetscToApplication()`, `ISGlobalToLocalMappingApply()` 758bb25748dSBarry Smith 759afcb2eb5SJed Brown @*/ 760afcb2eb5SJed Brown PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[]) 761afcb2eb5SJed Brown { 762cbc1caf0SMatthew G. Knepley PetscInt i,bs,Nmax; 76345b6f7e9SBarry Smith 76445b6f7e9SBarry Smith PetscFunctionBegin; 765cbc1caf0SMatthew G. Knepley PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 766cbc1caf0SMatthew G. Knepley bs = mapping->bs; 767cbc1caf0SMatthew G. Knepley Nmax = bs*mapping->n; 76845b6f7e9SBarry Smith if (bs == 1) { 769cbc1caf0SMatthew G. Knepley const PetscInt *idx = mapping->indices; 77045b6f7e9SBarry Smith for (i=0; i<N; i++) { 77145b6f7e9SBarry Smith if (in[i] < 0) { 77245b6f7e9SBarry Smith out[i] = in[i]; 77345b6f7e9SBarry Smith continue; 77445b6f7e9SBarry Smith } 77508401ef6SPierre Jolivet PetscCheck(in[i] < Nmax,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local index %" PetscInt_FMT " too large %" PetscInt_FMT " (max) at %" PetscInt_FMT,in[i],Nmax-1,i); 77645b6f7e9SBarry Smith out[i] = idx[in[i]]; 77745b6f7e9SBarry Smith } 77845b6f7e9SBarry Smith } else { 779cbc1caf0SMatthew G. Knepley const PetscInt *idx = mapping->indices; 78045b6f7e9SBarry Smith for (i=0; i<N; i++) { 78145b6f7e9SBarry Smith if (in[i] < 0) { 78245b6f7e9SBarry Smith out[i] = in[i]; 78345b6f7e9SBarry Smith continue; 78445b6f7e9SBarry Smith } 78508401ef6SPierre 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); 78645b6f7e9SBarry Smith out[i] = idx[in[i]/bs]*bs + (in[i] % bs); 78745b6f7e9SBarry Smith } 78845b6f7e9SBarry Smith } 78945b6f7e9SBarry Smith PetscFunctionReturn(0); 79045b6f7e9SBarry Smith } 79145b6f7e9SBarry Smith 79245b6f7e9SBarry Smith /*@ 7936006e8d2SBarry Smith ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local block numbering and converts them to the global block numbering 79445b6f7e9SBarry Smith 79545b6f7e9SBarry Smith Not collective 79645b6f7e9SBarry Smith 79745b6f7e9SBarry Smith Input Parameters: 79845b6f7e9SBarry Smith + mapping - the local to global mapping context 79945b6f7e9SBarry Smith . N - number of integers 8006006e8d2SBarry Smith - in - input indices in local block numbering 80145b6f7e9SBarry Smith 80245b6f7e9SBarry Smith Output Parameter: 8036006e8d2SBarry Smith . out - indices in global block numbering 80445b6f7e9SBarry Smith 80545b6f7e9SBarry Smith Notes: 80645b6f7e9SBarry Smith The in and out array parameters may be identical. 80745b6f7e9SBarry Smith 8086006e8d2SBarry Smith Example: 8096006e8d2SBarry 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 8106006e8d2SBarry Smith (the first block) would produce 0 and the mapping applied to 1 (the second block) would produce 3. 8116006e8d2SBarry Smith 81245b6f7e9SBarry Smith Level: advanced 81345b6f7e9SBarry Smith 814*c2e3fba1SPatrick Sanan .seealso: `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingDestroy()`, 815*c2e3fba1SPatrick Sanan `ISLocalToGlobalMappingApplyIS()`, `AOCreateBasic()`, `AOApplicationToPetsc()`, 816db781477SPatrick Sanan `AOPetscToApplication()`, `ISGlobalToLocalMappingApply()` 81745b6f7e9SBarry Smith 81845b6f7e9SBarry Smith @*/ 81945b6f7e9SBarry Smith PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[]) 82045b6f7e9SBarry Smith { 8218a1f772fSStefano Zampini PetscInt i,Nmax; 8228a1f772fSStefano Zampini const PetscInt *idx; 823d4bb536fSBarry Smith 824a0d79125SStefano Zampini PetscFunctionBegin; 825a0d79125SStefano Zampini PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 8268a1f772fSStefano Zampini Nmax = mapping->n; 8278a1f772fSStefano Zampini idx = mapping->indices; 828afcb2eb5SJed Brown for (i=0; i<N; i++) { 829afcb2eb5SJed Brown if (in[i] < 0) { 830afcb2eb5SJed Brown out[i] = in[i]; 831afcb2eb5SJed Brown continue; 832afcb2eb5SJed Brown } 83308401ef6SPierre 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); 834afcb2eb5SJed Brown out[i] = idx[in[i]]; 835afcb2eb5SJed Brown } 836afcb2eb5SJed Brown PetscFunctionReturn(0); 837afcb2eb5SJed Brown } 838d4bb536fSBarry Smith 8397e99dc12SLawrence Mitchell /*@ 840a997ad1aSLois Curfman McInnes ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers 841a997ad1aSLois Curfman McInnes specified with a global numbering. 842d4bb536fSBarry Smith 843b9cd556bSLois Curfman McInnes Not collective 844b9cd556bSLois Curfman McInnes 845d4bb536fSBarry Smith Input Parameters: 846b9cd556bSLois Curfman McInnes + mapping - mapping between local and global numbering 8470040bde1SJunchao Zhang . type - IS_GTOLM_MASK - maps global indices with no local value to -1 in the output list (i.e., mask them) 848d4bb536fSBarry Smith IS_GTOLM_DROP - drops the indices with no local value from the output list 849d4bb536fSBarry Smith . n - number of global indices to map 850b9cd556bSLois Curfman McInnes - idx - global indices to map 851d4bb536fSBarry Smith 852d4bb536fSBarry Smith Output Parameters: 853b9cd556bSLois Curfman McInnes + nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n) 854b9cd556bSLois Curfman McInnes - idxout - local index of each global index, one must pass in an array long enough 855e182c471SBarry Smith to hold all the indices. You can call ISGlobalToLocalMappingApply() with 8560298fd71SBarry Smith idxout == NULL to determine the required length (returned in nout) 857e182c471SBarry Smith and then allocate the required space and call ISGlobalToLocalMappingApply() 858e182c471SBarry Smith a second time to set the values. 859d4bb536fSBarry Smith 860b9cd556bSLois Curfman McInnes Notes: 8610298fd71SBarry Smith Either nout or idxout may be NULL. idx and idxout may be identical. 862d4bb536fSBarry Smith 8639a7b7924SJed Brown For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used; 864413f72f0SBarry 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. 865413f72f0SBarry Smith Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used. 8660f5bd95cSBarry Smith 867a997ad1aSLois Curfman McInnes Level: advanced 868a997ad1aSLois Curfman McInnes 86932fd6b96SBarry Smith Developer Note: The manual page states that idx and idxout may be identical but the calling 87032fd6b96SBarry Smith sequence declares idx as const so it cannot be the same as idxout. 87132fd6b96SBarry Smith 872db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingApply()`, `ISGlobalToLocalMappingApplyBlock()`, `ISLocalToGlobalMappingCreate()`, 873db781477SPatrick Sanan `ISLocalToGlobalMappingDestroy()` 874d4bb536fSBarry Smith @*/ 875413f72f0SBarry Smith PetscErrorCode ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[]) 876d4bb536fSBarry Smith { 8779d90f715SBarry Smith PetscFunctionBegin; 8789d90f715SBarry Smith PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 879413f72f0SBarry Smith if (!mapping->data) { 8809566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingSetUp(mapping)); 8819d90f715SBarry Smith } 8829566063dSJacob Faibussowitsch PetscCall((*mapping->ops->globaltolocalmappingapply)(mapping,type,n,idx,nout,idxout)); 8839d90f715SBarry Smith PetscFunctionReturn(0); 8849d90f715SBarry Smith } 8859d90f715SBarry Smith 886d4fe737eSStefano Zampini /*@ 887d4fe737eSStefano Zampini ISGlobalToLocalMappingApplyIS - Creates from an IS in the global numbering 888d4fe737eSStefano Zampini a new index set using the local numbering defined in an ISLocalToGlobalMapping 889d4fe737eSStefano Zampini context. 890d4fe737eSStefano Zampini 891d4fe737eSStefano Zampini Not collective 892d4fe737eSStefano Zampini 893d4fe737eSStefano Zampini Input Parameters: 894d4fe737eSStefano Zampini + mapping - mapping between local and global numbering 8950040bde1SJunchao Zhang . type - IS_GTOLM_MASK - maps global indices with no local value to -1 in the output list (i.e., mask them) 8962785b321SStefano Zampini IS_GTOLM_DROP - drops the indices with no local value from the output list 897d4fe737eSStefano Zampini - is - index set in global numbering 898d4fe737eSStefano Zampini 899d4fe737eSStefano Zampini Output Parameters: 900d4fe737eSStefano Zampini . newis - index set in local numbering 901d4fe737eSStefano Zampini 9024cb36875SStefano Zampini Notes: 9034cb36875SStefano Zampini The output IS will be sequential, as it encodes a purely local operation 9044cb36875SStefano Zampini 905d4fe737eSStefano Zampini Level: advanced 906d4fe737eSStefano Zampini 907db781477SPatrick Sanan .seealso: `ISGlobalToLocalMappingApply()`, `ISLocalToGlobalMappingCreate()`, 908db781477SPatrick Sanan `ISLocalToGlobalMappingDestroy()` 909d4fe737eSStefano Zampini @*/ 910413f72f0SBarry Smith PetscErrorCode ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,IS is,IS *newis) 911d4fe737eSStefano Zampini { 912d4fe737eSStefano Zampini PetscInt n,nout,*idxout; 913d4fe737eSStefano Zampini const PetscInt *idxin; 914d4fe737eSStefano Zampini 915d4fe737eSStefano Zampini PetscFunctionBegin; 916d4fe737eSStefano Zampini PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 917d4fe737eSStefano Zampini PetscValidHeaderSpecific(is,IS_CLASSID,3); 918d4fe737eSStefano Zampini PetscValidPointer(newis,4); 919d4fe737eSStefano Zampini 9209566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is,&n)); 9219566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is,&idxin)); 922d4fe737eSStefano Zampini if (type == IS_GTOLM_MASK) { 9239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&idxout)); 924d4fe737eSStefano Zampini } else { 9259566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,NULL)); 9269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nout,&idxout)); 927d4fe737eSStefano Zampini } 9289566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,idxout)); 9299566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is,&idxin)); 9309566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,nout,idxout,PETSC_OWN_POINTER,newis)); 931d4fe737eSStefano Zampini PetscFunctionReturn(0); 932d4fe737eSStefano Zampini } 933d4fe737eSStefano Zampini 9349d90f715SBarry Smith /*@ 9359d90f715SBarry Smith ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers 9369d90f715SBarry Smith specified with a block global numbering. 9379d90f715SBarry Smith 9389d90f715SBarry Smith Not collective 9399d90f715SBarry Smith 9409d90f715SBarry Smith Input Parameters: 9419d90f715SBarry Smith + mapping - mapping between local and global numbering 9420040bde1SJunchao Zhang . type - IS_GTOLM_MASK - maps global indices with no local value to -1 in the output list (i.e., mask them) 9439d90f715SBarry Smith IS_GTOLM_DROP - drops the indices with no local value from the output list 9449d90f715SBarry Smith . n - number of global indices to map 9459d90f715SBarry Smith - idx - global indices to map 9469d90f715SBarry Smith 9479d90f715SBarry Smith Output Parameters: 9489d90f715SBarry Smith + nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n) 9499d90f715SBarry Smith - idxout - local index of each global index, one must pass in an array long enough 9509d90f715SBarry Smith to hold all the indices. You can call ISGlobalToLocalMappingApplyBlock() with 9519d90f715SBarry Smith idxout == NULL to determine the required length (returned in nout) 9529d90f715SBarry Smith and then allocate the required space and call ISGlobalToLocalMappingApplyBlock() 9539d90f715SBarry Smith a second time to set the values. 9549d90f715SBarry Smith 9559d90f715SBarry Smith Notes: 9569d90f715SBarry Smith Either nout or idxout may be NULL. idx and idxout may be identical. 9579d90f715SBarry Smith 9589a7b7924SJed Brown For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used; 959413f72f0SBarry 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. 960413f72f0SBarry Smith Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used. 9619d90f715SBarry Smith 9629d90f715SBarry Smith Level: advanced 9639d90f715SBarry Smith 9649d90f715SBarry Smith Developer Note: The manual page states that idx and idxout may be identical but the calling 9659d90f715SBarry Smith sequence declares idx as const so it cannot be the same as idxout. 9669d90f715SBarry Smith 967db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingApply()`, `ISGlobalToLocalMappingApply()`, `ISLocalToGlobalMappingCreate()`, 968db781477SPatrick Sanan `ISLocalToGlobalMappingDestroy()` 9699d90f715SBarry Smith @*/ 970413f72f0SBarry Smith PetscErrorCode ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type, 9719d90f715SBarry Smith PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[]) 9729d90f715SBarry Smith { 9733a40ed3dSBarry Smith PetscFunctionBegin; 9740700a824SBarry Smith PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 975413f72f0SBarry Smith if (!mapping->data) { 9769566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingSetUp(mapping)); 977d4bb536fSBarry Smith } 9789566063dSJacob Faibussowitsch PetscCall((*mapping->ops->globaltolocalmappingapplyblock)(mapping,type,n,idx,nout,idxout)); 9793a40ed3dSBarry Smith PetscFunctionReturn(0); 980d4bb536fSBarry Smith } 98190f02eecSBarry Smith 98289d82c54SBarry Smith /*@C 9836a818285SBarry Smith ISLocalToGlobalMappingGetBlockInfo - Gets the neighbor information for each processor and 98489d82c54SBarry Smith each index shared by more than one processor 98589d82c54SBarry Smith 98689d82c54SBarry Smith Collective on ISLocalToGlobalMapping 98789d82c54SBarry Smith 988f899ff85SJose E. Roman Input Parameter: 98989d82c54SBarry Smith . mapping - the mapping from local to global indexing 99089d82c54SBarry Smith 991d8d19677SJose E. Roman Output Parameters: 99289d82c54SBarry Smith + nproc - number of processors that are connected to this one 99389d82c54SBarry Smith . proc - neighboring processors 99407b52d57SBarry Smith . numproc - number of indices for each subdomain (processor) 9953463a7baSJed Brown - indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering) 99689d82c54SBarry Smith 99789d82c54SBarry Smith Level: advanced 99889d82c54SBarry Smith 9992cfcea29SBarry Smith Fortran Usage: 10002cfcea29SBarry Smith $ ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by 10012cfcea29SBarry Smith $ ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc], 10022cfcea29SBarry Smith PetscInt indices[nproc][numprocmax],ierr) 10032cfcea29SBarry Smith There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and 10042cfcea29SBarry Smith indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough. 10052cfcea29SBarry Smith 1006db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 1007db781477SPatrick Sanan `ISLocalToGlobalMappingRestoreInfo()` 100889d82c54SBarry Smith @*/ 10096a818285SBarry Smith PetscErrorCode ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[]) 101089d82c54SBarry Smith { 1011268a049cSStefano Zampini PetscFunctionBegin; 1012268a049cSStefano Zampini PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 1013268a049cSStefano Zampini if (mapping->info_cached) { 1014268a049cSStefano Zampini *nproc = mapping->info_nproc; 1015268a049cSStefano Zampini *procs = mapping->info_procs; 1016268a049cSStefano Zampini *numprocs = mapping->info_numprocs; 1017268a049cSStefano Zampini *indices = mapping->info_indices; 1018268a049cSStefano Zampini } else { 10199566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockInfo_Private(mapping,nproc,procs,numprocs,indices)); 1020268a049cSStefano Zampini } 1021268a049cSStefano Zampini PetscFunctionReturn(0); 1022268a049cSStefano Zampini } 1023268a049cSStefano Zampini 1024268a049cSStefano Zampini static PetscErrorCode ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[]) 1025268a049cSStefano Zampini { 102697f1f81fSBarry Smith PetscMPIInt size,rank,tag1,tag2,tag3,*len,*source,imdex; 102732dcc486SBarry Smith PetscInt i,n = mapping->n,Ng,ng,max = 0,*lindices = mapping->indices; 102832dcc486SBarry Smith PetscInt *nprocs,*owner,nsends,*sends,j,*starts,nmax,nrecvs,*recvs,proc; 102997f1f81fSBarry Smith PetscInt cnt,scale,*ownedsenders,*nownedsenders,rstart,nowned; 103032dcc486SBarry Smith PetscInt node,nownedm,nt,*sends2,nsends2,*starts2,*lens2,*dest,nrecvs2,*starts3,*recvs2,k,*bprocs,*tmp; 103132dcc486SBarry Smith PetscInt first_procs,first_numprocs,*first_indices; 103289d82c54SBarry Smith MPI_Request *recv_waits,*send_waits; 103330dcb7c9SBarry Smith MPI_Status recv_status,*send_status,*recv_statuses; 1034ce94432eSBarry Smith MPI_Comm comm; 1035ace3abfcSBarry Smith PetscBool debug = PETSC_FALSE; 103689d82c54SBarry Smith 103789d82c54SBarry Smith PetscFunctionBegin; 10389566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)mapping,&comm)); 10399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 10409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 104124cf384cSBarry Smith if (size == 1) { 104224cf384cSBarry Smith *nproc = 0; 10430298fd71SBarry Smith *procs = NULL; 10449566063dSJacob Faibussowitsch PetscCall(PetscNew(numprocs)); 10451e2105dcSBarry Smith (*numprocs)[0] = 0; 10469566063dSJacob Faibussowitsch PetscCall(PetscNew(indices)); 10470298fd71SBarry Smith (*indices)[0] = NULL; 1048268a049cSStefano Zampini /* save info for reuse */ 1049268a049cSStefano Zampini mapping->info_nproc = *nproc; 1050268a049cSStefano Zampini mapping->info_procs = *procs; 1051268a049cSStefano Zampini mapping->info_numprocs = *numprocs; 1052268a049cSStefano Zampini mapping->info_indices = *indices; 1053268a049cSStefano Zampini mapping->info_cached = PETSC_TRUE; 105424cf384cSBarry Smith PetscFunctionReturn(0); 105524cf384cSBarry Smith } 105624cf384cSBarry Smith 10579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)mapping)->options,NULL,"-islocaltoglobalmappinggetinfo_debug",&debug,NULL)); 105807b52d57SBarry Smith 10593677ff5aSBarry Smith /* 10606a818285SBarry Smith Notes on ISLocalToGlobalMappingGetBlockInfo 10613677ff5aSBarry Smith 10623677ff5aSBarry Smith globally owned node - the nodes that have been assigned to this processor in global 10633677ff5aSBarry Smith numbering, just for this routine. 10643677ff5aSBarry Smith 10653677ff5aSBarry Smith nontrivial globally owned node - node assigned to this processor that is on a subdomain 10663677ff5aSBarry Smith boundary (i.e. is has more than one local owner) 10673677ff5aSBarry Smith 10683677ff5aSBarry Smith locally owned node - node that exists on this processors subdomain 10693677ff5aSBarry Smith 10703677ff5aSBarry Smith nontrivial locally owned node - node that is not in the interior (i.e. has more than one 10713677ff5aSBarry Smith local subdomain 10723677ff5aSBarry Smith */ 10739566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)mapping,&tag1)); 10749566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)mapping,&tag2)); 10759566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)mapping,&tag3)); 107689d82c54SBarry Smith 107789d82c54SBarry Smith for (i=0; i<n; i++) { 107889d82c54SBarry Smith if (lindices[i] > max) max = lindices[i]; 107989d82c54SBarry Smith } 10801c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&max,&Ng,1,MPIU_INT,MPI_MAX,comm)); 108178058e43SBarry Smith Ng++; 10829566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 10839566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 1084bc8ff85bSBarry Smith scale = Ng/size + 1; 1085a2e34c3dSBarry Smith ng = scale; if (rank == size-1) ng = Ng - scale*(size-1); ng = PetscMax(1,ng); 1086caba0dd0SBarry Smith rstart = scale*rank; 108789d82c54SBarry Smith 108889d82c54SBarry Smith /* determine ownership ranges of global indices */ 10899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2*size,&nprocs)); 10909566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(nprocs,2*size)); 109189d82c54SBarry Smith 109289d82c54SBarry Smith /* determine owners of each local node */ 10939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&owner)); 109489d82c54SBarry Smith for (i=0; i<n; i++) { 10953677ff5aSBarry Smith proc = lindices[i]/scale; /* processor that globally owns this index */ 109627c402fcSBarry Smith nprocs[2*proc+1] = 1; /* processor globally owns at least one of ours */ 10973677ff5aSBarry Smith owner[i] = proc; 109827c402fcSBarry Smith nprocs[2*proc]++; /* count of how many that processor globally owns of ours */ 109989d82c54SBarry Smith } 110027c402fcSBarry Smith nsends = 0; for (i=0; i<size; i++) nsends += nprocs[2*i+1]; 11019566063dSJacob Faibussowitsch PetscCall(PetscInfo(mapping,"Number of global owners for my local data %" PetscInt_FMT "\n",nsends)); 110289d82c54SBarry Smith 110389d82c54SBarry Smith /* inform other processors of number of messages and max length*/ 11049566063dSJacob Faibussowitsch PetscCall(PetscMaxSum(comm,nprocs,&nmax,&nrecvs)); 11059566063dSJacob Faibussowitsch PetscCall(PetscInfo(mapping,"Number of local owners for my global data %" PetscInt_FMT "\n",nrecvs)); 110689d82c54SBarry Smith 110789d82c54SBarry Smith /* post receives for owned rows */ 11089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((2*nrecvs+1)*(nmax+1),&recvs)); 11099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs+1,&recv_waits)); 111089d82c54SBarry Smith for (i=0; i<nrecvs; i++) { 11119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Irecv(recvs+2*nmax*i,2*nmax,MPIU_INT,MPI_ANY_SOURCE,tag1,comm,recv_waits+i)); 111289d82c54SBarry Smith } 111389d82c54SBarry Smith 111489d82c54SBarry Smith /* pack messages containing lists of local nodes to owners */ 11159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2*n+1,&sends)); 11169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size+1,&starts)); 111789d82c54SBarry Smith starts[0] = 0; 1118f6e5521dSKarl Rupp for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2]; 111989d82c54SBarry Smith for (i=0; i<n; i++) { 112089d82c54SBarry Smith sends[starts[owner[i]]++] = lindices[i]; 112130dcb7c9SBarry Smith sends[starts[owner[i]]++] = i; 112289d82c54SBarry Smith } 11239566063dSJacob Faibussowitsch PetscCall(PetscFree(owner)); 112489d82c54SBarry Smith starts[0] = 0; 1125f6e5521dSKarl Rupp for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2]; 112689d82c54SBarry Smith 112789d82c54SBarry Smith /* send the messages */ 11289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends+1,&send_waits)); 11299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends+1,&dest)); 113089d82c54SBarry Smith cnt = 0; 113189d82c54SBarry Smith for (i=0; i<size; i++) { 113227c402fcSBarry Smith if (nprocs[2*i]) { 11339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Isend(sends+starts[i],2*nprocs[2*i],MPIU_INT,i,tag1,comm,send_waits+cnt)); 113430dcb7c9SBarry Smith dest[cnt] = i; 113589d82c54SBarry Smith cnt++; 113689d82c54SBarry Smith } 113789d82c54SBarry Smith } 11389566063dSJacob Faibussowitsch PetscCall(PetscFree(starts)); 113989d82c54SBarry Smith 114089d82c54SBarry Smith /* wait on receives */ 11419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs+1,&source)); 11429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs+1,&len)); 114389d82c54SBarry Smith cnt = nrecvs; 11449566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(ng+1,&nownedsenders)); 114589d82c54SBarry Smith while (cnt) { 11469566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status)); 114789d82c54SBarry Smith /* unpack receives into our local space */ 11489566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&recv_status,MPIU_INT,&len[imdex])); 114989d82c54SBarry Smith source[imdex] = recv_status.MPI_SOURCE; 115030dcb7c9SBarry Smith len[imdex] = len[imdex]/2; 1151caba0dd0SBarry Smith /* count how many local owners for each of my global owned indices */ 115230dcb7c9SBarry Smith for (i=0; i<len[imdex]; i++) nownedsenders[recvs[2*imdex*nmax+2*i]-rstart]++; 115389d82c54SBarry Smith cnt--; 115489d82c54SBarry Smith } 11559566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_waits)); 115689d82c54SBarry Smith 115730dcb7c9SBarry Smith /* count how many globally owned indices are on an edge multiplied by how many processors own them. */ 1158bc8ff85bSBarry Smith nowned = 0; 1159bc8ff85bSBarry Smith nownedm = 0; 1160bc8ff85bSBarry Smith for (i=0; i<ng; i++) { 1161bc8ff85bSBarry Smith if (nownedsenders[i] > 1) {nownedm += nownedsenders[i]; nowned++;} 1162bc8ff85bSBarry Smith } 1163bc8ff85bSBarry Smith 1164bc8ff85bSBarry Smith /* create single array to contain rank of all local owners of each globally owned index */ 11659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nownedm+1,&ownedsenders)); 11669566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ng+1,&starts)); 1167bc8ff85bSBarry Smith starts[0] = 0; 1168bc8ff85bSBarry Smith for (i=1; i<ng; i++) { 1169bc8ff85bSBarry Smith if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1]; 1170bc8ff85bSBarry Smith else starts[i] = starts[i-1]; 1171bc8ff85bSBarry Smith } 1172bc8ff85bSBarry Smith 117330dcb7c9SBarry Smith /* for each nontrival globally owned node list all arriving processors */ 1174bc8ff85bSBarry Smith for (i=0; i<nrecvs; i++) { 1175bc8ff85bSBarry Smith for (j=0; j<len[i]; j++) { 117630dcb7c9SBarry Smith node = recvs[2*i*nmax+2*j]-rstart; 1177f6e5521dSKarl Rupp if (nownedsenders[node] > 1) ownedsenders[starts[node]++] = source[i]; 1178bc8ff85bSBarry Smith } 1179bc8ff85bSBarry Smith } 1180bc8ff85bSBarry Smith 118107b52d57SBarry Smith if (debug) { /* ----------------------------------- */ 118230dcb7c9SBarry Smith starts[0] = 0; 118330dcb7c9SBarry Smith for (i=1; i<ng; i++) { 118430dcb7c9SBarry Smith if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1]; 118530dcb7c9SBarry Smith else starts[i] = starts[i-1]; 118630dcb7c9SBarry Smith } 118730dcb7c9SBarry Smith for (i=0; i<ng; i++) { 118830dcb7c9SBarry Smith if (nownedsenders[i] > 1) { 11899566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm,"[%d] global node %" PetscInt_FMT " local owner processors: ",rank,i+rstart)); 119030dcb7c9SBarry Smith for (j=0; j<nownedsenders[i]; j++) { 11919566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm,"%" PetscInt_FMT " ",ownedsenders[starts[i]+j])); 119230dcb7c9SBarry Smith } 11939566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm,"\n")); 119430dcb7c9SBarry Smith } 119530dcb7c9SBarry Smith } 11969566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(comm,PETSC_STDOUT)); 119707b52d57SBarry Smith } /* ----------------------------------- */ 119830dcb7c9SBarry Smith 11993677ff5aSBarry Smith /* wait on original sends */ 12003a96401aSBarry Smith if (nsends) { 12019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends,&send_status)); 12029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nsends,send_waits,send_status)); 12039566063dSJacob Faibussowitsch PetscCall(PetscFree(send_status)); 12043a96401aSBarry Smith } 12059566063dSJacob Faibussowitsch PetscCall(PetscFree(send_waits)); 12069566063dSJacob Faibussowitsch PetscCall(PetscFree(sends)); 12079566063dSJacob Faibussowitsch PetscCall(PetscFree(nprocs)); 12083677ff5aSBarry Smith 12093677ff5aSBarry Smith /* pack messages to send back to local owners */ 121030dcb7c9SBarry Smith starts[0] = 0; 121130dcb7c9SBarry Smith for (i=1; i<ng; i++) { 121230dcb7c9SBarry Smith if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1]; 121330dcb7c9SBarry Smith else starts[i] = starts[i-1]; 121430dcb7c9SBarry Smith } 121530dcb7c9SBarry Smith nsends2 = nrecvs; 12169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends2+1,&nprocs)); /* length of each message */ 121730dcb7c9SBarry Smith for (i=0; i<nrecvs; i++) { 121830dcb7c9SBarry Smith nprocs[i] = 1; 121930dcb7c9SBarry Smith for (j=0; j<len[i]; j++) { 122030dcb7c9SBarry Smith node = recvs[2*i*nmax+2*j]-rstart; 1221f6e5521dSKarl Rupp if (nownedsenders[node] > 1) nprocs[i] += 2 + nownedsenders[node]; 122230dcb7c9SBarry Smith } 122330dcb7c9SBarry Smith } 1224f6e5521dSKarl Rupp nt = 0; 1225f6e5521dSKarl Rupp for (i=0; i<nsends2; i++) nt += nprocs[i]; 1226f6e5521dSKarl Rupp 12279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nt+1,&sends2)); 12289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends2+1,&starts2)); 1229f6e5521dSKarl Rupp 1230f6e5521dSKarl Rupp starts2[0] = 0; 1231f6e5521dSKarl Rupp for (i=1; i<nsends2; i++) starts2[i] = starts2[i-1] + nprocs[i-1]; 123230dcb7c9SBarry Smith /* 123330dcb7c9SBarry Smith Each message is 1 + nprocs[i] long, and consists of 123430dcb7c9SBarry Smith (0) the number of nodes being sent back 123530dcb7c9SBarry Smith (1) the local node number, 123630dcb7c9SBarry Smith (2) the number of processors sharing it, 123730dcb7c9SBarry Smith (3) the processors sharing it 123830dcb7c9SBarry Smith */ 123930dcb7c9SBarry Smith for (i=0; i<nsends2; i++) { 124030dcb7c9SBarry Smith cnt = 1; 124130dcb7c9SBarry Smith sends2[starts2[i]] = 0; 124230dcb7c9SBarry Smith for (j=0; j<len[i]; j++) { 124330dcb7c9SBarry Smith node = recvs[2*i*nmax+2*j]-rstart; 124430dcb7c9SBarry Smith if (nownedsenders[node] > 1) { 124530dcb7c9SBarry Smith sends2[starts2[i]]++; 124630dcb7c9SBarry Smith sends2[starts2[i]+cnt++] = recvs[2*i*nmax+2*j+1]; 124730dcb7c9SBarry Smith sends2[starts2[i]+cnt++] = nownedsenders[node]; 12489566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(&sends2[starts2[i]+cnt],&ownedsenders[starts[node]],nownedsenders[node])); 124930dcb7c9SBarry Smith cnt += nownedsenders[node]; 125030dcb7c9SBarry Smith } 125130dcb7c9SBarry Smith } 125230dcb7c9SBarry Smith } 125330dcb7c9SBarry Smith 125430dcb7c9SBarry Smith /* receive the message lengths */ 125530dcb7c9SBarry Smith nrecvs2 = nsends; 12569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2+1,&lens2)); 12579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2+1,&starts3)); 12589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2+1,&recv_waits)); 125930dcb7c9SBarry Smith for (i=0; i<nrecvs2; i++) { 12609566063dSJacob Faibussowitsch PetscCallMPI(MPI_Irecv(&lens2[i],1,MPIU_INT,dest[i],tag2,comm,recv_waits+i)); 126130dcb7c9SBarry Smith } 1262d44834fbSBarry Smith 12638a8e0b3aSBarry Smith /* send the message lengths */ 12648a8e0b3aSBarry Smith for (i=0; i<nsends2; i++) { 12659566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(&nprocs[i],1,MPIU_INT,source[i],tag2,comm)); 12668a8e0b3aSBarry Smith } 12678a8e0b3aSBarry Smith 1268d44834fbSBarry Smith /* wait on receives of lens */ 12690c468ba9SBarry Smith if (nrecvs2) { 12709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2,&recv_statuses)); 12719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrecvs2,recv_waits,recv_statuses)); 12729566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_statuses)); 12730c468ba9SBarry Smith } 12749566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_waits)); 1275d44834fbSBarry Smith 127630dcb7c9SBarry Smith starts3[0] = 0; 1277d44834fbSBarry Smith nt = 0; 127830dcb7c9SBarry Smith for (i=0; i<nrecvs2-1; i++) { 127930dcb7c9SBarry Smith starts3[i+1] = starts3[i] + lens2[i]; 1280d44834fbSBarry Smith nt += lens2[i]; 128130dcb7c9SBarry Smith } 128276466f69SStefano Zampini if (nrecvs2) nt += lens2[nrecvs2-1]; 1283d44834fbSBarry Smith 12849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nt+1,&recvs2)); 12859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2+1,&recv_waits)); 128652b72c4aSBarry Smith for (i=0; i<nrecvs2; i++) { 12879566063dSJacob Faibussowitsch PetscCallMPI(MPI_Irecv(recvs2+starts3[i],lens2[i],MPIU_INT,dest[i],tag3,comm,recv_waits+i)); 128830dcb7c9SBarry Smith } 128930dcb7c9SBarry Smith 129030dcb7c9SBarry Smith /* send the messages */ 12919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends2+1,&send_waits)); 129230dcb7c9SBarry Smith for (i=0; i<nsends2; i++) { 12939566063dSJacob Faibussowitsch PetscCallMPI(MPI_Isend(sends2+starts2[i],nprocs[i],MPIU_INT,source[i],tag3,comm,send_waits+i)); 129430dcb7c9SBarry Smith } 129530dcb7c9SBarry Smith 129630dcb7c9SBarry Smith /* wait on receives */ 12970c468ba9SBarry Smith if (nrecvs2) { 12989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrecvs2,&recv_statuses)); 12999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrecvs2,recv_waits,recv_statuses)); 13009566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_statuses)); 13010c468ba9SBarry Smith } 13029566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_waits)); 13039566063dSJacob Faibussowitsch PetscCall(PetscFree(nprocs)); 130430dcb7c9SBarry Smith 130507b52d57SBarry Smith if (debug) { /* ----------------------------------- */ 130630dcb7c9SBarry Smith cnt = 0; 130730dcb7c9SBarry Smith for (i=0; i<nrecvs2; i++) { 130830dcb7c9SBarry Smith nt = recvs2[cnt++]; 130930dcb7c9SBarry Smith for (j=0; j<nt; j++) { 13109566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm,"[%d] local node %" PetscInt_FMT " number of subdomains %" PetscInt_FMT ": ",rank,recvs2[cnt],recvs2[cnt+1])); 131130dcb7c9SBarry Smith for (k=0; k<recvs2[cnt+1]; k++) { 13129566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm,"%" PetscInt_FMT " ",recvs2[cnt+2+k])); 131330dcb7c9SBarry Smith } 131430dcb7c9SBarry Smith cnt += 2 + recvs2[cnt+1]; 13159566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm,"\n")); 131630dcb7c9SBarry Smith } 131730dcb7c9SBarry Smith } 13189566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(comm,PETSC_STDOUT)); 131907b52d57SBarry Smith } /* ----------------------------------- */ 132030dcb7c9SBarry Smith 132130dcb7c9SBarry Smith /* count number subdomains for each local node */ 13229566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(size,&nprocs)); 132330dcb7c9SBarry Smith cnt = 0; 132430dcb7c9SBarry Smith for (i=0; i<nrecvs2; i++) { 132530dcb7c9SBarry Smith nt = recvs2[cnt++]; 132630dcb7c9SBarry Smith for (j=0; j<nt; j++) { 1327f6e5521dSKarl Rupp for (k=0; k<recvs2[cnt+1]; k++) nprocs[recvs2[cnt+2+k]]++; 132830dcb7c9SBarry Smith cnt += 2 + recvs2[cnt+1]; 132930dcb7c9SBarry Smith } 133030dcb7c9SBarry Smith } 133130dcb7c9SBarry Smith nt = 0; for (i=0; i<size; i++) nt += (nprocs[i] > 0); 133230dcb7c9SBarry Smith *nproc = nt; 13339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nt+1,procs)); 13349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nt+1,numprocs)); 13359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nt+1,indices)); 13360298fd71SBarry Smith for (i=0;i<nt+1;i++) (*indices)[i]=NULL; 13379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size,&bprocs)); 133830dcb7c9SBarry Smith cnt = 0; 133930dcb7c9SBarry Smith for (i=0; i<size; i++) { 134030dcb7c9SBarry Smith if (nprocs[i] > 0) { 134130dcb7c9SBarry Smith bprocs[i] = cnt; 134230dcb7c9SBarry Smith (*procs)[cnt] = i; 134330dcb7c9SBarry Smith (*numprocs)[cnt] = nprocs[i]; 13449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nprocs[i],&(*indices)[cnt])); 134530dcb7c9SBarry Smith cnt++; 134630dcb7c9SBarry Smith } 134730dcb7c9SBarry Smith } 134830dcb7c9SBarry Smith 134930dcb7c9SBarry Smith /* make the list of subdomains for each nontrivial local node */ 13509566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(*numprocs,nt)); 135130dcb7c9SBarry Smith cnt = 0; 135230dcb7c9SBarry Smith for (i=0; i<nrecvs2; i++) { 135330dcb7c9SBarry Smith nt = recvs2[cnt++]; 135430dcb7c9SBarry Smith for (j=0; j<nt; j++) { 1355f6e5521dSKarl Rupp for (k=0; k<recvs2[cnt+1]; k++) (*indices)[bprocs[recvs2[cnt+2+k]]][(*numprocs)[bprocs[recvs2[cnt+2+k]]]++] = recvs2[cnt]; 135630dcb7c9SBarry Smith cnt += 2 + recvs2[cnt+1]; 135730dcb7c9SBarry Smith } 135830dcb7c9SBarry Smith } 13599566063dSJacob Faibussowitsch PetscCall(PetscFree(bprocs)); 13609566063dSJacob Faibussowitsch PetscCall(PetscFree(recvs2)); 136130dcb7c9SBarry Smith 136207b52d57SBarry Smith /* sort the node indexing by their global numbers */ 136307b52d57SBarry Smith nt = *nproc; 136407b52d57SBarry Smith for (i=0; i<nt; i++) { 13659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((*numprocs)[i],&tmp)); 1366f6e5521dSKarl Rupp for (j=0; j<(*numprocs)[i]; j++) tmp[j] = lindices[(*indices)[i][j]]; 13679566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithArray((*numprocs)[i],tmp,(*indices)[i])); 13689566063dSJacob Faibussowitsch PetscCall(PetscFree(tmp)); 136907b52d57SBarry Smith } 137007b52d57SBarry Smith 137107b52d57SBarry Smith if (debug) { /* ----------------------------------- */ 137230dcb7c9SBarry Smith nt = *nproc; 137330dcb7c9SBarry Smith for (i=0; i<nt; i++) { 13749566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm,"[%d] subdomain %" PetscInt_FMT " number of indices %" PetscInt_FMT ": ",rank,(*procs)[i],(*numprocs)[i])); 137530dcb7c9SBarry Smith for (j=0; j<(*numprocs)[i]; j++) { 13769566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm,"%" PetscInt_FMT " ",(*indices)[i][j])); 137730dcb7c9SBarry Smith } 13789566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm,"\n")); 137930dcb7c9SBarry Smith } 13809566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(comm,PETSC_STDOUT)); 138107b52d57SBarry Smith } /* ----------------------------------- */ 138230dcb7c9SBarry Smith 138330dcb7c9SBarry Smith /* wait on sends */ 138430dcb7c9SBarry Smith if (nsends2) { 13859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsends2,&send_status)); 13869566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nsends2,send_waits,send_status)); 13879566063dSJacob Faibussowitsch PetscCall(PetscFree(send_status)); 138830dcb7c9SBarry Smith } 138930dcb7c9SBarry Smith 13909566063dSJacob Faibussowitsch PetscCall(PetscFree(starts3)); 13919566063dSJacob Faibussowitsch PetscCall(PetscFree(dest)); 13929566063dSJacob Faibussowitsch PetscCall(PetscFree(send_waits)); 13933677ff5aSBarry Smith 13949566063dSJacob Faibussowitsch PetscCall(PetscFree(nownedsenders)); 13959566063dSJacob Faibussowitsch PetscCall(PetscFree(ownedsenders)); 13969566063dSJacob Faibussowitsch PetscCall(PetscFree(starts)); 13979566063dSJacob Faibussowitsch PetscCall(PetscFree(starts2)); 13989566063dSJacob Faibussowitsch PetscCall(PetscFree(lens2)); 139989d82c54SBarry Smith 14009566063dSJacob Faibussowitsch PetscCall(PetscFree(source)); 14019566063dSJacob Faibussowitsch PetscCall(PetscFree(len)); 14029566063dSJacob Faibussowitsch PetscCall(PetscFree(recvs)); 14039566063dSJacob Faibussowitsch PetscCall(PetscFree(nprocs)); 14049566063dSJacob Faibussowitsch PetscCall(PetscFree(sends2)); 140524cf384cSBarry Smith 140624cf384cSBarry Smith /* put the information about myself as the first entry in the list */ 140724cf384cSBarry Smith first_procs = (*procs)[0]; 140824cf384cSBarry Smith first_numprocs = (*numprocs)[0]; 140924cf384cSBarry Smith first_indices = (*indices)[0]; 141024cf384cSBarry Smith for (i=0; i<*nproc; i++) { 141124cf384cSBarry Smith if ((*procs)[i] == rank) { 141224cf384cSBarry Smith (*procs)[0] = (*procs)[i]; 141324cf384cSBarry Smith (*numprocs)[0] = (*numprocs)[i]; 141424cf384cSBarry Smith (*indices)[0] = (*indices)[i]; 141524cf384cSBarry Smith (*procs)[i] = first_procs; 141624cf384cSBarry Smith (*numprocs)[i] = first_numprocs; 141724cf384cSBarry Smith (*indices)[i] = first_indices; 141824cf384cSBarry Smith break; 141924cf384cSBarry Smith } 142024cf384cSBarry Smith } 1421268a049cSStefano Zampini 1422268a049cSStefano Zampini /* save info for reuse */ 1423268a049cSStefano Zampini mapping->info_nproc = *nproc; 1424268a049cSStefano Zampini mapping->info_procs = *procs; 1425268a049cSStefano Zampini mapping->info_numprocs = *numprocs; 1426268a049cSStefano Zampini mapping->info_indices = *indices; 1427268a049cSStefano Zampini mapping->info_cached = PETSC_TRUE; 142889d82c54SBarry Smith PetscFunctionReturn(0); 142989d82c54SBarry Smith } 143089d82c54SBarry Smith 14316a818285SBarry Smith /*@C 14326a818285SBarry Smith ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by ISLocalToGlobalMappingGetBlockInfo() 14336a818285SBarry Smith 14346a818285SBarry Smith Collective on ISLocalToGlobalMapping 14356a818285SBarry Smith 1436f899ff85SJose E. Roman Input Parameter: 14376a818285SBarry Smith . mapping - the mapping from local to global indexing 14386a818285SBarry Smith 1439d8d19677SJose E. Roman Output Parameters: 14406a818285SBarry Smith + nproc - number of processors that are connected to this one 14416a818285SBarry Smith . proc - neighboring processors 14426a818285SBarry Smith . numproc - number of indices for each processor 14436a818285SBarry Smith - indices - indices of local nodes shared with neighbor (sorted by global numbering) 14446a818285SBarry Smith 14456a818285SBarry Smith Level: advanced 14466a818285SBarry Smith 1447db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 1448db781477SPatrick Sanan `ISLocalToGlobalMappingGetInfo()` 14496a818285SBarry Smith @*/ 14506a818285SBarry Smith PetscErrorCode ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[]) 14516a818285SBarry Smith { 14526a818285SBarry Smith PetscFunctionBegin; 1453cbc1caf0SMatthew G. Knepley PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 1454268a049cSStefano Zampini if (mapping->info_free) { 14559566063dSJacob Faibussowitsch PetscCall(PetscFree(*numprocs)); 14566a818285SBarry Smith if (*indices) { 1457268a049cSStefano Zampini PetscInt i; 1458268a049cSStefano Zampini 14599566063dSJacob Faibussowitsch PetscCall(PetscFree((*indices)[0])); 14606a818285SBarry Smith for (i=1; i<*nproc; i++) { 14619566063dSJacob Faibussowitsch PetscCall(PetscFree((*indices)[i])); 14626a818285SBarry Smith } 14639566063dSJacob Faibussowitsch PetscCall(PetscFree(*indices)); 14646a818285SBarry Smith } 1465268a049cSStefano Zampini } 1466268a049cSStefano Zampini *nproc = 0; 1467268a049cSStefano Zampini *procs = NULL; 1468268a049cSStefano Zampini *numprocs = NULL; 1469268a049cSStefano Zampini *indices = NULL; 14706a818285SBarry Smith PetscFunctionReturn(0); 14716a818285SBarry Smith } 14726a818285SBarry Smith 14736a818285SBarry Smith /*@C 14746a818285SBarry Smith ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each processor and 14756a818285SBarry Smith each index shared by more than one processor 14766a818285SBarry Smith 14776a818285SBarry Smith Collective on ISLocalToGlobalMapping 14786a818285SBarry Smith 1479f899ff85SJose E. Roman Input Parameter: 14806a818285SBarry Smith . mapping - the mapping from local to global indexing 14816a818285SBarry Smith 1482d8d19677SJose E. Roman Output Parameters: 14836a818285SBarry Smith + nproc - number of processors that are connected to this one 14846a818285SBarry Smith . proc - neighboring processors 14856a818285SBarry Smith . numproc - number of indices for each subdomain (processor) 14866a818285SBarry Smith - indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering) 14876a818285SBarry Smith 14886a818285SBarry Smith Level: advanced 14896a818285SBarry Smith 14901bd0b88eSStefano Zampini Notes: The user needs to call ISLocalToGlobalMappingRestoreInfo when the data is no longer needed. 14911bd0b88eSStefano Zampini 14926a818285SBarry Smith Fortran Usage: 14936a818285SBarry Smith $ ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by 14946a818285SBarry Smith $ ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc], 14956a818285SBarry Smith PetscInt indices[nproc][numprocmax],ierr) 14966a818285SBarry Smith There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and 14976a818285SBarry Smith indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough. 14986a818285SBarry Smith 1499db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 1500db781477SPatrick Sanan `ISLocalToGlobalMappingRestoreInfo()` 15016a818285SBarry Smith @*/ 15026a818285SBarry Smith PetscErrorCode ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[]) 15036a818285SBarry Smith { 15048a1f772fSStefano Zampini PetscInt **bindices = NULL,*bnumprocs = NULL,bs,i,j,k; 15056a818285SBarry Smith 15066a818285SBarry Smith PetscFunctionBegin; 15076a818285SBarry Smith PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 15088a1f772fSStefano Zampini bs = mapping->bs; 15099566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockInfo(mapping,nproc,procs,&bnumprocs,&bindices)); 1510268a049cSStefano Zampini if (bs > 1) { /* we need to expand the cached info */ 15119566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(*nproc,&*indices)); 15129566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(*nproc,&*numprocs)); 15136a818285SBarry Smith for (i=0; i<*nproc; i++) { 15149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs*bnumprocs[i],&(*indices)[i])); 1515268a049cSStefano Zampini for (j=0; j<bnumprocs[i]; j++) { 15166a818285SBarry Smith for (k=0; k<bs; k++) { 15176a818285SBarry Smith (*indices)[i][j*bs+k] = bs*bindices[i][j] + k; 15186a818285SBarry Smith } 15196a818285SBarry Smith } 1520268a049cSStefano Zampini (*numprocs)[i] = bnumprocs[i]*bs; 15216a818285SBarry Smith } 1522268a049cSStefano Zampini mapping->info_free = PETSC_TRUE; 1523268a049cSStefano Zampini } else { 1524268a049cSStefano Zampini *numprocs = bnumprocs; 1525268a049cSStefano Zampini *indices = bindices; 15266a818285SBarry Smith } 15276a818285SBarry Smith PetscFunctionReturn(0); 15286a818285SBarry Smith } 15296a818285SBarry Smith 153007b52d57SBarry Smith /*@C 153107b52d57SBarry Smith ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by ISLocalToGlobalMappingGetInfo() 153289d82c54SBarry Smith 153307b52d57SBarry Smith Collective on ISLocalToGlobalMapping 153407b52d57SBarry Smith 1535f899ff85SJose E. Roman Input Parameter: 153607b52d57SBarry Smith . mapping - the mapping from local to global indexing 153707b52d57SBarry Smith 1538d8d19677SJose E. Roman Output Parameters: 153907b52d57SBarry Smith + nproc - number of processors that are connected to this one 154007b52d57SBarry Smith . proc - neighboring processors 154107b52d57SBarry Smith . numproc - number of indices for each processor 154207b52d57SBarry Smith - indices - indices of local nodes shared with neighbor (sorted by global numbering) 154307b52d57SBarry Smith 154407b52d57SBarry Smith Level: advanced 154507b52d57SBarry Smith 1546db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 1547db781477SPatrick Sanan `ISLocalToGlobalMappingGetInfo()` 154807b52d57SBarry Smith @*/ 15497087cfbeSBarry Smith PetscErrorCode ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[]) 155007b52d57SBarry Smith { 155107b52d57SBarry Smith PetscFunctionBegin; 15529566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreBlockInfo(mapping,nproc,procs,numprocs,indices)); 155307b52d57SBarry Smith PetscFunctionReturn(0); 155407b52d57SBarry Smith } 155586994e45SJed Brown 155686994e45SJed Brown /*@C 15571bd0b88eSStefano Zampini ISLocalToGlobalMappingGetNodeInfo - Gets the neighbor information for each node 15581bd0b88eSStefano Zampini 15591bd0b88eSStefano Zampini Collective on ISLocalToGlobalMapping 15601bd0b88eSStefano Zampini 1561f899ff85SJose E. Roman Input Parameter: 15621bd0b88eSStefano Zampini . mapping - the mapping from local to global indexing 15631bd0b88eSStefano Zampini 1564d8d19677SJose E. Roman Output Parameters: 15651bd0b88eSStefano Zampini + nnodes - number of local nodes (same ISLocalToGlobalMappingGetSize()) 15661bd0b88eSStefano Zampini . count - number of neighboring processors per node 15671bd0b88eSStefano Zampini - indices - indices of processes sharing the node (sorted) 15681bd0b88eSStefano Zampini 15691bd0b88eSStefano Zampini Level: advanced 15701bd0b88eSStefano Zampini 15711bd0b88eSStefano Zampini Notes: The user needs to call ISLocalToGlobalMappingRestoreInfo when the data is no longer needed. 15721bd0b88eSStefano Zampini 1573db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 1574db781477SPatrick Sanan `ISLocalToGlobalMappingGetInfo()`, `ISLocalToGlobalMappingRestoreNodeInfo()` 15751bd0b88eSStefano Zampini @*/ 15761bd0b88eSStefano Zampini PetscErrorCode ISLocalToGlobalMappingGetNodeInfo(ISLocalToGlobalMapping mapping,PetscInt *nnodes,PetscInt *count[],PetscInt **indices[]) 15771bd0b88eSStefano Zampini { 15781bd0b88eSStefano Zampini PetscInt n; 15791bd0b88eSStefano Zampini 15801bd0b88eSStefano Zampini PetscFunctionBegin; 15811bd0b88eSStefano Zampini PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 15829566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(mapping,&n)); 15831bd0b88eSStefano Zampini if (!mapping->info_nodec) { 15841bd0b88eSStefano Zampini PetscInt i,m,n_neigh,*neigh,*n_shared,**shared; 15851bd0b88eSStefano Zampini 15869566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n+1,&mapping->info_nodec,n,&mapping->info_nodei)); 15879566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetInfo(mapping,&n_neigh,&neigh,&n_shared,&shared)); 1588071fcb05SBarry Smith for (i=0;i<n;i++) { mapping->info_nodec[i] = 1;} 1589071fcb05SBarry Smith m = n; 1590071fcb05SBarry Smith mapping->info_nodec[n] = 0; 15911bd0b88eSStefano Zampini for (i=1;i<n_neigh;i++) { 15921bd0b88eSStefano Zampini PetscInt j; 15931bd0b88eSStefano Zampini 15941bd0b88eSStefano Zampini m += n_shared[i]; 15951bd0b88eSStefano Zampini for (j=0;j<n_shared[i];j++) mapping->info_nodec[shared[i][j]] += 1; 15961bd0b88eSStefano Zampini } 15979566063dSJacob Faibussowitsch if (n) PetscCall(PetscMalloc1(m,&mapping->info_nodei[0])); 15981bd0b88eSStefano Zampini for (i=1;i<n;i++) mapping->info_nodei[i] = mapping->info_nodei[i-1] + mapping->info_nodec[i-1]; 15999566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(mapping->info_nodec,n)); 16001bd0b88eSStefano Zampini for (i=0;i<n;i++) { mapping->info_nodec[i] = 1; mapping->info_nodei[i][0] = neigh[0]; } 16011bd0b88eSStefano Zampini for (i=1;i<n_neigh;i++) { 16021bd0b88eSStefano Zampini PetscInt j; 16031bd0b88eSStefano Zampini 16041bd0b88eSStefano Zampini for (j=0;j<n_shared[i];j++) { 16051bd0b88eSStefano Zampini PetscInt k = shared[i][j]; 16061bd0b88eSStefano Zampini 16071bd0b88eSStefano Zampini mapping->info_nodei[k][mapping->info_nodec[k]] = neigh[i]; 16081bd0b88eSStefano Zampini mapping->info_nodec[k] += 1; 16091bd0b88eSStefano Zampini } 16101bd0b88eSStefano Zampini } 16119566063dSJacob Faibussowitsch for (i=0;i<n;i++) PetscCall(PetscSortRemoveDupsInt(&mapping->info_nodec[i],mapping->info_nodei[i])); 16129566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreInfo(mapping,&n_neigh,&neigh,&n_shared,&shared)); 16131bd0b88eSStefano Zampini } 16141bd0b88eSStefano Zampini if (nnodes) *nnodes = n; 16151bd0b88eSStefano Zampini if (count) *count = mapping->info_nodec; 16161bd0b88eSStefano Zampini if (indices) *indices = mapping->info_nodei; 16171bd0b88eSStefano Zampini PetscFunctionReturn(0); 16181bd0b88eSStefano Zampini } 16191bd0b88eSStefano Zampini 16201bd0b88eSStefano Zampini /*@C 16211bd0b88eSStefano Zampini ISLocalToGlobalMappingRestoreNodeInfo - Frees the memory allocated by ISLocalToGlobalMappingGetNodeInfo() 16221bd0b88eSStefano Zampini 16231bd0b88eSStefano Zampini Collective on ISLocalToGlobalMapping 16241bd0b88eSStefano Zampini 1625f899ff85SJose E. Roman Input Parameter: 16261bd0b88eSStefano Zampini . mapping - the mapping from local to global indexing 16271bd0b88eSStefano Zampini 1628d8d19677SJose E. Roman Output Parameters: 16291bd0b88eSStefano Zampini + nnodes - number of local nodes 16301bd0b88eSStefano Zampini . count - number of neighboring processors per node 16311bd0b88eSStefano Zampini - indices - indices of processes sharing the node (sorted) 16321bd0b88eSStefano Zampini 16331bd0b88eSStefano Zampini Level: advanced 16341bd0b88eSStefano Zampini 1635db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingDestroy()`, `ISLocalToGlobalMappingCreateIS()`, `ISLocalToGlobalMappingCreate()`, 1636db781477SPatrick Sanan `ISLocalToGlobalMappingGetInfo()` 16371bd0b88eSStefano Zampini @*/ 16381bd0b88eSStefano Zampini PetscErrorCode ISLocalToGlobalMappingRestoreNodeInfo(ISLocalToGlobalMapping mapping,PetscInt *nnodes,PetscInt *count[],PetscInt **indices[]) 16391bd0b88eSStefano Zampini { 16401bd0b88eSStefano Zampini PetscFunctionBegin; 16411bd0b88eSStefano Zampini PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 16421bd0b88eSStefano Zampini if (nnodes) *nnodes = 0; 16431bd0b88eSStefano Zampini if (count) *count = NULL; 16441bd0b88eSStefano Zampini if (indices) *indices = NULL; 16451bd0b88eSStefano Zampini PetscFunctionReturn(0); 16461bd0b88eSStefano Zampini } 16471bd0b88eSStefano Zampini 16481bd0b88eSStefano Zampini /*@C 1649107e9a97SBarry Smith ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped 165086994e45SJed Brown 165186994e45SJed Brown Not Collective 165286994e45SJed Brown 16534165533cSJose E. Roman Input Parameter: 165486994e45SJed Brown . ltog - local to global mapping 165586994e45SJed Brown 16564165533cSJose E. Roman Output Parameter: 1657565245c5SBarry Smith . array - array of indices, the length of this array may be obtained with ISLocalToGlobalMappingGetSize() 165886994e45SJed Brown 165986994e45SJed Brown Level: advanced 166086994e45SJed Brown 166195452b02SPatrick Sanan Notes: 166295452b02SPatrick Sanan ISLocalToGlobalMappingGetSize() returns the length the this array 1663107e9a97SBarry Smith 1664db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingRestoreIndices()`, `ISLocalToGlobalMappingGetBlockIndices()`, `ISLocalToGlobalMappingRestoreBlockIndices()` 166586994e45SJed Brown @*/ 16667087cfbeSBarry Smith PetscErrorCode ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog,const PetscInt **array) 166786994e45SJed Brown { 166886994e45SJed Brown PetscFunctionBegin; 166986994e45SJed Brown PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1); 167086994e45SJed Brown PetscValidPointer(array,2); 167145b6f7e9SBarry Smith if (ltog->bs == 1) { 167286994e45SJed Brown *array = ltog->indices; 167345b6f7e9SBarry Smith } else { 167445b6f7e9SBarry Smith PetscInt *jj,k,i,j,n = ltog->n, bs = ltog->bs; 167545b6f7e9SBarry Smith const PetscInt *ii; 167645b6f7e9SBarry Smith 16779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs*n,&jj)); 167845b6f7e9SBarry Smith *array = jj; 167945b6f7e9SBarry Smith k = 0; 168045b6f7e9SBarry Smith ii = ltog->indices; 168145b6f7e9SBarry Smith for (i=0; i<n; i++) 168245b6f7e9SBarry Smith for (j=0; j<bs; j++) 168345b6f7e9SBarry Smith jj[k++] = bs*ii[i] + j; 168445b6f7e9SBarry Smith } 168586994e45SJed Brown PetscFunctionReturn(0); 168686994e45SJed Brown } 168786994e45SJed Brown 168886994e45SJed Brown /*@C 1689193a2b41SJulian Andrej ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with ISLocalToGlobalMappingGetIndices() 169086994e45SJed Brown 169186994e45SJed Brown Not Collective 169286994e45SJed Brown 16934165533cSJose E. Roman Input Parameters: 169486994e45SJed Brown + ltog - local to global mapping 169586994e45SJed Brown - array - array of indices 169686994e45SJed Brown 169786994e45SJed Brown Level: advanced 169886994e45SJed Brown 1699db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingGetIndices()` 170086994e45SJed Brown @*/ 17017087cfbeSBarry Smith PetscErrorCode ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog,const PetscInt **array) 170286994e45SJed Brown { 170386994e45SJed Brown PetscFunctionBegin; 170486994e45SJed Brown PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1); 170586994e45SJed Brown PetscValidPointer(array,2); 1706c9cc58a2SBarry Smith PetscCheck(ltog->bs != 1 || *array == ltog->indices,PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer"); 170745b6f7e9SBarry Smith 170845b6f7e9SBarry Smith if (ltog->bs > 1) { 17099566063dSJacob Faibussowitsch PetscCall(PetscFree(*(void**)array)); 171045b6f7e9SBarry Smith } 171145b6f7e9SBarry Smith PetscFunctionReturn(0); 171245b6f7e9SBarry Smith } 171345b6f7e9SBarry Smith 171445b6f7e9SBarry Smith /*@C 171545b6f7e9SBarry Smith ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block 171645b6f7e9SBarry Smith 171745b6f7e9SBarry Smith Not Collective 171845b6f7e9SBarry Smith 17194165533cSJose E. Roman Input Parameter: 172045b6f7e9SBarry Smith . ltog - local to global mapping 172145b6f7e9SBarry Smith 17224165533cSJose E. Roman Output Parameter: 172345b6f7e9SBarry Smith . array - array of indices 172445b6f7e9SBarry Smith 172545b6f7e9SBarry Smith Level: advanced 172645b6f7e9SBarry Smith 1727db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingRestoreBlockIndices()` 172845b6f7e9SBarry Smith @*/ 172945b6f7e9SBarry Smith PetscErrorCode ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array) 173045b6f7e9SBarry Smith { 173145b6f7e9SBarry Smith PetscFunctionBegin; 173245b6f7e9SBarry Smith PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1); 173345b6f7e9SBarry Smith PetscValidPointer(array,2); 173445b6f7e9SBarry Smith *array = ltog->indices; 173545b6f7e9SBarry Smith PetscFunctionReturn(0); 173645b6f7e9SBarry Smith } 173745b6f7e9SBarry Smith 173845b6f7e9SBarry Smith /*@C 173945b6f7e9SBarry Smith ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with ISLocalToGlobalMappingGetBlockIndices() 174045b6f7e9SBarry Smith 174145b6f7e9SBarry Smith Not Collective 174245b6f7e9SBarry Smith 17434165533cSJose E. Roman Input Parameters: 174445b6f7e9SBarry Smith + ltog - local to global mapping 174545b6f7e9SBarry Smith - array - array of indices 174645b6f7e9SBarry Smith 174745b6f7e9SBarry Smith Level: advanced 174845b6f7e9SBarry Smith 1749db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingApply()`, `ISLocalToGlobalMappingGetIndices()` 175045b6f7e9SBarry Smith @*/ 175145b6f7e9SBarry Smith PetscErrorCode ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array) 175245b6f7e9SBarry Smith { 175345b6f7e9SBarry Smith PetscFunctionBegin; 175445b6f7e9SBarry Smith PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1); 175545b6f7e9SBarry Smith PetscValidPointer(array,2); 175608401ef6SPierre Jolivet PetscCheck(*array == ltog->indices,PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer"); 17570298fd71SBarry Smith *array = NULL; 175886994e45SJed Brown PetscFunctionReturn(0); 175986994e45SJed Brown } 1760f7efa3c7SJed Brown 1761f7efa3c7SJed Brown /*@C 1762f7efa3c7SJed Brown ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings 1763f7efa3c7SJed Brown 1764f7efa3c7SJed Brown Not Collective 1765f7efa3c7SJed Brown 17664165533cSJose E. Roman Input Parameters: 1767f7efa3c7SJed Brown + comm - communicator for the new mapping, must contain the communicator of every mapping to concatenate 1768f7efa3c7SJed Brown . n - number of mappings to concatenate 1769f7efa3c7SJed Brown - ltogs - local to global mappings 1770f7efa3c7SJed Brown 17714165533cSJose E. Roman Output Parameter: 1772f7efa3c7SJed Brown . ltogcat - new mapping 1773f7efa3c7SJed Brown 17749d90f715SBarry Smith Note: this currently always returns a mapping with block size of 1 17759d90f715SBarry Smith 17769d90f715SBarry Smith Developer Note: If all the input mapping have the same block size we could easily handle that as a special case 17779d90f715SBarry Smith 1778f7efa3c7SJed Brown Level: advanced 1779f7efa3c7SJed Brown 1780db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()` 1781f7efa3c7SJed Brown @*/ 1782f7efa3c7SJed Brown PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm,PetscInt n,const ISLocalToGlobalMapping ltogs[],ISLocalToGlobalMapping *ltogcat) 1783f7efa3c7SJed Brown { 1784f7efa3c7SJed Brown PetscInt i,cnt,m,*idx; 1785f7efa3c7SJed Brown 1786f7efa3c7SJed Brown PetscFunctionBegin; 178708401ef6SPierre Jolivet PetscCheck(n >= 0,comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have a non-negative number of mappings, given %" PetscInt_FMT,n); 1788f7efa3c7SJed Brown if (n > 0) PetscValidPointer(ltogs,3); 1789f7efa3c7SJed Brown for (i=0; i<n; i++) PetscValidHeaderSpecific(ltogs[i],IS_LTOGM_CLASSID,3); 1790f7efa3c7SJed Brown PetscValidPointer(ltogcat,4); 1791f7efa3c7SJed Brown for (cnt=0,i=0; i<n; i++) { 17929566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(ltogs[i],&m)); 1793f7efa3c7SJed Brown cnt += m; 1794f7efa3c7SJed Brown } 17959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(cnt,&idx)); 1796f7efa3c7SJed Brown for (cnt=0,i=0; i<n; i++) { 1797f7efa3c7SJed Brown const PetscInt *subidx; 17989566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(ltogs[i],&m)); 17999566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(ltogs[i],&subidx)); 18009566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(&idx[cnt],subidx,m)); 18019566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(ltogs[i],&subidx)); 1802f7efa3c7SJed Brown cnt += m; 1803f7efa3c7SJed Brown } 18049566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm,1,cnt,idx,PETSC_OWN_POINTER,ltogcat)); 1805f7efa3c7SJed Brown PetscFunctionReturn(0); 1806f7efa3c7SJed Brown } 180704a59952SBarry Smith 1808413f72f0SBarry Smith /*MC 1809413f72f0SBarry Smith ISLOCALTOGLOBALMAPPINGBASIC - basic implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is 1810413f72f0SBarry Smith used this is good for only small and moderate size problems. 1811413f72f0SBarry Smith 1812413f72f0SBarry Smith Options Database Keys: 1813a2b725a8SWilliam Gropp . -islocaltoglobalmapping_type basic - select this method 1814413f72f0SBarry Smith 1815413f72f0SBarry Smith Level: beginner 1816413f72f0SBarry Smith 1817db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`, `ISLOCALTOGLOBALMAPPINGHASH` 1818413f72f0SBarry Smith M*/ 1819413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Basic(ISLocalToGlobalMapping ltog) 1820413f72f0SBarry Smith { 1821413f72f0SBarry Smith PetscFunctionBegin; 1822413f72f0SBarry Smith ltog->ops->globaltolocalmappingapply = ISGlobalToLocalMappingApply_Basic; 1823413f72f0SBarry Smith ltog->ops->globaltolocalmappingsetup = ISGlobalToLocalMappingSetUp_Basic; 1824413f72f0SBarry Smith ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Basic; 1825413f72f0SBarry Smith ltog->ops->destroy = ISLocalToGlobalMappingDestroy_Basic; 1826413f72f0SBarry Smith PetscFunctionReturn(0); 1827413f72f0SBarry Smith } 1828413f72f0SBarry Smith 1829413f72f0SBarry Smith /*MC 1830413f72f0SBarry Smith ISLOCALTOGLOBALMAPPINGHASH - hash implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is 1831ed56e8eaSBarry Smith used this is good for large memory problems. 1832413f72f0SBarry Smith 1833413f72f0SBarry Smith Options Database Keys: 1834a2b725a8SWilliam Gropp . -islocaltoglobalmapping_type hash - select this method 1835413f72f0SBarry Smith 183695452b02SPatrick Sanan Notes: 183795452b02SPatrick Sanan This is selected automatically for large problems if the user does not set the type. 1838ed56e8eaSBarry Smith 1839413f72f0SBarry Smith Level: beginner 1840413f72f0SBarry Smith 1841db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()`, `ISLOCALTOGLOBALMAPPINGHASH` 1842413f72f0SBarry Smith M*/ 1843413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Hash(ISLocalToGlobalMapping ltog) 1844413f72f0SBarry Smith { 1845413f72f0SBarry Smith PetscFunctionBegin; 1846413f72f0SBarry Smith ltog->ops->globaltolocalmappingapply = ISGlobalToLocalMappingApply_Hash; 1847413f72f0SBarry Smith ltog->ops->globaltolocalmappingsetup = ISGlobalToLocalMappingSetUp_Hash; 1848413f72f0SBarry Smith ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Hash; 1849413f72f0SBarry Smith ltog->ops->destroy = ISLocalToGlobalMappingDestroy_Hash; 1850413f72f0SBarry Smith PetscFunctionReturn(0); 1851413f72f0SBarry Smith } 1852413f72f0SBarry Smith 1853413f72f0SBarry Smith /*@C 1854413f72f0SBarry Smith ISLocalToGlobalMappingRegister - Adds a method for applying a global to local mapping with an ISLocalToGlobalMapping 1855413f72f0SBarry Smith 1856413f72f0SBarry Smith Not Collective 1857413f72f0SBarry Smith 1858413f72f0SBarry Smith Input Parameters: 1859413f72f0SBarry Smith + sname - name of a new method 1860413f72f0SBarry Smith - routine_create - routine to create method context 1861413f72f0SBarry Smith 1862413f72f0SBarry Smith Notes: 1863ed56e8eaSBarry Smith ISLocalToGlobalMappingRegister() may be called multiple times to add several user-defined mappings. 1864413f72f0SBarry Smith 1865413f72f0SBarry Smith Sample usage: 1866413f72f0SBarry Smith .vb 1867413f72f0SBarry Smith ISLocalToGlobalMappingRegister("my_mapper",MyCreate); 1868413f72f0SBarry Smith .ve 1869413f72f0SBarry Smith 1870ed56e8eaSBarry Smith Then, your mapping can be chosen with the procedural interface via 1871413f72f0SBarry Smith $ ISLocalToGlobalMappingSetType(ltog,"my_mapper") 1872413f72f0SBarry Smith or at runtime via the option 1873ed56e8eaSBarry Smith $ -islocaltoglobalmapping_type my_mapper 1874413f72f0SBarry Smith 1875413f72f0SBarry Smith Level: advanced 1876413f72f0SBarry Smith 1877db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingRegisterAll()`, `ISLocalToGlobalMappingRegisterDestroy()`, `ISLOCALTOGLOBALMAPPINGBASIC`, `ISLOCALTOGLOBALMAPPINGHASH` 1878413f72f0SBarry Smith 1879413f72f0SBarry Smith @*/ 1880413f72f0SBarry Smith PetscErrorCode ISLocalToGlobalMappingRegister(const char sname[],PetscErrorCode (*function)(ISLocalToGlobalMapping)) 1881413f72f0SBarry Smith { 1882413f72f0SBarry Smith PetscFunctionBegin; 18839566063dSJacob Faibussowitsch PetscCall(ISInitializePackage()); 18849566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&ISLocalToGlobalMappingList,sname,function)); 1885413f72f0SBarry Smith PetscFunctionReturn(0); 1886413f72f0SBarry Smith } 1887413f72f0SBarry Smith 1888413f72f0SBarry Smith /*@C 1889ed56e8eaSBarry Smith ISLocalToGlobalMappingSetType - Builds ISLocalToGlobalMapping for a particular global to local mapping approach. 1890413f72f0SBarry Smith 1891413f72f0SBarry Smith Logically Collective on ISLocalToGlobalMapping 1892413f72f0SBarry Smith 1893413f72f0SBarry Smith Input Parameters: 1894413f72f0SBarry Smith + ltog - the ISLocalToGlobalMapping object 1895413f72f0SBarry Smith - type - a known method 1896413f72f0SBarry Smith 1897413f72f0SBarry Smith Options Database Key: 1898ed56e8eaSBarry Smith . -islocaltoglobalmapping_type <method> - Sets the method; use -help for a list 1899413f72f0SBarry Smith of available methods (for instance, basic or hash) 1900413f72f0SBarry Smith 1901413f72f0SBarry Smith Notes: 1902413f72f0SBarry Smith See "petsc/include/petscis.h" for available methods 1903413f72f0SBarry Smith 1904413f72f0SBarry Smith Normally, it is best to use the ISLocalToGlobalMappingSetFromOptions() command and 1905413f72f0SBarry Smith then set the ISLocalToGlobalMapping type from the options database rather than by using 1906413f72f0SBarry Smith this routine. 1907413f72f0SBarry Smith 1908413f72f0SBarry Smith Level: intermediate 1909413f72f0SBarry Smith 1910413f72f0SBarry Smith Developer Note: ISLocalToGlobalMappingRegister() is used to add new types to ISLocalToGlobalMappingList from which they 1911413f72f0SBarry Smith are accessed by ISLocalToGlobalMappingSetType(). 1912413f72f0SBarry Smith 1913db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingRegister()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingGetType()` 1914413f72f0SBarry Smith @*/ 1915413f72f0SBarry Smith PetscErrorCode ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType type) 1916413f72f0SBarry Smith { 1917413f72f0SBarry Smith PetscBool match; 19185f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(ISLocalToGlobalMapping) = NULL; 1919413f72f0SBarry Smith 1920413f72f0SBarry Smith PetscFunctionBegin; 1921413f72f0SBarry Smith PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1); 1922a0d79125SStefano Zampini if (type) PetscValidCharPointer(type,2); 1923413f72f0SBarry Smith 19249566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)ltog,type,&match)); 1925413f72f0SBarry Smith if (match) PetscFunctionReturn(0); 1926413f72f0SBarry Smith 1927a0d79125SStefano Zampini /* L2G maps defer type setup at globaltolocal calls, allow passing NULL here */ 1928a0d79125SStefano Zampini if (type) { 19299566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(ISLocalToGlobalMappingList,type,&r)); 1930a0d79125SStefano Zampini PetscCheck(r,PetscObjectComm((PetscObject)ltog),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested ISLocalToGlobalMapping type %s",type); 1931a0d79125SStefano Zampini } 1932413f72f0SBarry Smith /* Destroy the previous private LTOG context */ 1933413f72f0SBarry Smith if (ltog->ops->destroy) { 19349566063dSJacob Faibussowitsch PetscCall((*ltog->ops->destroy)(ltog)); 1935413f72f0SBarry Smith ltog->ops->destroy = NULL; 1936413f72f0SBarry Smith } 19379566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)ltog,type)); 19389566063dSJacob Faibussowitsch if (r) PetscCall((*r)(ltog)); 1939a0d79125SStefano Zampini PetscFunctionReturn(0); 1940a0d79125SStefano Zampini } 1941a0d79125SStefano Zampini 1942a0d79125SStefano Zampini /*@C 1943a0d79125SStefano Zampini ISLocalToGlobalMappingGetType - Get the type of the l2g map 1944a0d79125SStefano Zampini 1945a0d79125SStefano Zampini Not Collective 1946a0d79125SStefano Zampini 1947a0d79125SStefano Zampini Input Parameter: 1948a0d79125SStefano Zampini . ltog - the ISLocalToGlobalMapping object 1949a0d79125SStefano Zampini 1950a0d79125SStefano Zampini Output Parameter: 1951a0d79125SStefano Zampini . type - the type 1952a0d79125SStefano Zampini 195349762cbcSSatish Balay Level: intermediate 195449762cbcSSatish Balay 1955db781477SPatrick Sanan .seealso: `ISLocalToGlobalMappingType`, `ISLocalToGlobalMappingRegister()`, `ISLocalToGlobalMappingCreate()`, `ISLocalToGlobalMappingSetType()` 1956a0d79125SStefano Zampini @*/ 1957a0d79125SStefano Zampini PetscErrorCode ISLocalToGlobalMappingGetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType *type) 1958a0d79125SStefano Zampini { 1959a0d79125SStefano Zampini PetscFunctionBegin; 1960a0d79125SStefano Zampini PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1); 1961a0d79125SStefano Zampini PetscValidPointer(type,2); 1962a0d79125SStefano Zampini *type = ((PetscObject)ltog)->type_name; 1963413f72f0SBarry Smith PetscFunctionReturn(0); 1964413f72f0SBarry Smith } 1965413f72f0SBarry Smith 1966413f72f0SBarry Smith PetscBool ISLocalToGlobalMappingRegisterAllCalled = PETSC_FALSE; 1967413f72f0SBarry Smith 1968413f72f0SBarry Smith /*@C 1969413f72f0SBarry Smith ISLocalToGlobalMappingRegisterAll - Registers all of the local to global mapping components in the IS package. 1970413f72f0SBarry Smith 1971413f72f0SBarry Smith Not Collective 1972413f72f0SBarry Smith 1973413f72f0SBarry Smith Level: advanced 1974413f72f0SBarry Smith 1975db781477SPatrick Sanan .seealso: `ISRegister()`, `ISLocalToGlobalRegister()` 1976413f72f0SBarry Smith @*/ 1977413f72f0SBarry Smith PetscErrorCode ISLocalToGlobalMappingRegisterAll(void) 1978413f72f0SBarry Smith { 1979413f72f0SBarry Smith PetscFunctionBegin; 1980413f72f0SBarry Smith if (ISLocalToGlobalMappingRegisterAllCalled) PetscFunctionReturn(0); 1981413f72f0SBarry Smith ISLocalToGlobalMappingRegisterAllCalled = PETSC_TRUE; 19829566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGBASIC, ISLocalToGlobalMappingCreate_Basic)); 19839566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGHASH, ISLocalToGlobalMappingCreate_Hash)); 1984413f72f0SBarry Smith PetscFunctionReturn(0); 1985413f72f0SBarry Smith } 1986